]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCPreprocessor.cxx
Updates to include new run types CALIBRATION_PULSER and STANDALONE_PULSER
[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 "AliTPCCalibPulser.h"
29 #include "AliTPCCalibCE.h"
30 #include "TFile.h"
31 #include "TTree.h"
32 #include "TGraph.h" 
33 #include "TEnv.h"
34 #include "TParameter.h"
35
36 #include <TTimeStamp.h>
37
38 const Int_t kValCutTemp = 100;               // discard temperatures > 100 degrees
39 const Int_t kDiffCutTemp = 5;                // discard temperature differences > 5 degrees
40 const TString kPedestalRunType = "PEDESTAL";  // pedestal run identifier
41 const TString kPulserRunType = "CALIBRATION_PULSER";   // pulser run identifier
42 const TString kPhysicsRunType = "PHYSICS";   // physics run identifier
43 const TString kStandAloneRunType = "STANDALONE"; // standalone run identifier
44 const TString kStandAlonePulserRunType = "STANDALONE_PULSER"; // standalone run identifier
45 const TString kCosmicRunType = "COSMIC"; // cosmic run identifier
46 const TString kLaserRunType = "LASER";   // laser run identifier
47 const TString kDaqRunType = "DAQ"; // DAQ run identifier
48 const TString kAmandaTemp = "TPC_PT_%d_TEMPERATURE"; // Amanda string for temperature entries
49 //const Double_t kFitFraction = 0.7;                 // Fraction of DCS sensor fits required              
50 const Double_t kFitFraction = -1.0;          // Don't require minimum number of fits in commissioning run 
51
52 //
53 // This class is the SHUTTLE preprocessor for the TPC detector.
54 //
55
56 ClassImp(AliTPCPreprocessor)
57
58 //______________________________________________________________________________________________
59 AliTPCPreprocessor::AliTPCPreprocessor(AliShuttleInterface* shuttle) :
60   AliPreprocessor("TPC",shuttle),
61   fConfEnv(0), fTemp(0), fHighVoltage(0), fHighVoltageStat(0), fConfigOK(kTRUE), fROC(0)
62 {
63   // constructor
64   fROC = AliTPCROC::Instance();
65
66   // define run types to be processed
67   
68   AddRunType(kPedestalRunType);
69   AddRunType(kPulserRunType);
70   AddRunType(kPhysicsRunType);
71   AddRunType(kStandAloneRunType);
72   AddRunType(kStandAlonePulserRunType);
73   AddRunType(kCosmicRunType);
74   AddRunType(kLaserRunType);
75   AddRunType(kDaqRunType);
76   
77 }
78 //______________________________________________________________________________________________
79  AliTPCPreprocessor::AliTPCPreprocessor(const AliTPCPreprocessor&  ) :
80    AliPreprocessor("TPC",0),
81    fConfEnv(0), fTemp(0), fHighVoltage(0), fHighVoltageStat(0), fConfigOK(kTRUE), fROC(0)
82  {
83
84    Fatal("AliTPCPreprocessor", "copy constructor not implemented");
85 //
86 // //  fTemp = new AliTPCSensorTempArray(*(org.fTemp));
87  }
88
89 //______________________________________________________________________________________________
90 AliTPCPreprocessor::~AliTPCPreprocessor()
91 {
92   // destructor
93
94   delete fTemp;
95   delete fHighVoltage;
96 }
97 //______________________________________________________________________________________________
98 AliTPCPreprocessor& AliTPCPreprocessor::operator = (const AliTPCPreprocessor& )
99 {
100   Fatal("operator =", "assignment operator not implemented");
101   return *this;
102 }
103
104
105 //______________________________________________________________________________________________
106 void AliTPCPreprocessor::Initialize(Int_t run, UInt_t startTime,
107         UInt_t endTime)
108 {
109   // Creates AliTestDataDCS object -- start maps half an hour beforre actual run start
110
111   UInt_t startTimeLocal = startTime-1800;
112
113   AliPreprocessor::Initialize(run, startTimeLocal, endTime);
114
115         AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s", run,
116                 TTimeStamp((time_t)startTime,0).AsString(),
117                 TTimeStamp((time_t)endTime,0).AsString()));
118
119   // Preprocessor configuration
120
121         AliCDBEntry* entry = GetFromOCDB("Config", "Preprocessor");
122         if (entry) fConfEnv = (TEnv*) entry->GetObject();
123         if ( fConfEnv==0 ) {
124            Log("AliTPCPreprocsessor: Preprocessor Config OCDB entry missing.\n");
125            fConfigOK = kFALSE;
126            return;
127         }
128
129   // Temperature sensors
130
131        TTree *confTree = 0;
132
133        TString tempConf = fConfEnv->GetValue("Temperature","ON");
134        tempConf.ToUpper();
135        if (tempConf != "OFF" ) {
136         entry = GetFromOCDB("Config", "Temperature");
137         if (entry) confTree = (TTree*) entry->GetObject();
138         if ( confTree==0 ) {
139            Log("AliTPCPreprocsessor: Temperature Config OCDB entry missing.\n");
140            fConfigOK = kFALSE;
141            return;
142         }
143         fTemp = new AliTPCSensorTempArray(startTimeLocal, fEndTime, confTree, kAmandaTemp);
144         fTemp->SetValCut(kValCutTemp);
145         fTemp->SetDiffCut(kDiffCutTemp);
146        }
147
148   // High voltage measurements
149
150       TString hvConf = fConfEnv->GetValue("HighVoltage","ON");
151       hvConf.ToUpper();
152       if (hvConf != "OFF" ) { 
153         confTree=0;
154         entry=0;
155         entry = GetFromOCDB("Config", "HighVoltage");
156         if (entry) confTree = (TTree*) entry->GetObject();
157         if ( confTree==0 ) {
158            Log("AliTPCPreprocsessor: High Voltage Config OCDB entry missing.\n");
159            fConfigOK = kFALSE;
160            return;
161         }
162         fHighVoltage = new AliDCSSensorArray(startTimeLocal, fEndTime, confTree);
163       }
164
165    // High voltage status values
166      
167       TString hvStatConf = fConfEnv->GetValue("HighVoltageStat","ON");
168       hvStatConf.ToUpper();
169       if (hvStatConf != "OFF" ) { 
170         confTree=0;
171         entry=0;
172         entry = GetFromOCDB("Config", "HighVoltageStat");
173         if (entry) confTree = (TTree*) entry->GetObject();
174         if ( confTree==0 ) {
175            Log("AliTPCPreprocsessor: High Voltage Status Config OCDB entry missing.\n");
176            fConfigOK = kFALSE;
177            return;
178         }
179         fHighVoltageStat = new AliDCSSensorArray(startTimeLocal, fEndTime, confTree);
180       }
181 }
182
183 //______________________________________________________________________________________________
184 UInt_t AliTPCPreprocessor::Process(TMap* dcsAliasMap)
185 {
186   // Fills data into TPC calibrations objects
187
188   // Amanda servers provide information directly through dcsAliasMap
189
190   
191   if (!fConfigOK) return 9;
192   UInt_t result = 0;
193   TObjArray *resultArray = new TObjArray();
194   TString errorHandling = fConfEnv->GetValue("ErrorHandling","ON");
195   errorHandling.ToUpper();
196   TObject * status;
197
198   UInt_t dcsResult=0;
199   if (errorHandling == "OFF" ) {
200     if (!dcsAliasMap) dcsResult=1;
201     if (dcsAliasMap->GetEntries() == 0 ) dcsResult=1;  
202     status = new TParameter<int>("dcsResult",dcsResult);
203     resultArray->Add(status);
204   } else {
205     if (!dcsAliasMap) return 9;
206     if (dcsAliasMap->GetEntries() == 0 ) return 9;
207   }
208
209
210   
211
212   TString runType = GetRunType();
213
214   // Temperature sensors are processed by AliTPCCalTemp
215
216   TString tempConf = fConfEnv->GetValue("Temperature","ON");
217   tempConf.ToUpper();
218   if (tempConf != "OFF" ) {
219     UInt_t tempResult = MapTemperature(dcsAliasMap);
220     result=tempResult;
221     status = new TParameter<int>("tempResult",tempResult);
222     resultArray->Add(status);
223   }
224
225   // High Voltage recordings
226
227
228   TString hvConf = fConfEnv->GetValue("HighVoltage","ON");
229   hvConf.ToUpper();
230   if (hvConf != "OFF" ) { 
231    UInt_t hvResult = MapHighVoltage(dcsAliasMap);
232    result+=hvResult;
233    status = new TParameter<int>("hvResult",hvResult);
234    resultArray->Add(status);
235  }
236
237   // Other calibration information will be retrieved through FXS files
238   //  examples:
239   //    TList* fileSourcesDAQ = GetFile(AliShuttleInterface::kDAQ, "pedestals");
240   //    const char* fileNamePed = GetFile(AliShuttleInterface::kDAQ, "pedestals", "LDC1");
241   //
242   //    TList* fileSourcesHLT = GetFile(AliShuttleInterface::kHLT, "calib");
243   //    const char* fileNameHLT = GetFile(AliShuttleInterface::kHLT, "calib", "LDC1");
244
245   // pedestal entries
246
247   if(runType == kPedestalRunType) {
248     Int_t numSources = 1;
249     Int_t pedestalSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
250     TString source = fConfEnv->GetValue("Pedestal","DAQ");
251     source.ToUpper();
252     if (source != "OFF" ) { 
253      if ( source == "HLT") pedestalSource[0] = AliShuttleInterface::kHLT;
254      if (!GetHLTStatus()) pedestalSource[0] = AliShuttleInterface::kDAQ;
255      if (source == "HLTDAQ" ) {
256          numSources=2;
257          pedestalSource[0] = AliShuttleInterface::kHLT;
258          pedestalSource[1] = AliShuttleInterface::kDAQ;
259      }
260      if (source == "DAQHLT" ) numSources=2;
261      UInt_t pedestalResult=0;
262      for (Int_t i=0; i<numSources; i++ ) {      
263        pedestalResult = ExtractPedestals(pedestalSource[i]);
264        if ( pedestalResult == 0 ) break;
265      }
266      result += pedestalResult;
267      status = new TParameter<int>("pedestalResult",pedestalResult);
268      resultArray->Add(status);
269     }
270   }
271
272   // pulser trigger processing
273
274   if(runType == kPulserRunType) {
275     Int_t numSources = 1;
276     Int_t pulserSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
277     TString source = fConfEnv->GetValue("Pulser","DAQ");
278     source.ToUpper();
279     if ( source != "OFF") { 
280      if ( source == "HLT") pulserSource[0] = AliShuttleInterface::kHLT;
281      if (!GetHLTStatus()) pulserSource[0] = AliShuttleInterface::kDAQ;
282      if (source == "HLTDAQ" ) {
283          numSources=2;
284          pulserSource[0] = AliShuttleInterface::kHLT;
285          pulserSource[1] = AliShuttleInterface::kDAQ;
286      }
287      if (source == "DAQHLT" ) numSources=2;
288      UInt_t pulserResult=0;
289      for (Int_t i=0; i<numSources; i++ ) {      
290        pulserResult = ExtractPulser(pulserSource[i]);
291        if ( pulserResult == 0 ) break;
292      }
293      result += pulserResult;
294      status = new TParameter<int>("pulserResult",pulserResult);
295      resultArray->Add(status);
296     }
297   }
298
299
300
301   // Central Electrode processing
302
303   if( runType == kPhysicsRunType || runType == kStandAloneRunType || 
304       runType == kDaqRunType ) {    
305
306 //   if (true) {                 // do CE processing for all run types
307     Int_t numSources = 1;
308     Int_t ceSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
309     TString source = fConfEnv->GetValue("CE","DAQ");
310     source.ToUpper();
311     if ( source != "OFF" ) { 
312      if ( source == "HLT") ceSource[0] = AliShuttleInterface::kHLT;
313      if (!GetHLTStatus()) ceSource[0] = AliShuttleInterface::kDAQ;
314      if (source == "HLTDAQ" ) {
315         numSources=2;
316         ceSource[0] = AliShuttleInterface::kHLT;
317         ceSource[1] = AliShuttleInterface::kDAQ;
318      }
319      if (source == "DAQHLT" ) numSources=2;
320      UInt_t ceResult=0;
321      for (Int_t i=0; i<numSources; i++ ) {      
322        ceResult = ExtractCE(ceSource[i]);
323        if ( ceResult == 0 ) break;
324      }
325      result += ceResult;
326      status = new TParameter<int>("ceResult",ceResult);
327      resultArray->Add(status);
328     }
329   }
330   
331   if (errorHandling == "OFF" ) {
332     AliCDBMetaData metaData;
333     metaData.SetBeamPeriod(0);
334     metaData.SetResponsible("Haavard Helstrup");
335     metaData.SetComment("Preprocessor AliTPC status.");
336     Store("Calib", "PreprocStatus", resultArray, &metaData, 0, kFALSE);
337     resultArray->Delete();
338     return 0;
339   } else { 
340     return result;
341   }
342 }
343 //______________________________________________________________________________________________
344 UInt_t AliTPCPreprocessor::MapTemperature(TMap* dcsAliasMap)
345 {
346
347    // extract DCS temperature maps. Perform fits to save space
348
349   UInt_t result=0;
350   TMap *map = fTemp->ExtractDCS(dcsAliasMap);
351   if (map) {
352     fTemp->MakeSplineFit(map);
353     Double_t fitFraction = 1.0*fTemp->NumFits()/fTemp->NumSensors(); 
354     if (fitFraction > kFitFraction ) {
355       AliInfo(Form("Temperature values extracted, fits performed.\n"));
356     } else { 
357       Log ("Too few temperature maps fitted. \n");
358       result = 9;
359     }
360   } else {
361     Log("No temperature map extracted. \n");
362     result=9;
363   }
364   delete map;
365   // Now store the final CDB file
366
367   if ( result == 0 ) {
368         AliCDBMetaData metaData;
369         metaData.SetBeamPeriod(0);
370         metaData.SetResponsible("Haavard Helstrup");
371         metaData.SetComment("Preprocessor AliTPC data base entries.");
372
373         Bool_t storeOK = Store("Calib", "Temperature", fTemp, &metaData, 0, kFALSE);
374         if ( !storeOK )  result=1;
375
376    }
377
378    return result;
379
380 }
381
382 //______________________________________________________________________________________________
383 UInt_t AliTPCPreprocessor::MapHighVoltage(TMap* dcsAliasMap)
384 {
385
386    // extract DCS HV maps. Perform fits to save space
387
388   UInt_t result=0;
389   TMap *map = fHighVoltage->ExtractDCS(dcsAliasMap);
390   if (map) {
391     fHighVoltage->MakeSplineFit(map);
392     Double_t fitFraction = 1.0*fHighVoltage->NumFits()/fHighVoltage->NumSensors(); 
393     if (fitFraction > kFitFraction ) {
394       AliInfo(Form("High voltage recordings extracted, fits performed.\n"));
395     } else { 
396       Log ("Too few high voltage recordings fitted. \n");
397       result = 9;
398     }
399   } else {
400     Log("No high voltage recordings extracted. \n");
401     result=9;
402   }
403   delete map;
404
405   TString hvStatConf = fConfEnv->GetValue("HighVoltageStat","ON");
406   hvStatConf.ToUpper();
407   if (hvStatConf != "OFF" ) { 
408     TMap *map2 = fHighVoltageStat->ExtractDCS(dcsAliasMap);
409     if (map2) {
410       fHighVoltageStat->StoreGraph(map2);
411     } else {
412        Log("No high voltage status recordings extracted. \n");
413       result=9;
414     }
415     delete map2;
416
417     // add status maps to high voltage sensor array
418
419     fHighVoltage->AddSensors(fHighVoltageStat);
420    }
421   // Now store the final CDB file
422
423   if ( result == 0 ) {
424         AliCDBMetaData metaData;
425         metaData.SetBeamPeriod(0);
426         metaData.SetResponsible("Haavard Helstrup");
427         metaData.SetComment("Preprocessor AliTPC data base entries.");
428
429         Bool_t storeOK = Store("Calib", "HighVoltage", fHighVoltage, &metaData, 0, kFALSE);
430         if ( !storeOK )  result=1;
431
432    }
433
434    return result;
435
436 }
437
438
439 //______________________________________________________________________________________________
440
441 UInt_t AliTPCPreprocessor::ExtractPedestals(Int_t sourceFXS)
442 {
443  //
444  //  Read pedestal file from file exchage server
445  //  Keep original entry from OCDB in case no new pedestals are available
446  //
447  AliTPCCalPad *calPadPed=0;
448  AliCDBEntry* entry = GetFromOCDB("Calib", "Pedestals");
449  if (entry) calPadPed = (AliTPCCalPad*)entry->GetObject();
450  if ( calPadPed==NULL ) {
451      Log("AliTPCPreprocsessor: No previous TPC pedestal entry available.\n");
452      calPadPed = new AliTPCCalPad("PedestalsMean","PedestalsMean");
453  }
454
455  AliTPCCalPad *calPadRMS=0;
456  entry = GetFromOCDB("Calib", "PadNoise");
457  if (entry) calPadRMS = (AliTPCCalPad*)entry->GetObject();
458  if ( calPadRMS==NULL ) {
459      Log("AliTPCPreprocsessor: No previous TPC noise entry available.\n");
460      calPadRMS = new AliTPCCalPad("PedestalsRMS","PedestalsRMS");
461  }
462
463
464  UInt_t result=0;
465
466  Int_t nSectors = fROC->GetNSectors();
467  TList* list = GetFileSources(sourceFXS,"pedestals");
468  
469  if (list && list->GetEntries()>0) {
470
471 //  loop through all files from LDCs
472
473     Bool_t changed=false;
474     UInt_t index = 0;
475     while (list->At(index)!=NULL) {
476      TObjString* fileNameEntry = (TObjString*) list->At(index);
477      if (fileNameEntry!=NULL) {
478         TString fileName = GetFile(sourceFXS, "pedestals",
479                                          fileNameEntry->GetString().Data());
480         TFile *f = TFile::Open(fileName);
481         if (!f) {
482           Log ("Error opening pedestal file.");
483           result =2;
484           break;
485         }
486         AliTPCCalibPedestal *calPed;
487         f->GetObject("tpcCalibPedestal",calPed);
488         if ( !calPed ) {
489           Log ("No pedestal calibration object in file.");
490           result = 2;
491           break;
492         }
493
494         //  replace entries for the sectors available in the present file
495
496         changed=true;
497         for (Int_t sector=0; sector<nSectors; sector++) {
498            AliTPCCalROC *rocPed=calPed->GetCalRocPedestal(sector, kFALSE);
499            if ( rocPed )  calPadPed->SetCalROC(rocPed,sector);
500            AliTPCCalROC *rocRMS=calPed->GetCalRocRMS(sector, kFALSE);
501            if ( rocRMS )  calPadRMS->SetCalROC(rocRMS,sector);
502         }
503         delete calPed; 
504         f->Close();
505       }
506      ++index;
507     }  // while(list)
508 //
509 //  Store updated pedestal entry to OCDB
510 //
511     if (changed) {
512      AliCDBMetaData metaData;
513      metaData.SetBeamPeriod(0);
514      metaData.SetResponsible("Haavard Helstrup");
515      metaData.SetComment("Preprocessor AliTPC data base entries."); 
516  
517      Bool_t storeOK = Store("Calib", "Pedestals", calPadPed, &metaData, 0, kTRUE);
518      if ( !storeOK ) ++result;
519      storeOK = Store("Calib", "PadNoise", calPadRMS, &metaData, 0, kTRUE);
520      if ( !storeOK ) ++result;
521     }
522   } else {
523     Log ("Error: no entries in input file list!");
524     result = 1;
525   }
526
527   return result;
528 }
529
530 //______________________________________________________________________________________________
531
532
533 UInt_t AliTPCPreprocessor::ExtractPulser(Int_t sourceFXS)
534 {
535  //
536  //  Read pulser calibration file from file exchage server
537  //  Keep original entry from OCDB in case no new pulser calibration is available
538  //
539  TObjArray    *pulserObjects=0;
540  AliTPCCalPad *pulserTmean=0;
541  AliTPCCalPad *pulserTrms=0;
542  AliTPCCalPad *pulserQmean=0;
543  AliCDBEntry* entry = GetFromOCDB("Calib", "Pulser");
544  if (entry) pulserObjects = (TObjArray*)entry->GetObject();
545  if ( pulserObjects==NULL ) {
546      Log("AliTPCPreprocsessor: No previous TPC pulser entry available.\n");
547      pulserObjects = new TObjArray;    
548  }
549
550  pulserTmean = (AliTPCCalPad*)pulserObjects->FindObject("PulserTmean");
551  if ( !pulserTmean ) {
552     pulserTmean = new AliTPCCalPad("PulserTmean","PulserTmean");
553     pulserObjects->Add(pulserTmean);
554  }
555  pulserTrms = (AliTPCCalPad*)pulserObjects->FindObject("PulserTrms");
556  if ( !pulserTrms )  { 
557     pulserTrms = new AliTPCCalPad("PulserTrms","PulserTrms");
558     pulserObjects->Add(pulserTrms);
559  }
560  pulserQmean = (AliTPCCalPad*)pulserObjects->FindObject("PulserQmean");
561  if ( !pulserQmean )  { 
562     pulserQmean = new AliTPCCalPad("PulserQmean","PulserQmean");
563     pulserObjects->Add(pulserQmean);
564  }
565
566
567  UInt_t result=0;
568
569  Int_t nSectors = fROC->GetNSectors();
570  TList* list = GetFileSources(sourceFXS,"pulser");
571  
572  if (list && list->GetEntries()>0) {
573
574 //  loop through all files from LDCs
575
576     Bool_t changed=false;
577     UInt_t index = 0;
578     while (list->At(index)!=NULL) {
579      TObjString* fileNameEntry = (TObjString*) list->At(index);
580      if (fileNameEntry!=NULL) {
581         TString fileName = GetFile(sourceFXS, "pulser",
582                                          fileNameEntry->GetString().Data());
583         TFile *f = TFile::Open(fileName);
584         if (!f) {
585           Log ("Error opening pulser file.");
586           result =2;
587           break;
588         }
589         AliTPCCalibPulser *calPulser;
590         f->GetObject("tpcCalibPulser",calPulser);
591         if ( !calPulser ) {
592           Log ("No pulser calibration object in file.");
593           result = 2;
594           break;
595         }
596
597         //  replace entries for the sectors available in the present file
598
599         changed=true;
600         for (Int_t sector=0; sector<nSectors; sector++) {
601            AliTPCCalROC *rocTmean=calPulser->GetCalRocT0(sector);
602            if ( rocTmean )  pulserTmean->SetCalROC(rocTmean,sector);
603            AliTPCCalROC *rocTrms=calPulser->GetCalRocRMS(sector);
604            if ( rocTrms )  pulserTrms->SetCalROC(rocTrms,sector);
605            AliTPCCalROC *rocQmean=calPulser->GetCalRocQ(sector);
606            if ( rocQmean )  pulserQmean->SetCalROC(rocQmean,sector);
607         }
608        delete calPulser;
609        f->Close();
610       }
611      ++index;
612     }  // while(list)
613 //
614 //  Store updated pedestal entry to OCDB
615 //
616     if (changed) {
617      AliCDBMetaData metaData;
618      metaData.SetBeamPeriod(0);
619      metaData.SetResponsible("Haavard Helstrup");
620      metaData.SetComment("Preprocessor AliTPC data base entries.");
621
622      Bool_t storeOK = Store("Calib", "Pulser", pulserObjects, &metaData, 0, kTRUE);
623      if ( !storeOK ) ++result;
624     }  
625   } else {
626     Log ("Error: no entries in input file list!");
627     result = 1;
628   }
629
630   return result;
631 }
632
633 UInt_t AliTPCPreprocessor::ExtractCE(Int_t sourceFXS)
634 {
635  //
636  //  Read Central Electrode file from file exchage server
637  //  Keep original entry from OCDB in case no new CE calibration is available
638  //
639  TObjArray    *ceObjects=0;
640  AliTPCCalPad *ceTmean=0;
641  AliTPCCalPad *ceTrms=0;
642  AliTPCCalPad *ceQmean=0;
643  TObjArray    *rocTtime=0;  
644  TObjArray    *rocQtime=0;  
645
646  AliCDBEntry* entry = GetFromOCDB("Calib", "CE");
647  if (entry) ceObjects = (TObjArray*)entry->GetObject();
648  if ( ceObjects==NULL ) {
649      Log("AliTPCPreprocsessor: No previous TPC central electrode entry available.\n");
650      ceObjects = new TObjArray;    
651  }
652
653  Int_t nSectors = fROC->GetNSectors();
654
655  ceTmean = (AliTPCCalPad*)ceObjects->FindObject("CETmean");
656  if ( !ceTmean ) {
657     ceTmean = new AliTPCCalPad("CETmean","CETmean");
658     ceObjects->Add(ceTmean);
659  }
660  ceTrms = (AliTPCCalPad*)ceObjects->FindObject("CETrms");
661  if ( !ceTrms )  { 
662     ceTrms = new AliTPCCalPad("CETrms","CETrms");
663     ceObjects->Add(ceTrms);
664  }
665  ceQmean = (AliTPCCalPad*)ceObjects->FindObject("CEQmean");
666  if ( !ceQmean )  { 
667     ceQmean = new AliTPCCalPad("CEQmean","CEQmean");
668     ceObjects->Add(ceQmean);
669  }
670  //!new from here please have a look!!!
671  rocTtime = (TObjArray*)ceObjects->FindObject("rocTtime");
672  if ( !rocTtime ) {
673      rocTtime = new TObjArray(nSectors);
674      rocTtime->SetName("rocTtime");
675      ceObjects->Add(rocTtime);
676  }
677  
678  rocQtime = (TObjArray*)ceObjects->FindObject("rocQtime");
679  if ( !rocQtime ) {
680      rocQtime = new TObjArray(nSectors);
681      rocQtime->SetName("rocQtime");
682      ceObjects->Add(rocQtime);
683  }
684
685
686  UInt_t result=0;
687
688  TList* list = GetFileSources(sourceFXS,"CE");
689  
690  if (list && list->GetEntries()>0) {
691
692 //  loop through all files from LDCs
693
694     UInt_t index = 0;
695     while (list->At(index)!=NULL) {
696      TObjString* fileNameEntry = (TObjString*) list->At(index);
697      if (fileNameEntry!=NULL) {
698         TString fileName = GetFile(sourceFXS, "CE",
699                                          fileNameEntry->GetString().Data());
700         TFile *f = TFile::Open(fileName);
701         if (!f) {
702           Log ("Error opening central electrode file.");
703           result =2;
704           break;
705         }
706         AliTPCCalibCE *calCE;
707         f->GetObject("tpcCalibCE",calCE);
708
709         //  replace entries for the sectors available in the present file
710
711         for (Int_t sector=0; sector<nSectors; sector++) {
712            AliTPCCalROC *rocTmean=calCE->GetCalRocT0(sector);
713            if ( rocTmean )  ceTmean->SetCalROC(rocTmean,sector);
714            AliTPCCalROC *rocTrms=calCE->GetCalRocRMS(sector);
715            if ( rocTrms )  ceTrms->SetCalROC(rocTrms,sector);
716            AliTPCCalROC *rocQmean=calCE->GetCalRocQ(sector);
717            if ( rocQmean )  ceQmean->SetCalROC(rocQmean,sector);
718            TGraph *grT=calCE->MakeGraphTimeCE(sector,0,2); // T time graph
719            if ( grT ) rocTtime->AddAt(grT,sector);         
720            TGraph *grQ=calCE->MakeGraphTimeCE(sector,0,3); // Q time graph
721            if ( grQ ) rocTtime->AddAt(grQ,sector);         
722         }
723        delete calCE;
724        f->Close();
725       }
726      ++index;
727     }  // while(list)
728 //
729 //  Store updated pedestal entry to OCDB
730 //
731     AliCDBMetaData metaData;
732     metaData.SetBeamPeriod(0);
733     metaData.SetResponsible("Haavard Helstrup");
734     metaData.SetComment("Preprocessor AliTPC data base entries.");
735
736     Bool_t storeOK = Store("Calib", "CE", ceObjects, &metaData, 0, kTRUE);
737     if ( !storeOK ) ++result;
738     
739   } else {
740     Log ("Error: no entries!");
741     result = 1;
742   }
743
744   return result;
745 }