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