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