]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCPreprocessor.cxx
New method to get the ratio between the expected and actual cluster shape. Will be...
[u/mrichter/AliRoot.git] / TPC / AliTPCPreprocessor.cxx
1 /**************************************************************************
2  * Copyright(c) 2007, 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 #include "AliTPCPreprocessor.h"
18 #include "AliShuttleInterface.h"
19
20 #include "AliCDBMetaData.h"
21 #include "AliDCSValue.h"
22 #include "AliLog.h"
23 #include "AliTPCSensorTempArray.h"
24 #include "AliTPCROC.h"
25 #include "AliTPCCalROC.h"
26 #include "AliTPCCalPad.h"
27 #include "AliTPCCalibPedestal.h"
28 #include "TFile.h"
29 #include "TTree.h"
30 #include "TEnv.h"
31
32 #include <TTimeStamp.h>
33
34 const Int_t kValCutTemp = 100;               // discard temperatures > 100 degrees
35 const Int_t kDiffCutTemp = 5;                // discard temperature differences > 5 degrees
36 const TString kPedestalRunType = "PEDESTAL_RUN";  // pedestal run identifier
37
38 //
39 // This class is the SHUTTLE preprocessor for the TPC detector.
40 // It contains several components, this far the part containing
41 // temperatures is implemented
42 //
43
44 ClassImp(AliTPCPreprocessor)
45
46 //______________________________________________________________________________________________
47 AliTPCPreprocessor::AliTPCPreprocessor(AliShuttleInterface* shuttle) :
48   AliPreprocessor("TPC",shuttle),
49   fConfEnv(0), fTemp(0), fPressure(0), fConfigOK(kTRUE), fROC(0)
50 {
51   // constructor
52   fROC = AliTPCROC::Instance();
53 }
54 //______________________________________________________________________________________________
55 // AliTPCPreprocessor::AliTPCPreprocessor(const AliTPCPreprocessor& org) :
56 //   AliPreprocessor(org),
57 //   fConfEnv(0), fTemp(0), fPressure(0), fConfigOK(kTRUE)
58 // {
59 //   // copy constructor not implemented
60 //   //   -- missing underlying copy constructor in AliPreprocessor
61 //
62 //   Fatal("AliTPCPreprocessor", "copy constructor not implemented");
63 //
64 // //  fTemp = new AliTPCSensorTempArray(*(org.fTemp));
65 // }
66
67 //______________________________________________________________________________________________
68 AliTPCPreprocessor::~AliTPCPreprocessor()
69 {
70   // destructor
71
72   delete fTemp;
73   delete fPressure;
74 }
75 //______________________________________________________________________________________________
76 AliTPCPreprocessor& AliTPCPreprocessor::operator = (const AliTPCPreprocessor& )
77 {
78   Fatal("operator =", "assignment operator not implemented");
79   return *this;
80 }
81
82
83 //______________________________________________________________________________________________
84 void AliTPCPreprocessor::Initialize(Int_t run, UInt_t startTime,
85         UInt_t endTime)
86 {
87   // Creates AliTestDataDCS object
88
89   AliPreprocessor::Initialize(run, startTime, endTime);
90
91         AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s", run,
92                 TTimeStamp(startTime).AsString(),
93                 TTimeStamp(endTime).AsString()));
94
95   // Preprocessor configuration
96
97         AliCDBEntry* entry = GetFromOCDB("Config", "Preprocessor");
98         if (entry) fConfEnv = (TEnv*) entry->GetObject();
99         if ( fConfEnv==0 ) {
100            Log("AliTPCPreprocsessor: Preprocessor Config OCDB entry missing.\n");
101            fConfigOK = kFALSE;
102            return;
103         }
104
105   // Temperature sensors
106
107         TTree *confTree = 0;
108         entry = GetFromOCDB("Config", "Temperature");
109         if (entry) confTree = (TTree*) entry->GetObject();
110         if ( confTree==0 ) {
111            Log("AliTPCPreprocsessor: Temperature Config OCDB entry missing.\n");
112            fConfigOK = kFALSE;
113            return;
114         }
115         fTemp = new AliTPCSensorTempArray(fStartTime, fEndTime, confTree);
116         fTemp->SetValCut(kValCutTemp);
117         fTemp->SetDiffCut(kDiffCutTemp);
118
119   // Pressure sensors
120
121         confTree=0;
122         entry=0;
123         entry = GetFromOCDB("Config", "Pressure");
124         if (entry) confTree = (TTree*) entry->GetObject();
125         if ( confTree==0 ) {
126            Log("AliTPCPreprocsessor: Pressure Config OCDB entry missing.\n");
127            fConfigOK = kFALSE;
128            return;
129         }
130         fPressure = new AliDCSSensorArray(fStartTime, fEndTime, confTree);
131
132 }
133
134 //______________________________________________________________________________________________
135 UInt_t AliTPCPreprocessor::Process(TMap* dcsAliasMap)
136 {
137   // Fills data into TPC calibrations objects
138
139   // Amanda servers provide information directly through dcsAliasMap
140
141   if (!dcsAliasMap) return 9;
142   if (dcsAliasMap->GetEntries() == 0 ) return 9;
143   if (!fConfigOK) return 9;
144
145   TString runType = GetRunType();
146
147   // Temperature sensors are processed by AliTPCCalTemp
148
149
150   UInt_t tempResult = MapTemperature(dcsAliasMap);
151   UInt_t result=tempResult;
152
153   // Pressure sensors
154
155   UInt_t pressureResult = MapPressure(dcsAliasMap);
156   result += pressureResult;
157
158   // Other calibration information will be retrieved through FXS files
159   //  examples:
160   //    TList* fileSourcesDAQ = GetFile(AliShuttleInterface::kDAQ, "pedestals");
161   //    const char* fileNamePed = GetFile(AliShuttleInterface::kDAQ, "pedestals", "LDC1");
162   //
163   //    TList* fileSourcesHLT = GetFile(AliShuttleInterface::kHLT, "calib");
164   //    const char* fileNameHLT = GetFile(AliShuttleInterface::kHLT, "calib", "LDC1");
165
166
167   if(runType == kPedestalRunType) {
168     Int_t pedestalSource = AliShuttleInterface::kDAQ;
169     TString source = fConfEnv->GetValue("Pedestal","DAQ");
170     source.ToUpper();
171     if ( source == "HLT" ) pedestalSource = AliShuttleInterface::kHLT;
172     UInt_t pedestalResult = ExtractPedestals(pedestalSource);
173     result += pedestalResult;
174
175   }
176
177
178   return result;
179 }
180 //______________________________________________________________________________________________
181 UInt_t AliTPCPreprocessor::MapTemperature(TMap* dcsAliasMap)
182 {
183
184    // extract DCS temperature maps. Perform fits to save space
185
186   UInt_t result=0;
187   TMap *map = fTemp->ExtractDCS(dcsAliasMap);
188   if (map) {
189     fTemp->MakeSplineFit(map);
190     AliInfo(Form("Temperature values extracted, fits performed.\n"));
191   } else {
192     Log("AliTPCPreprocsessor: no temperature map extracted. \n");
193     result=9;
194   }
195   delete map;
196   // Now store the final CDB file
197
198   if ( result == 0 ) {
199         AliCDBMetaData metaData;
200         metaData.SetBeamPeriod(0);
201         metaData.SetResponsible("Haavard Helstrup");
202         metaData.SetComment("Preprocessor AliTPC data base entries.");
203
204         Bool_t storeOK = Store("Calib", "Temperature", fTemp, &metaData, 0, kFALSE);
205         if ( !storeOK )  result=1;
206
207    }
208
209    return result;
210
211 }
212 //______________________________________________________________________________________________
213
214 UInt_t AliTPCPreprocessor::MapPressure(TMap* dcsAliasMap)
215 {
216
217    // extract DCS temperature maps. Perform fits to save space
218
219   UInt_t result=0;
220   TMap *map = fPressure->ExtractDCS(dcsAliasMap);
221   if (map) {
222     fPressure->MakeSplineFit(map);
223     AliInfo(Form("Pressure values extracted, fits performed.\n"));
224   } else {
225     Log("AliTPCPreprocsessor: no atmospheric pressure map extracted. \n");
226     result=9;
227   }
228   delete map;
229   // Now store the final CDB file
230
231   if ( result == 0 ) {
232         AliCDBMetaData metaData;
233         metaData.SetBeamPeriod(0);
234         metaData.SetResponsible("Haavard Helstrup");
235         metaData.SetComment("Preprocessor AliTPC data base entries.");
236
237         Bool_t storeOK = Store("Calib", "Pressure", fPressure, &metaData, 0, 0);
238         if ( !storeOK ) result=1;
239
240    }
241
242    return result;
243
244 }
245
246
247 //______________________________________________________________________________________________
248
249 UInt_t AliTPCPreprocessor::ExtractPedestals(Int_t sourceFXS)
250 {
251  //
252  //  Read pedestal file from file exchage server
253  //  Keep original entry from OCDB in case no new pedestals are available
254  //
255  AliTPCCalPad *calPadPed=0;
256  AliCDBEntry* entry = GetFromOCDB("Calib", "Pedestals");
257  if (entry) calPadPed = (AliTPCCalPad*)entry->GetObject();
258  if ( calPadPed==NULL ) {
259      Log("AliTPCPreprocsessor: No previous TPC pedestal entry available.\n");
260      calPadPed = new AliTPCCalPad("PedestalsMean","PedestalsMean");
261  }
262
263  AliTPCCalPad *calPadRMS=0;
264  entry = GetFromOCDB("Calib", "Noise");
265  if (entry) calPadRMS = (AliTPCCalPad*)entry->GetObject();
266  if ( calPadRMS==NULL ) {
267      Log("AliTPCPreprocsessor: No previous TPC noise entry available.\n");
268      calPadRMS = new AliTPCCalPad("PedestalsRMS","PedestalsRMS");
269  }
270
271
272  UInt_t result=0;
273
274  Int_t nSectors = fROC->GetNSectors();
275  TList* list = GetFileSources(sourceFXS,"pedestals");
276  if (list) {
277
278 //  loop through all files from LDCs
279
280     UInt_t index = 0;
281     while (list->At(index)!=NULL) {
282      TObjString* fileNameEntry = (TObjString*) list->At(index);
283      if (fileNameEntry!=NULL) {
284         TString fileName = GetFile(sourceFXS, "pedestals",
285                                          fileNameEntry->GetString().Data());
286         TFile *f = TFile::Open(fileName);
287         AliTPCCalibPedestal *calPed;
288         f->GetObject("AliTPCCalibPedestal",calPed);
289
290         //  replace entries for the sectors available in the present file
291
292         for (Int_t sector=0; sector<nSectors; sector++) {
293            AliTPCCalROC *rocPed=calPed->GetCalRocPedestal(sector, kFALSE);
294            if ( rocPed )  calPadPed->SetCalROC(rocPed,sector);
295            AliTPCCalROC *rocRMS=calPed->GetCalRocRMS(sector, kFALSE);
296            if ( rocRMS )  calPadRMS->SetCalROC(rocRMS,sector);
297         }
298       }
299      ++index;
300     }  // while(list)
301 //
302 //  Store updated pedestal entry to OCDB
303 //
304     AliCDBMetaData metaData;
305     metaData.SetBeamPeriod(0);
306     metaData.SetResponsible("Haavard Helstrup");
307     metaData.SetComment("Preprocessor AliTPC data base entries.");
308
309     Bool_t storeOK = Store("Calib", "Pedestals", calPadPed, &metaData, 0, kTRUE);
310     if ( !storeOK ) ++result;
311     storeOK = Store("Calib", "PadNoise", calPadRMS, &metaData, 0, kTRUE);
312     if ( !storeOK ) ++result;
313
314   }
315
316   return result;
317 }
318
319 //______________________________________________________________________________________________
320
321