]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - STEER/AliDCSSensor.cxx
fix codding violation in AliTRDseedV1 class
[u/mrichter/AliRoot.git] / STEER / AliDCSSensor.cxx
... / ...
CommitLineData
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
17////////////////////////////////////////////////////////////////////////////////
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// //
23////////////////////////////////////////////////////////////////////////////////
24
25
26#include "AliDCSSensor.h"
27#include "TDatime.h"
28ClassImp(AliDCSSensor)
29
30const Double_t kSecInHour = 3600.; // seconds in one hour
31
32
33
34AliDCSSensor::AliDCSSensor():
35 fId(),
36 fIdDCS(0),
37 fStringID(),
38 fStartTime(0),
39 fEndTime(0),
40 fGraph(0),
41 fFit(0),
42 fX(0),
43 fY(0),
44 fZ(0)
45{
46 //
47 // Standard constructor
48 //
49}
50
51AliDCSSensor::AliDCSSensor(const AliDCSSensor& source) :
52 TNamed(source),
53 fId(source.fId),
54 fIdDCS(source.fIdDCS),
55 fStringID(source.fStringID),
56 fStartTime(source.fStartTime),
57 fEndTime(source.fEndTime),
58 fGraph(0),
59 fFit(0),
60 fX(source.fX),
61 fY(source.fY),
62 fZ(source.fZ)
63//
64// Copy constructor
65//
66{
67 if (source.fGraph) fGraph = (TGraph*)source.fGraph->Clone();
68 if (source.fFit) fFit = (AliSplineFit*)source.fFit->Clone();
69}
70
71AliDCSSensor& AliDCSSensor::operator=(const AliDCSSensor& source){
72//
73// assignment operator
74//
75 if (&source == this) return *this;
76 new (this) AliDCSSensor(source);
77
78 return *this;
79}
80
81//_____________________________________________________________________________
82Double_t AliDCSSensor::GetValue(UInt_t timeSec)
83{
84 //
85 // Get DCS value for actual sensor
86 // timeSec given as offset from start-of-map measured in seconds
87 // *NOTE* In the current TPC setup, start-of-map is defined as the
88 // first measured point for each sensor. This will be different
89 // for each sensor in the array. If you want to get a value at the
90 // same absolute time, use AliDCSSensor::GetValue(TTimeStamp time)
91 // or AliDCSSensorArray::GetValue (UInt_t timeSec, Int_t sensor)
92 // which measure offsets with respect to the (global) start-of-run
93 //
94 Bool_t inside=kTRUE;
95 return Eval(TTimeStamp((time_t)(fStartTime+timeSec),0),inside);
96}
97//_____________________________________________________________________________
98Double_t AliDCSSensor::GetValue(TTimeStamp time)
99{
100 // Get DCS value for actual sensor
101 // time given as absolute TTimeStamp
102 //
103 Bool_t inside=kTRUE;
104 return Eval(time, inside);
105}
106
107//_____________________________________________________________________________
108
109Double_t AliDCSSensor::Eval(const TTimeStamp& time, Bool_t& inside) const
110{
111 //
112 // Return DCS value at given time
113 // The value is calculated from the AliSplineFit, if a fit is not available
114 // the most recent reading from the Graph of DCS points is returned (if
115 // the graph is present)
116 // If time < start of map return value at start of map, inside = false
117 // If time > end of map return value at end of map, inside = false
118
119 UInt_t timeSec = time.GetSec();
120 UInt_t diff = timeSec-fStartTime;
121 inside = true;
122
123 if ( timeSec < fStartTime ) {
124 inside=false;
125 diff=0;
126 }
127 if ( timeSec > fEndTime ) {
128 inside=false;
129 diff = fEndTime-fStartTime;
130 }
131
132 Double_t timeHour = diff/kSecInHour;
133 if ( fFit ) {
134 return fFit->Eval(timeHour);
135 } else {
136 if ( fGraph ) {
137 return EvalGraph(timeHour);
138 } else {
139 return -99;
140 }
141 }
142}
143//_____________________________________________________________________________
144
145Double_t AliDCSSensor::EvalGraph(const TTimeStamp& time, Bool_t& inside) const
146{
147 //
148 // Return DCS value from graph of DCS points (i.e return last reading before
149 // the time specified by TTimeStamp
150 // If time < start of map return value at start of map, inside = false
151 // If time > end of map return value at end of map, inside = false
152
153 UInt_t timeSec = time.GetSec();
154 UInt_t diff = timeSec-fStartTime;
155 inside = true;
156
157 if ( timeSec < fStartTime ) {
158 inside=false;
159 diff=0;
160 }
161 if ( timeSec > fEndTime ) {
162 inside=false;
163 diff = fEndTime-fStartTime;
164 }
165
166 Double_t timeHour = diff/kSecInHour;
167 if ( fGraph ) {
168 return EvalGraph(timeHour);
169 } else {
170 return -99;
171 }
172}
173//_____________________________________________________________________________
174Double_t AliDCSSensor::EvalGraph(const Double_t& timeHour) const
175{
176 //
177 // Extract last value in graph observed before time given by timeHour
178 //
179
180 // return -99 if point specified is before beginning of graph
181 Double_t x=0; Double_t y=0;
182 fGraph->GetPoint(0,x,y);
183 if ( timeHour < x ) return -99;
184
185 // return previous point when first time > timeHour is observed
186
187 Int_t npoints = fGraph->GetN();
188 for (Int_t i=1; i<npoints; i++) {
189 fGraph->GetPoint(i,x,y);
190 if ( timeHour < x ) {
191 fGraph->GetPoint(i-1,x,y);
192 return y;
193 }
194 }
195
196 // return last point if all times are < timeHour
197 return y;
198}
199
200
201//_____________________________________________________________________________
202TGraph* AliDCSSensor::MakeGraph(Int_t nPoints, Bool_t debug) const
203{
204 //
205 // Make graph from start time to end time of DCS values
206 //
207
208
209
210 UInt_t stepTime = (fEndTime-fStartTime)/nPoints;
211
212 if (debug==kTRUE) {
213 printf ("Start time %d, End time %d, step time %d\n",
214 fStartTime,fEndTime,stepTime);
215 TTimeStamp t((time_t)fStartTime,0); t.Print();
216 TTimeStamp t2((time_t)fEndTime,0); t2.Print();
217 }
218
219 if ( !fFit ) return 0;
220
221 Double_t *x = new Double_t[nPoints+1];
222 Double_t *y = new Double_t[nPoints+1];
223 for (Int_t ip=0; ip<nPoints; ip++) {
224 x[ip] = (time_t)(fStartTime+ip*stepTime);
225 y[ip] = fFit->Eval(ip*stepTime/kSecInHour);
226 if (debug==kTRUE) {
227 TTimeStamp t3((time_t)x[ip],0);
228 printf ("x=%f, y=%f ",x[ip],y[ip]);
229 t3.Print();
230 }
231 }
232
233 TGraph *graph = new TGraph(nPoints,x,y);
234 delete [] x;
235 delete [] y;
236
237 graph->GetXaxis()->SetTimeDisplay(1);
238 graph->GetXaxis()->SetLabelOffset(0.02);
239 graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
240
241 return graph;
242}
243
244//_____________________________________________________________________________
245
246TClonesArray * AliDCSSensor::ReadTree(TTree* tree) {
247 //
248 // read values from ascii file
249 //
250
251 Int_t nentries = tree->GetEntries();
252
253 char stringId[100];
254 Int_t num=0;
255 Int_t idDCS=0;
256 Double_t x=0;
257 Double_t y=0;
258 Double_t z=0;
259
260 tree->SetBranchAddress("StringID",&stringId);
261 tree->SetBranchAddress("IdDCS",&idDCS);
262 tree->SetBranchAddress("Num",&num);
263 tree->SetBranchAddress("X",&x);
264 tree->SetBranchAddress("Y",&y);
265 tree->SetBranchAddress("Z",&z);
266
267 // firstSensor = (Int_t)tree->GetMinimum("ECha");
268 // lastSensor = (Int_t)tree->GetMaximum("ECha");
269
270 TClonesArray * array = new TClonesArray("AliDCSSensor",nentries);
271 printf ("nentries = %d\n",nentries);
272
273 for (Int_t isensor=0; isensor<nentries; isensor++){
274 AliDCSSensor * sens = new ((*array)[isensor])AliDCSSensor;
275 tree->GetEntry(isensor);
276 sens->SetId(isensor);
277 sens->SetIdDCS(idDCS);
278 sens->SetStringID(TString(stringId));
279 sens->SetX(x);
280 sens->SetY(y);
281 sens->SetZ(z);
282
283 }
284 return array;
285}
286