]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCPreprocessor.cxx
New version of the TPC preprocessor. New preprocessor configuration directory (Haavard)
[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", "Temperature");
98         if (entry) fConfEnv = (TEnv*) entry->GetObject();
99         if ( fConfEnv==0 ) {
100            AliWarning(Form("Preprocessor Config OCDB entry missing.\n"));
101            Log("AliTPCPreprocsessor: Preprocessor Config OCDB entry missing.\n");
102         }
103
104   // Temperature sensors
105
106         TTree *confTree = 0;
107         entry = GetFromOCDB("Config", "Temperature");
108         if (entry) confTree = (TTree*) entry->GetObject();
109         if ( confTree==0 ) {
110            AliError(Form("Temperature Config OCDB entry missing.\n"));
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            AliError(Form("Pressure Config OCDB entry missing.\n"));
127            Log("AliTPCPreprocsessor: Pressure Config OCDB entry missing.\n");
128            fConfigOK = kFALSE;
129            return;
130         }
131         fPressure = new AliDCSSensorArray(fStartTime, fEndTime, confTree);
132
133 }
134
135 //______________________________________________________________________________________________
136 UInt_t AliTPCPreprocessor::Process(TMap* dcsAliasMap)
137 {
138   // Fills data into TPC calibrations objects
139
140   // Amanda servers provide information directly through dcsAliasMap
141
142   if (!dcsAliasMap) return 9;
143   if (dcsAliasMap->GetEntries() == 0 ) return 9;
144   if (!fConfigOK) return 9;
145
146   TString runType = GetRunType();
147
148   // Temperature sensors are processed by AliTPCCalTemp
149
150
151   UInt_t tempResult = MapTemperature(dcsAliasMap);
152   UInt_t result=tempResult;
153
154   // Pressure sensors
155
156   UInt_t pressureResult = MapPressure(dcsAliasMap);
157   result += pressureResult;
158
159   // Other calibration information will be retrieved through FXS files
160   //  examples:
161   //    TList* fileSourcesDAQ = GetFile(AliShuttleInterface::kDAQ, "pedestals");
162   //    const char* fileNamePed = GetFile(AliShuttleInterface::kDAQ, "pedestals", "LDC1");
163   //
164   //    TList* fileSourcesHLT = GetFile(AliShuttleInterface::kHLT, "calib");
165   //    const char* fileNameHLT = GetFile(AliShuttleInterface::kHLT, "calib", "LDC1");
166
167
168   if(runType == kPedestalRunType) {
169     Int_t pedestalSource = AliShuttleInterface::kDAQ;
170     TString source = fConfEnv->GetValue("Pedestal","DAQ");
171     source.ToUpper();
172     if ( source == "HLT" ) pedestalSource = AliShuttleInterface::kHLT;
173     UInt_t pedestalResult = ExtractPedestals(pedestalSource);
174     result += pedestalResult;
175
176   }
177
178
179   return result;
180 }
181 //______________________________________________________________________________________________
182 UInt_t AliTPCPreprocessor::MapTemperature(TMap* dcsAliasMap)
183 {
184
185    // extract DCS temperature maps. Perform fits to save space
186
187   UInt_t result=0;
188   TMap *map = fTemp->ExtractDCS(dcsAliasMap);
189   if (map) {
190     fTemp->MakeSplineFit(map);
191     AliInfo(Form("Temperature values extracted, fits performed.\n"));
192   } else {
193     AliError(Form("No temperature map extracted.\n"));
194     Log("AliTPCPreprocsessor: no temperature map extracted. \n");
195     result=9;
196   }
197   delete map;
198   // Now store the final CDB file
199
200   if ( result == 0 ) {
201         AliCDBMetaData metaData;
202         metaData.SetBeamPeriod(0);
203         metaData.SetResponsible("Haavard Helstrup");
204         metaData.SetComment("Preprocessor AliTPC data base entries.");
205
206         Bool_t storeOK = Store("Calib", "Temperature", fTemp, &metaData, 0, kFALSE);
207         if ( !storeOK )  result=1;
208
209    }
210
211    return result;
212
213 }
214 //______________________________________________________________________________________________
215
216 UInt_t AliTPCPreprocessor::MapPressure(TMap* dcsAliasMap)
217 {
218
219    // extract DCS temperature maps. Perform fits to save space
220
221   UInt_t result=0;
222   TMap *map = fPressure->ExtractDCS(dcsAliasMap);
223   if (map) {
224     fPressure->MakeSplineFit(map);
225     AliInfo(Form("Pressure values extracted, fits performed.\n"));
226   } else {
227     AliError(Form("No atmospheric pressure map extracted.\n"));
228     Log("AliTPCPreprocsessor: no atmospheric pressure map extracted. \n");
229     result=9;
230   }
231   delete map;
232   // Now store the final CDB file
233
234   if ( result == 0 ) {
235         AliCDBMetaData metaData;
236         metaData.SetBeamPeriod(0);
237         metaData.SetResponsible("Haavard Helstrup");
238         metaData.SetComment("Preprocessor AliTPC data base entries.");
239
240         Bool_t storeOK = Store("Calib", "Pressure", fPressure, &metaData, 0, 0);
241         if ( !storeOK ) result=1;
242
243    }
244
245    return result;
246
247 }
248
249
250 //______________________________________________________________________________________________
251
252 UInt_t AliTPCPreprocessor::ExtractPedestals(Int_t sourceFXS)
253 {
254  //
255  //  Read pedestal file from file exchage server
256  //  Keep original entry from OCDB in case no new pedestals are available
257  //
258  AliTPCCalPad *calPadPed=0;
259  AliCDBEntry* entry = GetFromOCDB("Calib", "Pedestals");
260  if (entry) calPadPed = (AliTPCCalPad*)entry->GetObject();
261  if ( calPadPed==NULL ) {
262      AliWarning(Form("No previous TPC pedestal entry available.\n"));
263      Log("AliTPCPreprocsessor: No previous TPC pedestal entry available.\n");
264      calPadPed = new AliTPCCalPad("PedestalsMean","PedestalsMean");
265  }
266
267  AliTPCCalPad *calPadRMS=0;
268  entry = GetFromOCDB("Calib", "Noise");
269  if (entry) calPadRMS = (AliTPCCalPad*)entry->GetObject();
270  if ( calPadRMS==NULL ) {
271      AliWarning(Form("No previous TPC noise entry available.\n"));
272      Log("AliTPCPreprocsessor: No previous TPC noise entry available.\n");
273      calPadRMS = new AliTPCCalPad("PedestalsRMS","PedestalsRMS");
274  }
275
276
277  UInt_t result=0;
278
279  Int_t nSectors = fROC->GetNSectors();
280  TList* list = GetFileSources(sourceFXS,"pedestals");
281  if (list) {
282
283 //  loop through all files from LDCs
284
285     UInt_t index = 0;
286     while (list->At(index)!=NULL) {
287      TObjString* fileNameEntry = (TObjString*) list->At(index);
288      if (fileNameEntry!=NULL) {
289         TString fileName = GetFile(sourceFXS, "pedestals",
290                                          fileNameEntry->GetString().Data());
291         TFile *f = TFile::Open(fileName);
292         AliTPCCalibPedestal *calPed;
293         f->GetObject("AliTPCCalibPedestal",calPed);
294
295         //  replace entries for the sectors available in the present file
296
297         for (Int_t sector=0; sector<nSectors; sector++) {
298            AliTPCCalROC *rocPed=calPed->GetCalRocPedestal(sector, kFALSE);
299            if ( rocPed )  calPadPed->SetCalROC(rocPed,sector);
300            AliTPCCalROC *rocRMS=calPed->GetCalRocRMS(sector, kFALSE);
301            if ( rocRMS )  calPadRMS->SetCalROC(rocRMS,sector);
302         }
303       }
304      ++index;
305     }  // while(list)
306 //
307 //  Store updated pedestal entry to OCDB
308 //
309     AliCDBMetaData metaData;
310     metaData.SetBeamPeriod(0);
311     metaData.SetResponsible("Haavard Helstrup");
312     metaData.SetComment("Preprocessor AliTPC data base entries.");
313
314     Bool_t storeOK = Store("Calib", "Pedestals", calPadPed, &metaData, 0, kTRUE);
315     if ( !storeOK ) ++result;
316     storeOK = Store("Calib", "PadNoise", calPadRMS, &metaData, 0, kTRUE);
317     if ( !storeOK ) ++result;
318
319   }
320
321   return result;
322 }
323
324 //______________________________________________________________________________________________
325
326