1 /**************************************************************************
2 * Copyright(c) 2006-07, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 ////////////////////////////////////////////////////////////////////////////////
19 // Class describing TPC temperature sensors (including pointers to graphs/fits//
20 // Authors: Marian Ivanov, Haavard Helstrup and Martin Siska //
22 ////////////////////////////////////////////////////////////////////////////////
25 #include "AliDCSSensor.h"
27 ClassImp(AliDCSSensor)
29 const Double_t kSecInHour = 3600.; // seconds in one hour
33 AliDCSSensor::AliDCSSensor():
46 // Standard constructor
50 AliDCSSensor::AliDCSSensor(const AliDCSSensor& source) :
53 fIdDCS(source.fIdDCS),
54 fStringID(source.fStringID),
55 fStartTime(source.fStartTime),
56 fEndTime(source.fEndTime),
57 fGraph(source.fGraph),
67 AliDCSSensor& AliDCSSensor::operator=(const AliDCSSensor& source){
69 // assignment operator
71 if (&source == this) return *this;
72 new (this) AliDCSSensor(source);
77 //_____________________________________________________________________________
78 Double_t AliDCSSensor::GetValue(UInt_t timeSec)
81 // Get temperature value for actual sensor
82 // timeSec given as offset from start-of-map measured in seconds
83 // *NOTE* In the current TPC setup, start-of-map is defined as the
84 // first measured point for each sensor. This will be different
85 // for each sensor in the array. If you want to get a value at the
86 // same absolute time, use AliDCSSensor::GetValue(TTimeStamp time)
87 // or AliDCSSensorArray::GetValue (UInt_t timeSec, Int_t sensor)
88 // which measure offsets with respect to the (global) start-of-run
91 return Eval(TTimeStamp((time_t)(fStartTime+timeSec),0),inside);
93 //_____________________________________________________________________________
94 Double_t AliDCSSensor::GetValue(TTimeStamp time)
96 // Get temperature value for actual sensor
97 // time given as absolute TTimeStamp
100 return Eval(time, inside);
103 //_____________________________________________________________________________
105 Double_t AliDCSSensor::Eval(const TTimeStamp& time, Bool_t inside) const
108 // Return temperature at given time
109 // If time < start of map return value at start of map, inside = false
110 // If time > end of map return value at end of map, inside = false
112 UInt_t timeSec = time.GetSec();
113 UInt_t diff = timeSec-fStartTime;
116 if ( timeSec < fStartTime ) {
120 if ( timeSec > fEndTime ) {
122 diff = fEndTime-fStartTime;
125 Double_t timeHour = diff/kSecInHour;
127 return fFit->Eval(timeHour);
130 return EvalGraph(timeHour);
136 //_____________________________________________________________________________
137 Double_t AliDCSSensor::EvalGraph(const Double_t& timeHour) const
140 // Extract last value in graph observed before time given by timeHour
143 // return -99 if point specified is before beginning of graph
144 Double_t x=0; Double_t y=0;
145 fGraph->GetPoint(0,x,y);
146 if ( timeHour < x ) return -99;
148 // return previous point when first time > timeHour is observed
150 Int_t npoints = fGraph->GetN();
151 for (Int_t i=1; i<npoints; i++) {
152 fGraph->GetPoint(i,x,y);
153 if ( timeHour < x ) {
154 fGraph->GetPoint(i-1,x,y);
159 // return last point if all times are < timeHour
164 //_____________________________________________________________________________
165 TGraph* AliDCSSensor::MakeGraph(Int_t nPoints, Bool_t debug) const
168 // Make graph from start time to end time of DCS values
173 UInt_t stepTime = (fEndTime-fStartTime)/nPoints;
176 printf ("Start time %d, End time %d, step time %d\n",
177 fStartTime,fEndTime,stepTime);
178 TTimeStamp t((time_t)fStartTime,0); t.Print();
179 TTimeStamp t2((time_t)fEndTime,0); t2.Print();
182 if ( !fFit ) return 0;
184 Double_t *x = new Double_t[nPoints+1];
185 Double_t *y = new Double_t[nPoints+1];
186 for (Int_t ip=0; ip<nPoints; ip++) {
187 x[ip] = (time_t)(fStartTime+ip*stepTime);
188 y[ip] = fFit->Eval(ip*stepTime/kSecInHour);
190 TTimeStamp t3((time_t)x[ip],0);
191 printf ("x=%f, y=%f ",x[ip],y[ip]);
196 TGraph *graph = new TGraph(nPoints,x,y);
200 graph->GetXaxis()->SetTimeDisplay(1);
201 graph->GetXaxis()->SetLabelOffset(0.02);
202 graph->GetXaxis()->SetTimeFormat("#splitline{%d/%m}{%H:%M}");
207 //_____________________________________________________________________________
209 TClonesArray * AliDCSSensor::ReadTree(TTree* tree) {
211 // read values from ascii file
214 Int_t nentries = tree->GetEntries();
223 tree->SetBranchAddress("StringID",&stringId);
224 tree->SetBranchAddress("IdDCS",&idDCS);
225 tree->SetBranchAddress("Num",&num);
226 tree->SetBranchAddress("X",&x);
227 tree->SetBranchAddress("Y",&y);
228 tree->SetBranchAddress("Z",&z);
230 // firstSensor = (Int_t)tree->GetMinimum("ECha");
231 // lastSensor = (Int_t)tree->GetMaximum("ECha");
233 TClonesArray * array = new TClonesArray("AliDCSSensor",nentries);
234 printf ("nentries = %d\n",nentries);
236 for (Int_t isensor=0; isensor<nentries; isensor++){
237 AliDCSSensor * sens = new ((*array)[isensor])AliDCSSensor;
238 tree->GetEntry(isensor);
239 sens->SetId(isensor);
240 sens->SetIdDCS(idDCS);
241 sens->SetStringID(TString(stringId));