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