AliRoot Core  3dc7879 (3dc7879)
AliDCSSensor.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 2006-07, 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 
18 // //
19 // Class describing time dependent values read from DCS sensors //
20 // (including pointers to graphs/fits) //
21 // Authors: Marian Ivanov, Haavard Helstrup and Martin Siska //
22 // //
24 
25 
26 #include "AliDCSSensor.h"
27 #include "TDatime.h"
28 #include "TCanvas.h"
29 ClassImp(AliDCSSensor)
30 
31 const Double_t kSecInHour = 3600.; // seconds in one hour
32 
33 
34 
36  fId(),
37  fIdDCS(0),
38  fStringID(),
39  fStartTime(0),
40  fEndTime(0),
41  fGraph(0),
42  fFit(0),
43  fX(0),
44  fY(0),
45  fZ(0)
46 {
47  //
48  // Standard constructor
49  //
50 }
51 
53  TNamed(source),
54  fId(source.fId),
55  fIdDCS(source.fIdDCS),
56  fStringID(source.fStringID),
57  fStartTime(source.fStartTime),
58  fEndTime(source.fEndTime),
59  fGraph(0),
60  fFit(0),
61  fX(source.fX),
62  fY(source.fY),
63  fZ(source.fZ)
64 //
65 // Copy constructor
66 //
67 {
68  if (source.fGraph) fGraph = (TGraph*)source.fGraph->Clone();
69  if (source.fFit) fFit = (AliSplineFit*)source.fFit->Clone();
70 }
71 
73  //
74  // Destructor
75  //
76  if(fGraph)
77  delete fGraph;
78  fGraph=0;
79 
80 }
81 
83 //
84 // assignment operator
85 //
86  if (&source == this) return *this;
87  new (this) AliDCSSensor(source);
88 
89  return *this;
90 }
91 
92 
93 void AliDCSSensor::Print(const Option_t* option) const{
94  //
95  // print function
96  //
97  TString opt = option; opt.ToLower();
98  printf("%s:%s\n",GetTitle(), GetName());
99  printf("%s\n",fStringID.Data());
100 
101 }
102 
103 void AliDCSSensor::Draw(Option_t* option) {
104  //
105  // draw function - to viusalize sensor
106  // Unfortuantelly - it make a memory leak as function Draw does not return the object pointer
107  //
108  TCanvas * canvas = new TCanvas((fStringID+option).Data(), (fStringID+option).Data());
109  if (fGraph){
110  // transform points to time in s
111  Int_t npoints = fGraph->GetN();
112  for (Int_t i=0; i<npoints; i++){
113  fGraph->GetX()[i]=fGraph->GetX()[i]*3600+fStartTime;
114  }
115  fGraph->Draw("alp");
116  return;
117  }
118  canvas->cd();
119  TGraph * graph = MakeGraph(100); // memory leak - we can not modify the content - const method
120  if (graph){
121  graph->Draw(option); //
122  }
123 }
124 
125 
126 
127 //_____________________________________________________________________________
128 Double_t AliDCSSensor::GetValue(UInt_t timeSec)
129 {
130  //
131  // Get DCS value for actual sensor
132  // timeSec given as offset from start-of-map measured in seconds
133  // *NOTE* In the current TPC setup, start-of-map is defined as the
134  // first measured point for each sensor. This will be different
135  // for each sensor in the array. If you want to get a value at the
136  // same absolute time, use AliDCSSensor::GetValue(TTimeStamp time)
137  // or AliDCSSensorArray::GetValue (UInt_t timeSec, Int_t sensor)
138  // which measure offsets with respect to the (global) start-of-run
139  //
140  Bool_t inside=kTRUE;
141  return Eval(TTimeStamp((time_t)(fStartTime+timeSec),0),inside);
142 }
143 //_____________________________________________________________________________
144 Double_t AliDCSSensor::GetValue(TTimeStamp time)
145 {
146  // Get DCS value for actual sensor
147  // time given as absolute TTimeStamp
148  //
149  Bool_t inside=kTRUE;
150  return Eval(time, inside);
151 }
152 
153 //_____________________________________________________________________________
154 
155 Double_t AliDCSSensor::Eval(const TTimeStamp& time, Bool_t& inside) const
156 {
157  //
158  // Return DCS value at given time
159  // The value is calculated from the AliSplineFit, if a fit is not available
160  // the most recent reading from the Graph of DCS points is returned (if
161  // the graph is present)
162  // If time < start of map return value at start of map, inside = false
163  // If time > end of map return value at end of map, inside = false
164 
165  UInt_t timeSec = time.GetSec();
166  UInt_t diff = timeSec-fStartTime;
167  inside = true;
168 
169  if ( timeSec < fStartTime ) {
170  inside=false;
171  diff=0;
172  }
173  if ( timeSec > fEndTime ) {
174  inside=false;
175  diff = fEndTime-fStartTime;
176  }
177 
178  Double_t timeHour = diff/kSecInHour;
179  if ( fFit ) {
180  return fFit->Eval(timeHour);
181  } else {
182  if ( fGraph ) {
183  return EvalGraph(timeHour);
184  } else {
185  return -99;
186  }
187  }
188 }
189 //_____________________________________________________________________________
190 
191 Double_t AliDCSSensor::EvalGraph(const TTimeStamp& time, Bool_t& inside) const
192 {
193  //
194  // Return DCS value from graph of DCS points (i.e return last reading before
195  // the time specified by TTimeStamp
196  // If time < start of map return value at start of map, inside = false
197  // If time > end of map return value at end of map, inside = false
198 
199  UInt_t timeSec = time.GetSec();
200  UInt_t diff = timeSec-fStartTime;
201  inside = true;
202 
203  if ( timeSec < fStartTime ) {
204  inside=false;
205  diff=0;
206  }
207  if ( timeSec > fEndTime ) {
208  inside=false;
209  diff = fEndTime-fStartTime;
210  }
211 
212  Double_t timeHour = diff/kSecInHour;
213  if ( fGraph ) {
214  return EvalGraph(timeHour);
215  } else {
216  return -99;
217  }
218 }
219 //_____________________________________________________________________________
220 Double_t AliDCSSensor::EvalGraph(const Double_t& timeHour) const
221 {
222  //
223  // Extract last value in graph observed before time given by timeHour
224  //
225 
226  // return -99 if point specified is before beginning of graph
227  Double_t x=0; Double_t y=0;
228  fGraph->GetPoint(0,x,y);
229  if ( timeHour < x ) return -99;
230 
231  // return previous point when first time > timeHour is observed
232 
233  Int_t npoints = fGraph->GetN();
234  for (Int_t i=1; i<npoints; i++) {
235  fGraph->GetPoint(i,x,y);
236  if ( timeHour < x ) {
237  fGraph->GetPoint(i-1,x,y);
238  return y;
239  }
240  }
241 
242  // return last point if all times are < timeHour
243  return y;
244 }
245 
246 
247 //_____________________________________________________________________________
248 TGraph* AliDCSSensor::MakeGraph(Int_t nPoints, Bool_t debug) const
249 {
250  //
251  // Make graph from start time to end time of DCS values
252  //
253 
254 
255 
256  UInt_t stepTime = (fEndTime-fStartTime)/nPoints;
257 
258  if (debug==kTRUE) {
259  printf ("Start time %d, End time %d, step time %d\n",
260  fStartTime,fEndTime,stepTime);
261  TTimeStamp t((time_t)fStartTime,0); t.Print();
262  TTimeStamp t2((time_t)fEndTime,0); t2.Print();
263  }
264 
265  if ( !fFit ) return 0;
266 
267  Double_t *x = new Double_t[nPoints+1];
268  Double_t *y = new Double_t[nPoints+1];
269  for (Int_t ip=0; ip<nPoints; ip++) {
270  x[ip] = (time_t)(fStartTime+ip*stepTime);
271  y[ip] = fFit->Eval(ip*stepTime/kSecInHour);
272  if (debug==kTRUE) {
273  TTimeStamp t3((time_t)x[ip],0);
274  printf ("x=%f, y=%f ",x[ip],y[ip]);
275  t3.Print();
276  }
277  }
278 
279  TGraph *graph = new TGraph(nPoints,x,y);
280  delete [] x;
281  delete [] y;
282 
283  graph->GetXaxis()->SetTimeDisplay(1);
284  graph->GetXaxis()->SetLabelOffset(0.02);
285  graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
286 
287  return graph;
288 }
289 
290 //_____________________________________________________________________________
291 
292 TClonesArray * AliDCSSensor::ReadTree(TTree* tree) {
293  //
294  // read values from ascii file
295  //
296 
297  Int_t nentries = tree->GetEntries();
298 
299  char stringId[100];
300  Int_t num=0;
301  Int_t idDCS=0;
302  Double_t x=0;
303  Double_t y=0;
304  Double_t z=0;
305 
306  tree->SetBranchAddress("StringID",&stringId);
307  tree->SetBranchAddress("IdDCS",&idDCS);
308  tree->SetBranchAddress("Num",&num);
309  tree->SetBranchAddress("X",&x);
310  tree->SetBranchAddress("Y",&y);
311  tree->SetBranchAddress("Z",&z);
312 
313  // firstSensor = (Int_t)tree->GetMinimum("ECha");
314  // lastSensor = (Int_t)tree->GetMaximum("ECha");
315 
316  TClonesArray * array = new TClonesArray("AliDCSSensor",nentries);
317  printf ("nentries = %d\n",nentries);
318 
319  for (Int_t isensor=0; isensor<nentries; isensor++){
320  AliDCSSensor * sens = new ((*array)[isensor])AliDCSSensor;
321  tree->GetEntry(isensor);
322  sens->SetId(isensor);
323  sens->SetIdDCS(idDCS);
324  sens->SetStringID(TString(stringId));
325  sens->SetX(x);
326  sens->SetY(y);
327  sens->SetZ(z);
328 
329  }
330  return array;
331 }
332 
printf("Chi2/npoints = %f\n", TMath::Sqrt(chi2/npoints))
virtual void Draw(Option_t *option="")
Double_t EvalGraph(const TTimeStamp &time, Bool_t &inside) const
void SetStringID(const TString &stringID)
Definition: AliDCSSensor.h:58
UInt_t fStartTime
Definition: AliDCSSensor.h:83
npoints
Definition: driftITSTPC.C:85
AliSplineFit * fFit
Definition: AliDCSSensor.h:86
TObjArray * array
Definition: AnalyzeLaser.C:12
void SetZ(Double_t z)
Definition: AliDCSSensor.h:62
TTree * tree
Double_t Eval(const TTimeStamp &time, Bool_t &inside) const
void SetIdDCS(Int_t iddcs)
Definition: AliDCSSensor.h:57
static TClonesArray * ReadTree(TTree *tree)
void SetX(Double_t x)
Definition: AliDCSSensor.h:60
virtual ~AliDCSSensor()
AliDCSSensor & operator=(const AliDCSSensor &source)
const Double_t kSecInHour
void SetY(Double_t y)
Definition: AliDCSSensor.h:61
Double_t Eval(Double_t x, Int_t deriv=0) const
Double_t GetValue(UInt_t timeSec)
virtual void Print(const Option_t *option="") const
TString fStringID
Definition: AliDCSSensor.h:82
void SetId(Int_t id)
Definition: AliDCSSensor.h:56
UInt_t fEndTime
Definition: AliDCSSensor.h:84
TGraph * MakeGraph(Int_t nPoints=100, Bool_t debug=kFALSE) const
Int_t debug
TGraph * fGraph
Definition: AliDCSSensor.h:85