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