]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCPreprocessor.cxx
AddRunType calls added in the constructor (F.Prino)
[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), 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), 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
164 //______________________________________________________________________________________________
165 UInt_t AliTPCPreprocessor::Process(TMap* dcsAliasMap)
166 {
167   // Fills data into TPC calibrations objects
168
169   // Amanda servers provide information directly through dcsAliasMap
170
171   
172   if (!fConfigOK) return 9;
173   UInt_t result = 0;
174   TObjArray *resultArray = new TObjArray();
175   TString errorHandling = fConfEnv->GetValue("ErrorHandling","ON");
176   errorHandling.ToUpper();
177   TObject * status;
178
179   UInt_t dcsResult=0;
180   if (errorHandling == "OFF" ) {
181     if (!dcsAliasMap) dcsResult=1;
182     if (dcsAliasMap->GetEntries() == 0 ) dcsResult=1;  
183     status = new TParameter<int>("dcsResult",dcsResult);
184     resultArray->Add(status);
185   } else {
186     if (!dcsAliasMap) return 9;
187     if (dcsAliasMap->GetEntries() == 0 ) return 9;
188   }
189
190
191   
192
193   TString runType = GetRunType();
194
195   // Temperature sensors are processed by AliTPCCalTemp
196
197   TString tempConf = fConfEnv->GetValue("Temperature","ON");
198   tempConf.ToUpper();
199   if (tempConf != "OFF" ) {
200     UInt_t tempResult = MapTemperature(dcsAliasMap);
201     result=tempResult;
202     status = new TParameter<int>("tempResult",tempResult);
203     resultArray->Add(status);
204   }
205
206   // High Voltage recordings
207
208
209   TString hvConf = fConfEnv->GetValue("HighVoltage","ON");
210   hvConf.ToUpper();
211   if (hvConf != "OFF" ) { 
212    UInt_t hvResult = MapHighVoltage(dcsAliasMap);
213    result+=hvResult;
214    status = new TParameter<int>("hvResult",hvResult);
215    resultArray->Add(status);
216  }
217
218   // Other calibration information will be retrieved through FXS files
219   //  examples:
220   //    TList* fileSourcesDAQ = GetFile(AliShuttleInterface::kDAQ, "pedestals");
221   //    const char* fileNamePed = GetFile(AliShuttleInterface::kDAQ, "pedestals", "LDC1");
222   //
223   //    TList* fileSourcesHLT = GetFile(AliShuttleInterface::kHLT, "calib");
224   //    const char* fileNameHLT = GetFile(AliShuttleInterface::kHLT, "calib", "LDC1");
225
226   // pedestal entries
227
228   if(runType == kPedestalRunType) {
229     Int_t numSources = 1;
230     Int_t pedestalSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
231     TString source = fConfEnv->GetValue("Pedestal","DAQ");
232     source.ToUpper();
233     if (source != "OFF" ) { 
234      if ( source == "HLT") pedestalSource[0] = AliShuttleInterface::kHLT;
235      if (!GetHLTStatus()) pedestalSource[0] = AliShuttleInterface::kDAQ;
236      if (source == "HLTDAQ" ) {
237          numSources=2;
238          pedestalSource[0] = AliShuttleInterface::kHLT;
239          pedestalSource[1] = AliShuttleInterface::kDAQ;
240      }
241      if (source == "DAQHLT" ) numSources=2;
242      UInt_t pedestalResult=0;
243      for (Int_t i=0; i<numSources; i++ ) {      
244        UInt_t pedestalResult = ExtractPedestals(pedestalSource[i]);
245        if ( pedestalResult == 0 ) break;
246      }
247      result += pedestalResult;
248      status = new TParameter<int>("pedestalResult",pedestalResult);
249      resultArray->Add(status);
250     }
251   }
252
253   // pulser trigger processing
254
255   if(runType == kPulserRunType) {
256     Int_t numSources = 1;
257     Int_t pulserSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
258     TString source = fConfEnv->GetValue("Pulser","DAQ");
259     source.ToUpper();
260     if ( source != "OFF") { 
261      if ( source == "HLT") pulserSource[0] = AliShuttleInterface::kHLT;
262      if (!GetHLTStatus()) pulserSource[0] = AliShuttleInterface::kDAQ;
263      if (source == "HLTDAQ" ) {
264          numSources=2;
265          pulserSource[0] = AliShuttleInterface::kHLT;
266          pulserSource[1] = AliShuttleInterface::kDAQ;
267      }
268      if (source == "DAQHLT" ) numSources=2;
269      UInt_t pulserResult=0;
270      for (Int_t i=0; i<numSources; i++ ) {      
271        pulserResult = ExtractPulser(pulserSource[i]);
272        if ( pulserResult == 0 ) break;
273      }
274      result += pulserResult;
275      status = new TParameter<int>("pulserResult",pulserResult);
276      resultArray->Add(status);
277     }
278   }
279
280
281
282   // Central Electrode processing
283
284 //  if( runType == kPhysicsRunType || runType == kStandAloneRunType || 
285 //      runType == kDaqRunType ) {    
286
287    if (true) {                 // do CE processing for all run types
288     Int_t numSources = 1;
289     Int_t ceSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
290     TString source = fConfEnv->GetValue("CE","DAQ");
291     source.ToUpper();
292     if ( source != "OFF" ) { 
293      if ( source == "HLT") ceSource[0] = AliShuttleInterface::kHLT;
294      if (!GetHLTStatus()) ceSource[0] = AliShuttleInterface::kDAQ;
295      if (source == "HLTDAQ" ) {
296         numSources=2;
297         ceSource[0] = AliShuttleInterface::kHLT;
298         ceSource[1] = AliShuttleInterface::kDAQ;
299      }
300      if (source == "DAQHLT" ) numSources=2;
301      UInt_t ceResult=0;
302      for (Int_t i=0; i<numSources; i++ ) {      
303        ceResult = ExtractCE(ceSource[i]);
304        if ( ceResult == 0 ) break;
305      }
306      result += ceResult;
307      status = new TParameter<int>("ceResult",ceResult);
308      resultArray->Add(status);
309     }
310   }
311   
312   if (errorHandling == "OFF" ) {
313     AliCDBMetaData metaData;
314     metaData.SetBeamPeriod(0);
315     metaData.SetResponsible("Haavard Helstrup");
316     metaData.SetComment("Preprocessor AliTPC status.");
317     Store("Calib", "PreprocStatus", resultArray, &metaData, 0, kFALSE);
318     resultArray->Delete();
319     return 0;
320   } else { 
321     return result;
322   }
323 }
324 //______________________________________________________________________________________________
325 UInt_t AliTPCPreprocessor::MapTemperature(TMap* dcsAliasMap)
326 {
327
328    // extract DCS temperature maps. Perform fits to save space
329
330   UInt_t result=0;
331   TMap *map = fTemp->ExtractDCS(dcsAliasMap);
332   if (map) {
333     fTemp->MakeSplineFit(map);
334     Double_t fitFraction = 1.0*fTemp->NumFits()/fTemp->NumSensors(); 
335     if (fitFraction > kFitFraction ) {
336       AliInfo(Form("Temperature values extracted, fits performed.\n"));
337     } else { 
338       Log ("Too few temperature maps fitted. \n");
339       result = 9;
340     }
341   } else {
342     Log("No temperature map extracted. \n");
343     result=9;
344   }
345   delete map;
346   // Now store the final CDB file
347
348   if ( result == 0 ) {
349         AliCDBMetaData metaData;
350         metaData.SetBeamPeriod(0);
351         metaData.SetResponsible("Haavard Helstrup");
352         metaData.SetComment("Preprocessor AliTPC data base entries.");
353
354         Bool_t storeOK = Store("Calib", "Temperature", fTemp, &metaData, 0, kFALSE);
355         if ( !storeOK )  result=1;
356
357    }
358
359    return result;
360
361 }
362
363 //______________________________________________________________________________________________
364 UInt_t AliTPCPreprocessor::MapHighVoltage(TMap* dcsAliasMap)
365 {
366
367    // extract DCS HV maps. Perform fits to save space
368
369   UInt_t result=0;
370   TMap *map = fHighVoltage->ExtractDCS(dcsAliasMap);
371   if (map) {
372     fHighVoltage->MakeSplineFit(map);
373     Double_t fitFraction = 1.0*fHighVoltage->NumFits()/fHighVoltage->NumSensors(); 
374     if (fitFraction > kFitFraction ) {
375       AliInfo(Form("High voltage recordings extracted, fits performed.\n"));
376     } else { 
377       Log ("Too few high voltage recordings fitted. \n");
378       result = 9;
379     }
380   } else {
381     Log("No high voltage recordings extracted. \n");
382     result=9;
383   }
384   delete map;
385   // Now store the final CDB file
386
387   if ( result == 0 ) {
388         AliCDBMetaData metaData;
389         metaData.SetBeamPeriod(0);
390         metaData.SetResponsible("Haavard Helstrup");
391         metaData.SetComment("Preprocessor AliTPC data base entries.");
392
393         Bool_t storeOK = Store("Calib", "HighVoltage", fHighVoltage, &metaData, 0, kFALSE);
394         if ( !storeOK )  result=1;
395
396    }
397
398    return result;
399
400 }
401
402
403 //______________________________________________________________________________________________
404
405 UInt_t AliTPCPreprocessor::ExtractPedestals(Int_t sourceFXS)
406 {
407  //
408  //  Read pedestal file from file exchage server
409  //  Keep original entry from OCDB in case no new pedestals are available
410  //
411  AliTPCCalPad *calPadPed=0;
412  AliCDBEntry* entry = GetFromOCDB("Calib", "Pedestals");
413  if (entry) calPadPed = (AliTPCCalPad*)entry->GetObject();
414  if ( calPadPed==NULL ) {
415      Log("AliTPCPreprocsessor: No previous TPC pedestal entry available.\n");
416      calPadPed = new AliTPCCalPad("PedestalsMean","PedestalsMean");
417  }
418
419  AliTPCCalPad *calPadRMS=0;
420  entry = GetFromOCDB("Calib", "PadNoise");
421  if (entry) calPadRMS = (AliTPCCalPad*)entry->GetObject();
422  if ( calPadRMS==NULL ) {
423      Log("AliTPCPreprocsessor: No previous TPC noise entry available.\n");
424      calPadRMS = new AliTPCCalPad("PedestalsRMS","PedestalsRMS");
425  }
426
427
428  UInt_t result=0;
429
430  Int_t nSectors = fROC->GetNSectors();
431  TList* list = GetFileSources(sourceFXS,"pedestals");
432  
433  if (list && list->GetEntries()>0) {
434
435 //  loop through all files from LDCs
436
437     Bool_t changed=false;
438     UInt_t index = 0;
439     while (list->At(index)!=NULL) {
440      TObjString* fileNameEntry = (TObjString*) list->At(index);
441      if (fileNameEntry!=NULL) {
442         TString fileName = GetFile(sourceFXS, "pedestals",
443                                          fileNameEntry->GetString().Data());
444         TFile *f = TFile::Open(fileName);
445         if (!f) {
446           Log ("Error opening pedestal file.");
447           result =2;
448           break;
449         }
450         AliTPCCalibPedestal *calPed;
451         f->GetObject("tpcCalibPedestal",calPed);
452         if ( !calPed ) {
453           Log ("No pedestal calibration object in file.");
454           result = 2;
455           break;
456         }
457
458         //  replace entries for the sectors available in the present file
459
460         changed=true;
461         for (Int_t sector=0; sector<nSectors; sector++) {
462            AliTPCCalROC *rocPed=calPed->GetCalRocPedestal(sector, kFALSE);
463            if ( rocPed )  calPadPed->SetCalROC(rocPed,sector);
464            AliTPCCalROC *rocRMS=calPed->GetCalRocRMS(sector, kFALSE);
465            if ( rocRMS )  calPadRMS->SetCalROC(rocRMS,sector);
466         }
467       }
468      ++index;
469     }  // while(list)
470 //
471 //  Store updated pedestal entry to OCDB
472 //
473     if (changed) {
474      AliCDBMetaData metaData;
475      metaData.SetBeamPeriod(0);
476      metaData.SetResponsible("Haavard Helstrup");
477      metaData.SetComment("Preprocessor AliTPC data base entries."); 
478  
479      Bool_t storeOK = Store("Calib", "Pedestals", calPadPed, &metaData, 0, kTRUE);
480      if ( !storeOK ) ++result;
481      storeOK = Store("Calib", "PadNoise", calPadRMS, &metaData, 0, kTRUE);
482      if ( !storeOK ) ++result;
483     }
484   } else {
485     Log ("Error: no entries in input file list!");
486     result = 1;
487   }
488
489   return result;
490 }
491
492 //______________________________________________________________________________________________
493
494
495 UInt_t AliTPCPreprocessor::ExtractPulser(Int_t sourceFXS)
496 {
497  //
498  //  Read pulser calibration file from file exchage server
499  //  Keep original entry from OCDB in case no new pulser calibration is available
500  //
501  TObjArray    *pulserObjects=0;
502  AliTPCCalPad *pulserTmean=0;
503  AliTPCCalPad *pulserTrms=0;
504  AliTPCCalPad *pulserQmean=0;
505  AliCDBEntry* entry = GetFromOCDB("Calib", "Pulser");
506  if (entry) pulserObjects = (TObjArray*)entry->GetObject();
507  if ( pulserObjects==NULL ) {
508      Log("AliTPCPreprocsessor: No previous TPC pulser entry available.\n");
509      pulserObjects = new TObjArray;    
510  }
511
512  pulserTmean = (AliTPCCalPad*)pulserObjects->FindObject("PulserTmean");
513  if ( !pulserTmean ) {
514     pulserTmean = new AliTPCCalPad("PulserTmean","PulserTmean");
515     pulserObjects->Add(pulserTmean);
516  }
517  pulserTrms = (AliTPCCalPad*)pulserObjects->FindObject("PulserTrms");
518  if ( !pulserTrms )  { 
519     pulserTrms = new AliTPCCalPad("PulserTrms","PulserTrms");
520     pulserObjects->Add(pulserTrms);
521  }
522  pulserQmean = (AliTPCCalPad*)pulserObjects->FindObject("PulserQmean");
523  if ( !pulserQmean )  { 
524     pulserQmean = new AliTPCCalPad("PulserQmean","PulserQmean");
525     pulserObjects->Add(pulserQmean);
526  }
527
528
529  UInt_t result=0;
530
531  Int_t nSectors = fROC->GetNSectors();
532  TList* list = GetFileSources(sourceFXS,"pulser");
533  
534  if (list && list->GetEntries()>0) {
535
536 //  loop through all files from LDCs
537
538     Bool_t changed=false;
539     UInt_t index = 0;
540     while (list->At(index)!=NULL) {
541      TObjString* fileNameEntry = (TObjString*) list->At(index);
542      if (fileNameEntry!=NULL) {
543         TString fileName = GetFile(sourceFXS, "pulser",
544                                          fileNameEntry->GetString().Data());
545         TFile *f = TFile::Open(fileName);
546         if (!f) {
547           Log ("Error opening pulser file.");
548           result =2;
549           break;
550         }
551         AliTPCCalibPulser *calPulser;
552         f->GetObject("tpcCalibPulser",calPulser);
553         if ( !calPulser ) {
554           Log ("No pulser calibration object in file.");
555           result = 2;
556           break;
557         }
558
559         //  replace entries for the sectors available in the present file
560
561         changed=true;
562         for (Int_t sector=0; sector<nSectors; sector++) {
563            AliTPCCalROC *rocTmean=calPulser->GetCalRocT0(sector);
564            if ( rocTmean )  pulserTmean->SetCalROC(rocTmean,sector);
565            AliTPCCalROC *rocTrms=calPulser->GetCalRocRMS(sector);
566            if ( rocTrms )  pulserTrms->SetCalROC(rocTrms,sector);
567            AliTPCCalROC *rocQmean=calPulser->GetCalRocQ(sector);
568            if ( rocQmean )  pulserQmean->SetCalROC(rocQmean,sector);
569         }
570       }
571      ++index;
572     }  // while(list)
573 //
574 //  Store updated pedestal entry to OCDB
575 //
576     if (changed) {
577      AliCDBMetaData metaData;
578      metaData.SetBeamPeriod(0);
579      metaData.SetResponsible("Haavard Helstrup");
580      metaData.SetComment("Preprocessor AliTPC data base entries.");
581
582      Bool_t storeOK = Store("Calib", "Pulser", pulserObjects, &metaData, 0, kTRUE);
583      if ( !storeOK ) ++result;
584     }  
585   } else {
586     Log ("Error: no entries in input file list!");
587     result = 1;
588   }
589
590   return result;
591 }
592
593 UInt_t AliTPCPreprocessor::ExtractCE(Int_t sourceFXS)
594 {
595  //
596  //  Read Central Electrode file from file exchage server
597  //  Keep original entry from OCDB in case no new CE calibration is available
598  //
599  TObjArray    *ceObjects=0;
600  AliTPCCalPad *ceTmean=0;
601  AliTPCCalPad *ceTrms=0;
602  AliTPCCalPad *ceQmean=0;
603  TObjArray    *rocTtime=0;  
604  TObjArray    *rocQtime=0;  
605
606  AliCDBEntry* entry = GetFromOCDB("Calib", "CE");
607  if (entry) ceObjects = (TObjArray*)entry->GetObject();
608  if ( ceObjects==NULL ) {
609      Log("AliTPCPreprocsessor: No previous TPC central electrode entry available.\n");
610      ceObjects = new TObjArray;    
611  }
612
613  Int_t nSectors = fROC->GetNSectors();
614
615  ceTmean = (AliTPCCalPad*)ceObjects->FindObject("CETmean");
616  if ( !ceTmean ) {
617     ceTmean = new AliTPCCalPad("CETmean","CETmean");
618     ceObjects->Add(ceTmean);
619  }
620  ceTrms = (AliTPCCalPad*)ceObjects->FindObject("CETrms");
621  if ( !ceTrms )  { 
622     ceTrms = new AliTPCCalPad("CETrms","CETrms");
623     ceObjects->Add(ceTrms);
624  }
625  ceQmean = (AliTPCCalPad*)ceObjects->FindObject("CEQmean");
626  if ( !ceQmean )  { 
627     ceQmean = new AliTPCCalPad("CEQmean","CEQmean");
628     ceObjects->Add(ceQmean);
629  }
630  //!new from here please have a look!!!
631  rocTtime = (TObjArray*)ceObjects->FindObject("rocTtime");
632  if ( !rocTtime ) {
633      rocTtime = new TObjArray(nSectors);
634      rocTtime->SetName("rocTtime");
635      ceObjects->Add(rocTtime);
636  }
637  
638  rocQtime = (TObjArray*)ceObjects->FindObject("rocQtime");
639  if ( !rocQtime ) {
640      rocQtime = new TObjArray(nSectors);
641      rocQtime->SetName("rocQtime");
642      ceObjects->Add(rocQtime);
643  }
644
645
646  UInt_t result=0;
647
648  TList* list = GetFileSources(sourceFXS,"CE");
649  
650  if (list && list->GetEntries()>0) {
651
652 //  loop through all files from LDCs
653
654     UInt_t index = 0;
655     while (list->At(index)!=NULL) {
656      TObjString* fileNameEntry = (TObjString*) list->At(index);
657      if (fileNameEntry!=NULL) {
658         TString fileName = GetFile(sourceFXS, "CE",
659                                          fileNameEntry->GetString().Data());
660         TFile *f = TFile::Open(fileName);
661         if (!f) {
662           Log ("Error opening central electrode file.");
663           result =2;
664           break;
665         }
666         AliTPCCalibCE *calCE;
667         f->GetObject("tpcCalibCE",calCE);
668
669         //  replace entries for the sectors available in the present file
670
671         for (Int_t sector=0; sector<nSectors; sector++) {
672            AliTPCCalROC *rocTmean=calCE->GetCalRocT0(sector);
673            if ( rocTmean )  ceTmean->SetCalROC(rocTmean,sector);
674            AliTPCCalROC *rocTrms=calCE->GetCalRocRMS(sector);
675            if ( rocTrms )  ceTrms->SetCalROC(rocTrms,sector);
676            AliTPCCalROC *rocQmean=calCE->GetCalRocQ(sector);
677            if ( rocQmean )  ceQmean->SetCalROC(rocQmean,sector);
678            TGraph *grT=calCE->MakeGraphTimeCE(sector,0,2); // T time graph
679            if ( grT ) rocTtime->AddAt(grT,sector);         
680            TGraph *grQ=calCE->MakeGraphTimeCE(sector,0,3); // Q time graph
681            if ( grQ ) rocTtime->AddAt(grQ,sector);         
682         }
683       }
684      ++index;
685     }  // while(list)
686 //
687 //  Store updated pedestal entry to OCDB
688 //
689     AliCDBMetaData metaData;
690     metaData.SetBeamPeriod(0);
691     metaData.SetResponsible("Haavard Helstrup");
692     metaData.SetComment("Preprocessor AliTPC data base entries.");
693
694     Bool_t storeOK = Store("Calib", "CE", ceObjects, &metaData, 0, kTRUE);
695     if ( !storeOK ) ++result;
696     
697   } else {
698     Log ("Error: no entries!");
699     result = 1;
700   }
701
702   return result;
703 }