]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliDCSSensorArray.cxx
Read Config entry from OCDB (Haavard)
[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
26 ClassImp(AliDCSSensorArray)
27
28 const Double_t kSecInHour = 3600.; // seconds in one hour
29
30 //_____________________________________________________________________________
31 AliDCSSensorArray::AliDCSSensorArray():TNamed(), 
32   fMinGraph(10),
33   fMinPoints(10),
34   fIter(10),
35   fMaxDelta(0.0),
36   fFitReq(2),
37   fValCut(-1),
38   fDiffCut(-1),
39   fStartTime (2000,1,1,0,0,0),
40   fEndTime   (2000,1,1,0,0,0),
41   fSensors(0)
42 {
43   //
44   // AliDCSSensorArray default constructor
45   //
46
47 }
48 //_____________________________________________________________________________
49 AliDCSSensorArray::AliDCSSensorArray(TClonesArray *arr):TNamed(), 
50   fMinGraph(10),
51   fMinPoints(10),
52   fIter(10),
53   fMaxDelta(0.0),
54   fFitReq(2),
55   fValCut(-1),
56   fDiffCut(-1),
57   fStartTime (2000,1,1,0,0,0),
58   fEndTime   (2000,1,1,0,0,0),
59   fSensors(arr)
60 {
61   //
62   // AliDCSSensorArray special constructor taking TClonesArray from ReadList
63   //
64
65 }
66 //_____________________________________________________________________________
67 AliDCSSensorArray::AliDCSSensorArray(Int_t run, const char* dbEntry) : 
68   TNamed(),
69   fMinGraph(10),
70   fMinPoints(10),
71   fIter(10),
72   fMaxDelta(0.0),
73   fFitReq(2),
74   fValCut(-1),
75   fDiffCut(-1),
76   fStartTime (2000,1,1,0,0,0),
77   fEndTime   (2000,1,1,0,0,0),
78   fSensors(0)
79 {
80   //
81   // Read configuration from OCDB
82   //
83   
84   AliCDBEntry *entry = AliCDBManager::Instance()->Get(dbEntry,run); 
85   TTree *tree = (TTree*) entry->GetObject();
86   fSensors = AliDCSSensor::ReadTree(tree);
87 }
88 //_____________________________________________________________________________
89 AliDCSSensorArray::AliDCSSensorArray(UInt_t startTime, UInt_t endTime,
90                        TTree* confTree) :
91   TNamed(),
92   fMinGraph(10),
93   fMinPoints(10),
94   fIter(10),
95   fMaxDelta(0.0),
96   fFitReq(2),
97   fValCut(-1),
98   fDiffCut(-1),
99   fStartTime (2000,1,1,0,0,0),
100   fEndTime   (2000,1,1,0,0,0),
101   fSensors(0)
102
103 {
104   //
105   // AliDCSSensorArray constructor for Shuttle preprocessor 
106   //  (confTree read from OCDB)
107   //
108   fSensors = AliDCSSensor::ReadTree(confTree);
109   fSensors->BypassStreamer(kFALSE);
110   fStartTime = TTimeStamp(startTime);
111   fEndTime   = TTimeStamp(endTime);
112 }
113
114
115
116 //_____________________________________________________________________________
117 AliDCSSensorArray::AliDCSSensorArray(const AliDCSSensorArray &c):TNamed(c),
118   fMinGraph(c.fMinGraph),
119   fMinPoints(c.fMinPoints),
120   fIter(c.fIter),
121   fMaxDelta(c.fMaxDelta),
122   fFitReq(c.fFitReq),
123   fValCut(c.fValCut),
124   fDiffCut(c.fDiffCut),
125   fStartTime (c.fStartTime),
126   fEndTime   (c.fEndTime),
127   fSensors(0)
128
129 {
130   //
131   // AliDCSSensorArray copy constructor
132   //
133
134   ((AliDCSSensorArray &) c).Copy(*this);
135
136 }
137
138 ///_____________________________________________________________________________
139 AliDCSSensorArray::~AliDCSSensorArray()
140 {
141   //
142   // AliDCSSensorArray destructor
143   //
144   fSensors->Delete();
145   delete fSensors;
146
147 }
148
149 //_____________________________________________________________________________
150 AliDCSSensorArray &AliDCSSensorArray::operator=(const AliDCSSensorArray &c)
151 {
152   //
153   // Assignment operator
154   //
155
156   if (this != &c) ((AliDCSSensorArray &) c).Copy(*this);
157   return *this;
158
159 }
160
161 //_____________________________________________________________________________
162 void AliDCSSensorArray::Copy(TObject &c) const
163 {
164   //
165   // Copy function
166   //
167
168   TObject::Copy(c);
169 }
170 //_____________________________________________________________________________
171 void AliDCSSensorArray::SetGraph(TMap *map) 
172 {
173   // 
174   // Read graphs from DCS maps 
175   //
176   Int_t nsensors = fSensors->GetEntries();
177   for ( Int_t isensor=0; isensor<nsensors; isensor++) {
178     AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
179     TString stringID = entry->GetStringID();    
180     TGraph *gr = (TGraph*)map->GetValue(stringID.Data());
181     if ( gr !=0 ) {
182        entry->SetGraph((TGraph*)gr->Clone());
183     } else { 
184        entry->SetGraph(0);
185     }
186   } 
187 }  
188 //_____________________________________________________________________________
189 void AliDCSSensorArray::MakeSplineFit(TMap *map, Bool_t keepMap) 
190 {
191   // 
192   // Make spline fits 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 || gr->GetN() < fMinGraph) { 
200       entry->SetFit(0);
201       continue;    
202     }
203     AliSplineFit *fit = new AliSplineFit();
204     fit->InitKnots(gr,fMinPoints,fIter,fMaxDelta);
205     fit->SplineFit(fFitReq);
206     entry->SetStartTime(fStartTime);
207     entry->SetEndTime(fEndTime);
208     fit->Cleanup();
209     entry->SetFit(fit);
210     if (keepMap) entry->SetGraph(gr);
211   } 
212 }  
213
214
215 //_____________________________________________________________________________
216 Double_t AliDCSSensorArray::GetValue(UInt_t timeSec, Int_t sensor) 
217 {
218   // 
219   // Return sensor value at time timeSec (obtained from fitted function)
220   //  timeSec = time in seconds from start of run
221   //
222   AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(sensor);
223   return entry->GetValue(timeSec);
224 }
225     
226
227 //_____________________________________________________________________________
228 TMap* AliDCSSensorArray::ExtractDCS(TMap *dcsMap) 
229 {
230  //
231  // Extract temperature graphs from DCS maps
232  //
233  TMap *values = new TMap;
234  TObjArray * valueSet;
235  Int_t nsensors = fSensors->GetEntries();
236  for ( Int_t isensor=0; isensor<nsensors; isensor++) {
237    AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
238    TString stringID = entry->GetStringID();    
239    TPair *pair = (TPair*)dcsMap->FindObject(stringID.Data());
240    if ( pair ) {                            // only try to read values 
241                                             // if DCS object available
242      valueSet = (TObjArray*)pair->Value();
243      TGraph *graph = MakeGraph(valueSet);
244      values->Add(new TObjString(stringID.Data()),graph);
245    }
246  }
247  return values;
248 }
249
250 //_____________________________________________________________________________
251 TGraph* AliDCSSensorArray::MakeGraph(TObjArray* valueSet){
252   //
253   // Make graph of temperature values read from DCS map
254   //   (spline fit parameters will subsequently be obtained from this graph) 
255   //
256   Int_t nentries = valueSet->GetEntriesFast(); 
257   if ( nentries == 0 ) return 0;
258   
259   Float_t *x = new Float_t[nentries];
260   Float_t *y = new Float_t[nentries];
261   Int_t time0=0;
262   Int_t out=0;
263   Int_t skipped=0;
264   AliDCSValue *val = (AliDCSValue *)valueSet->At(0);
265   AliDCSValue::Type type = val->GetType();
266   if ( type == AliDCSValue::kInvalid || type == AliDCSValue::kBool ) return 0;
267   Float_t value;
268   for (Int_t i=0; i<nentries; i++){
269     val = (AliDCSValue *)valueSet->At(i);
270     if (!val) continue;
271     if (time0==0){
272       time0=val->GetTimeStamp();
273     }
274     switch ( type )
275     { 
276       case AliDCSValue::kFloat:
277         value = val->GetFloat();
278         break;
279       case AliDCSValue::kChar:
280         value = static_cast<Float_t>(val->GetChar());
281         break;
282       case AliDCSValue::kInt:
283         value = static_cast<Float_t>(val->GetInt());
284         break;
285       case AliDCSValue::kUInt:
286         value = static_cast<Float_t>(val->GetUInt());
287         break;
288       default:
289         continue;
290     }
291     if (fValCut>0 && TMath::Abs(value)>fValCut) continue;   // refuse values greater than cut
292     if (fDiffCut>0 ) {
293       if ( out>0 && skipped<10 && TMath::Abs(value-y[out-1])>fDiffCut) {
294         skipped++;                               // refuse values changing 
295         continue;                                // by > cut  in one time step
296       }                                          
297       skipped=0;
298     }                                         
299     if (val->GetTimeStamp()-time0>1000000) continue;
300     x[out] = (val->GetTimeStamp()-time0)/kSecInHour; // give times in fractions of hours 
301     y[out] = val->GetFloat();
302     out++;
303     
304   }
305   TGraph * graph = new TGraph(out,x,y);
306   delete [] x;
307   delete [] y;
308   return graph;
309 }
310   
311 //_____________________________________________________________________________
312 AliDCSSensor* AliDCSSensorArray::GetSensor(Int_t IdDCS) 
313 {
314  //
315  //  Return sensor information for sensor specified by IdDCS
316  //
317  Int_t nsensors = fSensors->GetEntries();
318  for (Int_t isensor=0; isensor<nsensors; isensor++) {
319    AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
320    if (entry->GetIdDCS() == IdDCS) return entry;
321  }
322  return 0;
323 }
324 //_____________________________________________________________________________
325 AliDCSSensor* AliDCSSensorArray::GetSensor(const TString& stringID) 
326 {
327  //
328  //  Return sensor information for sensor specified by IdDCS
329  //
330  Int_t nsensors = fSensors->GetEntries();
331  for (Int_t isensor=0; isensor<nsensors; isensor++) {
332    AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
333    if (entry->GetStringID() == stringID) return entry;
334  }
335  return 0;
336 }
337 //_____________________________________________________________________________
338 AliDCSSensor* AliDCSSensorArray::GetSensor(Double_t x, Double_t y, Double_t z) 
339 {
340  //
341  //  Return sensor closest to given position
342  //
343  Int_t nsensors = fSensors->GetEntries();
344  Double_t dist2min=1e99;
345  Double_t xs,ys,zs,dist2;
346  Int_t ind=-1;
347  for (Int_t isensor=0; isensor<nsensors; isensor++) {
348    AliDCSSensor *entry = (AliDCSSensor*)fSensors->At(isensor);
349    xs = entry->GetX();
350    ys = entry->GetY();
351    zs = entry->GetZ();
352    dist2 = (x-xs)*(x-xs) + (y-ys)*(y-ys) + (z-zs)*(z-zs);
353    if (dist2 < dist2min) {
354       ind=isensor;
355       dist2min = dist2;
356    } 
357  }
358  if ( ind >= 0 ) {
359     return (AliDCSSensor*)fSensors->At(ind);
360  } else {
361     return 0;
362  }
363 }
364
365 AliDCSSensor* AliDCSSensorArray::GetSensorNum(Int_t ind) 
366 {
367  //
368  //  Return sensor given by array index
369  //
370  return (AliDCSSensor*)fSensors->At(ind);
371 }
372
373 Int_t AliDCSSensorArray::GetFirstIdDCS() const
374 {
375  //
376  //  Return DCS Id of first sensor
377  //
378  if ( fSensors != 0 ) {
379     return ((AliDCSSensor*)fSensors->At(0))->GetIdDCS();
380  } else {
381     return 0;
382  }
383 }
384
385 Int_t AliDCSSensorArray::GetLastIdDCS() const 
386 {
387  //
388  //  Return DCS Id of last sensor
389  //
390  if ( fSensors != 0 ) {
391     Int_t last = fSensors->GetEntries();
392     return ((AliDCSSensor*)fSensors->At(last-1))->GetIdDCS();
393  } else {
394     return 0;
395  }
396 }
397 void AliDCSSensorArray::ClearGraph()
398 {
399   //
400   // Delete DCS graphs from all sensors in array
401   //
402    
403    Int_t nsensors = fSensors->GetEntries();
404    for ( Int_t isensor=0; isensor<nsensors; isensor++) {
405      AliDCSSensor *sensor = (AliDCSSensor*)fSensors->At(isensor);
406      TGraph *gr = sensor->GetGraph();
407      if ( gr != 0 ) {
408        delete gr;
409        gr = 0;
410      }
411      sensor->SetGraph(0);
412    }    
413 }  
414 void AliDCSSensorArray::ClearFit()
415 {
416   //
417   // Delete spline fits from all sensors in array
418   //
419    
420    Int_t nsensors = fSensors->GetEntries();
421    for ( Int_t isensor=0; isensor<nsensors; isensor++) {
422      AliDCSSensor *sensor = (AliDCSSensor*)fSensors->At(isensor);
423      AliSplineFit *fit = sensor->GetFit();
424      if ( fit != 0 ) {
425        delete fit;
426        fit = 0;
427      }
428      sensor->SetFit(0);
429    }    
430 }