]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliDCSSensorArray.cxx
Bug fix in copy constructor and assignement operator (Matthias + me)
[u/mrichter/AliRoot.git] / STEER / AliDCSSensorArray.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 //  Calibration class for DCS sensors                                        //
20 //  Authors: Marian Ivanov and Haavard Helstrup                              //
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23
24 #include "AliDCSSensorArray.h"
25 #include "AliLog.h"
26
27 ClassImp(AliDCSSensorArray)
28
29 const Double_t kSecInHour = 3600.; // seconds in one hour
30 const UInt_t   kMinMapTime = 60;   // don't fit maps shorter than one minute
31
32 //_____________________________________________________________________________
33 AliDCSSensorArray::AliDCSSensorArray():TNamed(), 
34   fMinGraph(10),
35   fMinPoints(10),
36   fIter(10),
37   fMaxDelta(0.0),
38   fFitReq(2),
39   fValCut(-1),
40   fDiffCut(-1),
41   fStartTime (2000,1,1,0,0,0),
42   fEndTime   (2000,1,1,0,0,0),
43   fSensors(0)
44 {
45   //
46   // AliDCSSensorArray default constructor
47   //
48
49 }
50 //_____________________________________________________________________________
51 AliDCSSensorArray::AliDCSSensorArray(TClonesArray *arr):TNamed(),
52   fMinGraph(10),
53   fMinPoints(10),
54   fIter(10),
55   fMaxDelta(0.0),
56   fFitReq(2),
57   fValCut(-1),
58   fDiffCut(-1),
59   fStartTime (2000,1,1,0,0,0),
60   fEndTime   (2000,1,1,0,0,0),
61   fSensors(arr)
62 {
63   //
64   // AliDCSSensorArray special constructor taking TClonesArray from ReadList
65   //
66
67 }
68 //_____________________________________________________________________________
69 AliDCSSensorArray::AliDCSSensorArray(Int_t run, const char* dbEntry) :
70   TNamed(),
71   fMinGraph(10),
72   fMinPoints(10),
73   fIter(10),
74   fMaxDelta(0.0),
75   fFitReq(2),
76   fValCut(-1),
77   fDiffCut(-1),
78   fStartTime (2000,1,1,0,0,0),
79   fEndTime   (2000,1,1,0,0,0),
80   fSensors(0)
81 {
82   //
83   // Read configuration from OCDB
84   //
85
86   AliCDBEntry *entry = AliCDBManager::Instance()->Get(dbEntry,run);
87   TTree *tree = (TTree*) entry->GetObject();
88   fSensors = AliDCSSensor::ReadTree(tree);
89 }
90 //_____________________________________________________________________________
91 AliDCSSensorArray::AliDCSSensorArray(UInt_t startTime, UInt_t endTime,
92                        TTree* confTree) :
93   TNamed(),
94   fMinGraph(10),
95   fMinPoints(10),
96   fIter(10),
97   fMaxDelta(0.0),
98   fFitReq(2),
99   fValCut(-1),
100   fDiffCut(-1),
101   fStartTime (2000,1,1,0,0,0),
102   fEndTime   (2000,1,1,0,0,0),
103   fSensors(0)
104
105 {
106   //
107   // AliDCSSensorArray constructor for Shuttle preprocessor
108   //  (confTree read from OCDB)
109   //
110   fSensors = AliDCSSensor::ReadTree(confTree);
111   fSensors->BypassStreamer(kFALSE);
112   fStartTime = TTimeStamp((time_t)startTime,0);
113   fEndTime   = TTimeStamp((time_t)endTime,0);
114 }
115
116
117 //_____________________________________________________________________________
118 AliDCSSensorArray::AliDCSSensorArray(UInt_t startTime, UInt_t endTime,
119                        TClonesArray *sensors) :
120   TNamed(),
121   fMinGraph(10),
122   fMinPoints(10),
123   fIter(10),
124   fMaxDelta(0.0),
125   fFitReq(2),
126   fValCut(-1),
127   fDiffCut(-1),
128   fStartTime (2000,1,1,0,0,0),
129   fEndTime   (2000,1,1,0,0,0),
130   fSensors(sensors)
131
132 {
133   //
134   // AliDCSSensorArray constructor for Shuttle preprocessor
135   //  (TClonesArray of AliDCSSensor objects)
136   //
137   fStartTime = TTimeStamp((time_t)startTime,0);
138   fEndTime   = TTimeStamp((time_t)endTime,0);
139 }
140
141 //_____________________________________________________________________________
142 AliDCSSensorArray::AliDCSSensorArray(const AliDCSSensorArray &c):TNamed(c),
143   fMinGraph(c.fMinGraph),
144   fMinPoints(c.fMinPoints),
145   fIter(c.fIter),
146   fMaxDelta(c.fMaxDelta),
147   fFitReq(c.fFitReq),
148   fValCut(c.fValCut),
149   fDiffCut(c.fDiffCut),
150   fStartTime (c.fStartTime),
151   fEndTime   (c.fEndTime),
152   fSensors(0)
153
154 {
155   //
156   // AliDCSSensorArray copy constructor
157   //
158
159   fSensors = (TClonesArray*)c.fSensors->Clone();
160 }
161
162 ///_____________________________________________________________________________
163 AliDCSSensorArray::~AliDCSSensorArray()
164 {
165   //
166   // AliDCSSensorArray destructor
167   //
168   fSensors->Delete();
169   delete fSensors;
170
171 }
172
173 //_____________________________________________________________________________
174 AliDCSSensorArray &AliDCSSensorArray::operator=(const AliDCSSensorArray &c)
175 {
176   //
177   // Assignment operator
178   //
179   if (this != &c) {
180      fSensors->Delete();
181      new (this) AliDCSSensorArray(c);
182      fSensors = (TClonesArray*)c.fSensors->Clone();
183   }
184   return *this;
185 }
186
187 //____________________________________________________________________________
188
189 void AliDCSSensorArray::SetGraph(TMap *map)
190 {
191   //
192   // Read graphs from DCS maps
193   //
194   Int_t nsensors = fSensors->GetEntries();
195   for ( Int_t isensor=0; isensor<nsensors; isensor++) {
196     AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
197     TString stringID = entry->GetStringID();
198     TGraph *gr = (TGraph*)map->GetValue(stringID.Data());
199     if ( gr !=0 ) {
200        entry->SetGraph((TGraph*)gr->Clone());
201     } else {
202        entry->SetGraph(0);
203     }
204   }
205 }
206 //_____________________________________________________________________________
207 void AliDCSSensorArray::MakeSplineFit(TMap *map, Bool_t keepMap)
208 {
209   //
210   // Make spline fits from DCS maps
211   //
212   Int_t nsensors = fSensors->GetEntries();
213   for ( Int_t isensor=0; isensor<nsensors; isensor++) {
214     AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
215     TString stringID = entry->GetStringID();
216     TGraph *gr = (TGraph*)map->GetValue(stringID.Data());
217     if (!gr ) {
218       entry->SetFit(0);
219       entry->SetGraph(0);
220       AliWarning(Form("sensor %s: no input graph",stringID.Data()));
221       continue;
222     }
223     UInt_t timeDiff = entry->GetEndTime() - entry->GetStartTime();
224     if ( timeDiff < kMinMapTime ) {
225       AliWarning(Form("sensor %s: map length < 60 s, DCS graph kept.",stringID.Data()));
226       entry->SetGraph((TGraph*)gr->Clone());
227     } else {
228       AliSplineFit *fit = new AliSplineFit();
229       fit->SetMinPoints(fMinGraph);
230       fit->InitKnots(gr,fMinPoints,fIter,fMaxDelta);
231       fit->SplineFit(fFitReq);
232       fit->Cleanup();
233       if (fit) {
234         entry->SetFit(fit);
235       } else {
236         AliWarning(Form("sensor %s: no fit performed, DCS graph kept.",stringID.Data()));
237         entry->SetGraph((TGraph*)gr->Clone());
238       }
239     }
240     if (keepMap) entry->SetGraph((TGraph*)gr->Clone());
241   }
242 }
243 //_____________________________________________________________________________
244 void AliDCSSensorArray::MakeSplineFitAddPoints(TMap *map)
245 {
246   //
247   // Make spline fits from DCS maps
248   //
249   Int_t nsensors = fSensors->GetEntries();
250   for ( Int_t isensor=0; isensor<nsensors; isensor++) {
251     AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
252
253   // fetch old points from existing graph
254
255     TGraph *gr = entry->GetGraph();
256     if (!gr) {
257       gr = new TGraph();
258       entry->SetGraph(gr);
259     } 
260     TString stringID = entry->GetStringID();
261
262   // fetch new points from DCS map
263   
264     TGraph *grAdd = (TGraph*)map->GetValue(stringID.Data());
265     if (!grAdd ) return;
266
267   // add new points to end of graph
268   
269     Int_t nPointsOld=gr->GetN();
270     Int_t nPointsAdd=grAdd->GetN();
271     gr->Expand(nPointsOld+nPointsAdd);
272     gr->Set(nPointsOld+nPointsAdd);
273     Double_t *addX=grAdd->GetX();
274     Double_t *addY=grAdd->GetY();
275     for (Int_t i=0;i<nPointsAdd;i++) {
276       gr->SetPoint(nPointsOld+i,addX[i],addY[i]);
277     }
278  
279    // make fit to complete graph
280    
281     AliSplineFit *fit = new AliSplineFit();
282     fit->SetMinPoints(fMinGraph);
283     fit->InitKnots(gr,fMinPoints,fIter,fMaxDelta);
284     fit->SplineFit(fFitReq);
285     fit->Cleanup();
286     if (fit) {
287       AliSplineFit *oldFit = entry->GetFit();
288       if (oldFit) delete oldFit;
289       entry->SetFit(fit);
290     } else {
291       AliWarning(Form("sensor %s: no new fit performed. If available, old fit kept.",stringID.Data()));
292     }
293   }
294 }
295
296 //_____________________________________________________________________________
297 Int_t AliDCSSensorArray::NumFits() const 
298 {
299  //
300  // Return number of sensors where a succesful fit has been made
301  //
302   Int_t nfit=0;
303   Int_t nsensors = fSensors->GetEntries();
304   for ( Int_t isensor=0; isensor<nsensors; isensor++) {
305     AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
306     if (entry->GetFit()) nfit++;
307   }    
308   return nfit;
309 }
310 //_____________________________________________________________________________
311 Double_t AliDCSSensorArray::GetValue(UInt_t timeSec, Int_t sensor) 
312 {
313   //
314   // Return sensor value at time timeSec (obtained from fitted function)
315   //  timeSec = time in seconds from start of run
316   //
317
318   AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(sensor);
319   return entry->GetValue(TTimeStamp((time_t)fStartTime.GetSec()+timeSec,0));
320 }
321
322
323 //_____________________________________________________________________________
324 TMap* AliDCSSensorArray::ExtractDCS(TMap *dcsMap, Bool_t keepStart)
325 {
326  //
327  // Extract temperature graphs from DCS maps
328  //
329  TMap *values = new TMap;
330  TObjArray * valueSet;
331  //
332  // Keep global start/end times
333  //    to avoid extrapolations, the fits will only be valid from first 
334  //    measured point to last measured point. This is consistent with hardware,
335  //    as there would be a new measured point if the value changed.
336  
337  TTimeStamp startTime=fStartTime;
338  TTimeStamp endTime=fEndTime;
339  
340  Int_t nsensors = fSensors->GetEntries();
341  for ( Int_t isensor=0; isensor<nsensors; isensor++) {
342    AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
343    TString stringID = entry->GetStringID();
344    TPair *pair = (TPair*)dcsMap->FindObject(stringID.Data());
345    if ( pair ) {                            // only try to read values
346                                             // if DCS object available
347      valueSet = (TObjArray*)pair->Value();
348      TGraph *graph = MakeGraph(valueSet,keepStart);   // MakeGraph sets start/end time
349                                             // per sensor
350      values->Add(new TObjString(stringID.Data()),graph);
351      entry->SetStartTime(fStartTime);
352      entry->SetEndTime(fEndTime);
353    }
354  }
355  // Reset global start/end time 
356  //    ..... yes, I know this won't get a prize for structured programming..:-)
357
358  fStartTime=startTime;
359  fEndTime=endTime;
360  return values;
361 }
362
363 //_____________________________________________________________________________
364 TGraph* AliDCSSensorArray::MakeGraph(TObjArray* valueSet, Bool_t keepStart){
365   //
366   // Make graph of temperature values read from DCS map
367   //   (spline fit parameters will subsequently be obtained from this graph) 
368   //
369   Int_t nentries = valueSet->GetEntriesFast(); 
370   if ( nentries == 0 ) return 0;
371   
372   Float_t *x = new Float_t[nentries];
373   Float_t *y = new Float_t[nentries];
374   Int_t time0=0;
375   TTimeStamp firstTime(0);
376   TTimeStamp lastTime(0);
377   if (keepStart) { 
378      firstTime = fStartTime;
379      time0 = firstTime.GetSec();
380   }
381   Int_t out=0;
382   Int_t skipped=0;
383   AliDCSValue *val = (AliDCSValue *)valueSet->At(0);
384   AliDCSValue::Type type = val->GetType();
385   if ( type == AliDCSValue::kInvalid || type == AliDCSValue::kBool ) return 0;
386   Float_t value;
387   for (Int_t i=0; i<nentries; i++){
388     val = (AliDCSValue *)valueSet->At(i);
389     if (!val) continue;
390     if (time0==0){
391       time0=val->GetTimeStamp();
392       firstTime= TTimeStamp((time_t)val->GetTimeStamp(),0);
393       lastTime=TTimeStamp((time_t)val->GetTimeStamp(),0);
394      }
395     switch ( type )
396     { 
397       case AliDCSValue::kFloat:
398         value = val->GetFloat();
399         break;
400       case AliDCSValue::kChar:
401         value = static_cast<Float_t>(val->GetChar());
402         break;
403       case AliDCSValue::kInt:
404         value = static_cast<Float_t>(val->GetInt());
405         break;
406       case AliDCSValue::kUInt:
407         value = static_cast<Float_t>(val->GetUInt());
408         break;
409       default:
410         continue;
411     }
412     if (fValCut>0 && TMath::Abs(value)>fValCut) continue;   // refuse values greater than cut
413     if (fDiffCut>0 ) {
414       if ( out>0 && skipped<10 && TMath::Abs(value-y[out-1])>fDiffCut) {
415         skipped++;                               // refuse values changing 
416         continue;                                // by > cut  in one time step
417       }                                          
418       skipped=0;
419     }                                         
420     if (val->GetTimeStamp()-time0>1000000) continue;
421     lastTime=TTimeStamp((time_t)val->GetTimeStamp(),0);
422     x[out] = (val->GetTimeStamp()-time0)/kSecInHour; // give times in fractions of hours 
423     y[out] = val->GetFloat();
424     out++;    
425   }
426   if (!keepStart) fStartTime=firstTime;
427   fEndTime=lastTime;
428   TGraph * graph = new TGraph(out,x,y);
429   delete [] x;
430   delete [] y;
431   return graph;
432 }
433   
434 //_____________________________________________________________________________
435 AliDCSSensor* AliDCSSensorArray::GetSensor(Int_t IdDCS) 
436 {
437  //
438  //  Return sensor information for sensor specified by IdDCS
439  //
440  Int_t nsensors = fSensors->GetEntries();
441  for (Int_t isensor=0; isensor<nsensors; isensor++) {
442    AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
443    if (entry->GetIdDCS() == IdDCS) return entry;
444  }
445  return 0;
446 }
447 //_____________________________________________________________________________
448 AliDCSSensor* AliDCSSensorArray::GetSensor(const TString& stringID)
449 {
450  //
451  //  Return sensor information for sensor specified by stringID
452  //
453  Int_t nsensors = fSensors->GetEntries();
454  for (Int_t isensor=0; isensor<nsensors; isensor++) {
455    AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
456    if (entry->GetStringID() == stringID) return entry;
457  }
458  return 0;
459 }
460 //_____________________________________________________________________________
461 AliDCSSensor* AliDCSSensorArray::GetSensor(Double_t x, Double_t y, Double_t z)
462 {
463  //
464  //  Return sensor closest to given position
465  //
466  Int_t nsensors = fSensors->GetEntries();
467  Double_t dist2min=1e99;
468  Double_t xs,ys,zs,dist2;
469  Int_t ind=-1;
470  for (Int_t isensor=0; isensor<nsensors; isensor++) {
471    AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
472    xs = entry->GetX();
473    ys = entry->GetY();
474    zs = entry->GetZ();
475    dist2 = (x-xs)*(x-xs) + (y-ys)*(y-ys) + (z-zs)*(z-zs);
476    if (dist2 < dist2min) {
477       ind=isensor;
478       dist2min = dist2;
479    }
480  }
481  if ( ind >= 0 ) {
482     return (AliDCSSensor*)fSensors->At(ind);
483  } else {
484     return 0;
485  }
486 }
487 //_____________________________________________________________________________
488 AliDCSSensor* AliDCSSensorArray::GetSensorNum(Int_t ind)
489 {
490  //
491  //  Return sensor given by array index
492  //
493  return (AliDCSSensor*)fSensors->At(ind);
494 }
495
496 //_____________________________________________________________________________
497 Int_t AliDCSSensorArray::SetSensor(const TString& stringID,
498                           const  AliDCSSensor& sensor)
499 {
500  //
501  //  Update sensor information for sensor specified by stringID
502  //
503  Int_t nsensors = fSensors->GetEntries();
504  for (Int_t isensor=0; isensor<nsensors; isensor++) {
505    AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
506    if (entry->GetStringID() == stringID) 
507      {
508       new ((*fSensors)[isensor])AliDCSSensor(sensor);
509       return isensor;
510      }
511  }
512  return -1;
513 }
514 //_____________________________________________________________________________
515 void AliDCSSensorArray::SetSensorNum(const Int_t ind, const AliDCSSensor& sensor)
516 {
517  //
518  //  Update sensor information for sensor at index ind
519  //
520    new ((*fSensors)[ind])AliDCSSensor(sensor);
521    return;
522 }
523 //_____________________________________________________________________________
524 void AliDCSSensorArray::RemoveSensorNum(Int_t ind)
525 {
526  //
527  //  Return sensor given by array index
528  //
529
530   delete fSensors->RemoveAt(ind);
531   fSensors->Compress();
532 }
533 //_____________________________________________________________________________
534 void AliDCSSensorArray::RemoveSensor(Int_t IdDCS)
535 {
536  //
537  //  Deletes Sensor by given IdDCS
538  //
539
540   Int_t nsensors = fSensors->GetEntries();
541   for (Int_t isensor=0; isensor<nsensors; isensor++) { // loop over sensors
542     AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
543     if (entry->GetIdDCS()==IdDCS) {
544       delete fSensors->RemoveAt(isensor);
545       break;
546     }
547   }
548   fSensors->Compress();
549 }
550 //_____________________________________________________________________________
551 TArrayI AliDCSSensorArray::OutsideThreshold(Double_t threshold, UInt_t timeSec, Bool_t below) const
552 {
553  //
554  // Return sensors with values outside threshold at time defined by second
555  // parameter
556  // By default sensors with values below threshold are listed, if third
557  // parameter is set to kFALSE sensors with values above threshold are listed
558  //
559   Int_t nsensors = fSensors->GetEntries();
560   TArrayI array(nsensors);
561   Int_t outside=0;
562   for (Int_t isensor=0; isensor<nsensors; isensor++) { // loop over sensors
563     AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
564     Double_t val=entry->GetValue(timeSec);
565     if (below) {
566       if (val<threshold) array[outside++] = entry->GetIdDCS();
567     } else {
568       if (val>threshold) array[outside++] = entry->GetIdDCS();
569     }    
570   }
571   array.Set(outside);
572   return array;
573 }
574  
575 //_____________________________________________________________________________
576 Int_t AliDCSSensorArray::GetFirstIdDCS() const
577 {
578  //
579  //  Return DCS Id of first sensor
580  //
581  if ( fSensors != 0 ) {
582     return ((AliDCSSensor*)fSensors->At(0))->GetIdDCS();
583  } else {
584     return 0;
585  }
586 }
587
588 //_____________________________________________________________________________
589 Int_t AliDCSSensorArray::GetLastIdDCS() const 
590 {
591  //
592  //  Return DCS Id of last sensor
593  //
594  if ( fSensors != 0 ) {
595     Int_t last = fSensors->GetEntries();
596     return ((AliDCSSensor*)fSensors->At(last-1))->GetIdDCS();
597  } else {
598     return 0;
599  }
600 }
601 //_____________________________________________________________________________
602 void AliDCSSensorArray::ClearGraph()
603 {
604   //
605   // Delete DCS graphs from all sensors in array
606   //
607    
608    Int_t nsensors = fSensors->GetEntries();
609    for ( Int_t isensor=0; isensor<nsensors; isensor++) {
610      AliDCSSensor *sensor = (AliDCSSensor*)fSensors->At(isensor);
611      TGraph *gr = sensor->GetGraph();
612      if ( gr != 0 ) {
613        delete gr;
614        gr = 0;
615      }
616      sensor->SetGraph(0);
617    }
618 }
619 //_____________________________________________________________________________
620 void AliDCSSensorArray::ClearFit()
621 {
622   //
623   // Delete spline fits from all sensors in array
624   //
625
626    Int_t nsensors = fSensors->GetEntries();
627    for ( Int_t isensor=0; isensor<nsensors; isensor++) {
628      AliDCSSensor *sensor = (AliDCSSensor*)fSensors->At(isensor);
629      AliSplineFit *fit = sensor->GetFit();
630      if ( fit != 0 ) {
631        delete fit;
632        fit = 0;
633      }
634      sensor->SetFit(0);
635    }
636 }
637 //_____________________________________________________________________________
638 void AliDCSSensorArray::AddSensors(AliDCSSensorArray *newSensors)
639 {
640   //
641   // add sensors from two sensor arrays
642   //
643   
644   Int_t numNew = newSensors->NumSensors();
645   Int_t numOld = fSensors->GetEntries();
646   fSensors->Expand(numOld+numNew);
647   for (Int_t i=0;i<numNew;i++) {
648     AliDCSSensor *sens = newSensors->GetSensorNum(i);
649     new ((*fSensors)[numOld+i]) AliDCSSensor(*sens);
650   }
651 }  
652