]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/Base/AliTPCPreprocessor.cxx
c289a7d385e795367dd2cb355c9557d938b7d7bd
[u/mrichter/AliRoot.git] / TPC / Base / 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 "AliTPCCalibRaw.h"
31 #include "AliTPCdataQA.h"
32 #include "ARVersion.h"
33 #include "TFile.h"
34 #include "TTree.h"
35 #include "TGraph.h" 
36 #include "TEnv.h"
37 #include "TParameter.h"
38
39 #include <TTimeStamp.h>
40
41 const Int_t kValCutTemp = 100;               // discard temperatures > 100 degrees
42 const Int_t kDiffCutTemp = 5;                // discard temperature differences > 5 degrees
43 const Double_t kHighVoltageDifference = 1e-4; // don't record High Voltage points 
44                                              // differing by less than 1e-4 from
45                                              // previous point.
46 const Double_t kGasCompositionDifference = 1e-3; // don't record Gas Composition points 
47                                              // differing by less than 1e-3 from
48                                              // previous point
49 const TString kPedestalRunType = "PEDESTAL";  // pedestal run identifier
50 const TString kPulserRunType = "PULSER";     // pulser run identifier
51 const TString kPhysicsRunType = "PHYSICS";   // physics run identifier
52 const TString kCosmicRunType = "COSMIC";     // cosmic run identifier
53 const TString kLaserRunType = "LASER";       // laser run identifier
54 const TString kDaqRunType = "DAQ"; // DAQ run identifier
55 const TString kAmandaTemp = "TPC_PT_%d_TEMPERATURE"; // Amanda string for temperature entries
56 const TString kAmandaDDL = "DDL%d";   // Amanda string for list of active DDLs
57 const Int_t  kNumDDL = 216;           // number of TPC DDLs
58 const Int_t  kFirstDDL = 768;         // identifier of first DDL
59 //const Double_t kFitFraction = 0.7;                 // Fraction of DCS sensor fits required              
60 const Double_t kFitFraction = -1.0;          // Don't require minimum number of fits in commissioning run 
61 const Int_t   kNumPressureSensors = 3;    // number of pressure sensors
62 const char* kPressureSensorNames[kNumPressureSensors] = {
63                    "CavernAtmosPressure",
64                    "CavernAtmosPressure2",
65                    "SurfaceAtmosPressure" };
66 const Int_t  kMinCESectors = 32;      // minimum number of sectors (per side)
67                                       // to accept CE calibration
68       
69
70 //
71 // This class is the SHUTTLE preprocessor for the TPC detector.
72 //
73
74 ClassImp(AliTPCPreprocessor)
75
76 //______________________________________________________________________________________________
77 AliTPCPreprocessor::AliTPCPreprocessor(AliShuttleInterface* shuttle) :
78   AliPreprocessor("TPC",shuttle),
79   fConfEnv(0), fTemp(0), fHighVoltage(0), fHighVoltageStat(0), fGoofie(0),
80   fPressure(0), fGasComposition(0), fConfigOK(kTRUE), fROC(0)
81 {
82   // constructor
83   fROC = AliTPCROC::Instance();
84
85   // define run types to be processed
86   
87   AddRunType(kPedestalRunType);
88   AddRunType(kPulserRunType);
89   AddRunType(kPhysicsRunType);
90   AddRunType(kCosmicRunType);
91   AddRunType(kLaserRunType);
92   AddRunType(kDaqRunType);
93   
94 }
95 //______________________________________________________________________________________________
96  AliTPCPreprocessor::AliTPCPreprocessor(const AliTPCPreprocessor&  ) :
97    AliPreprocessor("TPC",0),
98    fConfEnv(0), fTemp(0), fHighVoltage(0), fHighVoltageStat(0), fGoofie(0),
99    fPressure(0), fGasComposition(0), fConfigOK(kTRUE), fROC(0)
100  {
101
102    Fatal("AliTPCPreprocessor", "copy constructor not implemented");
103 //
104 // //  fTemp = new AliTPCSensorTempArray(*(org.fTemp));
105  }
106
107 //______________________________________________________________________________________________
108 AliTPCPreprocessor::~AliTPCPreprocessor()
109 {
110   // destructor
111
112   delete fTemp;
113   delete fHighVoltage;
114   delete fHighVoltageStat;
115   delete fGoofie;
116   delete fPressure;
117 }
118 //______________________________________________________________________________________________
119 AliTPCPreprocessor& AliTPCPreprocessor::operator = (const AliTPCPreprocessor& )
120 {
121   Fatal("operator =", "assignment operator not implemented");
122   return *this;
123 }
124
125
126 //______________________________________________________________________________________________
127 void AliTPCPreprocessor::Initialize(Int_t run, UInt_t startTime,
128         UInt_t endTime)
129 {
130
131   AliPreprocessor::Initialize(run, startTime, endTime);
132
133         AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s", run,
134                 TTimeStamp((time_t)startTime,0).AsString(),
135                 TTimeStamp((time_t)endTime,0).AsString()));
136
137   // Preprocessor configuration
138
139         AliCDBEntry* entry = GetFromOCDB("Config", "Preprocessor");
140         if (entry) fConfEnv = (TEnv*) entry->GetObject();
141         if ( fConfEnv==0 ) {
142            Log("AliTPCPreprocsessor: Preprocessor Config OCDB entry missing.\n");
143            fConfigOK = kFALSE;
144            return;
145         }
146
147   // Temperature sensors
148
149        TTree *confTree = 0;
150
151        TString tempConf = fConfEnv->GetValue("Temperature","ON");
152        tempConf.ToUpper();
153        if (tempConf != "OFF" ) {
154         entry = GetFromOCDB("Config", "Temperature");
155         if (entry) confTree = (TTree*) entry->GetObject();
156         if ( confTree==0 ) {
157            Log("AliTPCPreprocsessor: Temperature Config OCDB entry missing.\n");
158            fConfigOK = kFALSE;
159            return;
160         }
161         fTemp = new AliTPCSensorTempArray(startTime, endTime, confTree, kAmandaTemp);
162         fTemp->SetValCut(kValCutTemp);
163         fTemp->SetDiffCut(kDiffCutTemp);
164        }
165
166   // High voltage measurements
167
168       TString hvConf = fConfEnv->GetValue("HighVoltage","ON");
169       hvConf.ToUpper();
170       if (hvConf != "OFF" ) { 
171         confTree=0;
172         entry=0;
173         entry = GetFromOCDB("Config", "HighVoltage");
174         if (entry) confTree = (TTree*) entry->GetObject();
175         if ( confTree==0 ) {
176            Log("AliTPCPreprocsessor: High Voltage Config OCDB entry missing.\n");
177            fConfigOK = kFALSE;
178            return;
179         }
180         time_t timeStart = (time_t)(((TString)GetRunParameter("DAQ_time_start")).Atoi());
181         time_t timeEnd = (time_t)(((TString)GetRunParameter("DAQ_time_end")).Atoi());
182         fHighVoltage = new AliDCSSensorArray (UInt_t(timeStart), 
183                                             UInt_t(timeEnd), confTree);
184       }
185
186    // High voltage status values
187      
188       TString hvStatConf = fConfEnv->GetValue("HighVoltageStat","ON");
189       hvStatConf.ToUpper();
190       if (hvStatConf != "OFF" ) { 
191         confTree=0;
192         entry=0;
193         entry = GetFromOCDB("Config", "HighVoltageStat");
194         if (entry) confTree = (TTree*) entry->GetObject();
195         if ( confTree==0 ) {
196            Log("AliTPCPreprocsessor: High Voltage Status Config OCDB entry missing.\n");
197            fConfigOK = kFALSE;
198            return;
199         }
200         fHighVoltageStat = new AliDCSSensorArray(startTime, endTime, confTree);
201       }
202
203    // Goofie values
204      
205       TString goofieConf = fConfEnv->GetValue("Goofie","ON");
206       goofieConf.ToUpper();
207       if (goofieConf != "OFF" ) { 
208         confTree=0;
209         entry=0;
210         entry = GetFromOCDB("Config", "Goofie");
211         if (entry) confTree = (TTree*) entry->GetObject();
212         if ( confTree==0 ) {
213            Log("AliTPCPreprocsessor: Goofie Config OCDB entry missing.\n");
214            fConfigOK = kFALSE;
215            return;
216         }
217         fGoofie = new AliDCSSensorArray(startTime, endTime, confTree);
218       }
219
220
221   // Gas composition measurements
222
223       TString gasConf = fConfEnv->GetValue("GasComposition","OFF");
224       gasConf.ToUpper();
225       if (gasConf != "OFF" ) { 
226         confTree=0;
227         entry=0;
228         entry = GetFromOCDB("Config", "GasComposition");
229         if (entry) confTree = (TTree*) entry->GetObject();
230         if ( confTree==0 ) {
231            Log("AliTPCPreprocsessor: Gas Composition Config OCDB entry missing.\n");
232            fConfigOK = kFALSE;
233            return;
234         }
235         time_t timeStart = (time_t)(((TString)GetRunParameter("DAQ_time_start")).Atoi());
236         time_t timeEnd = (time_t)(((TString)GetRunParameter("DAQ_time_end")).Atoi());
237         fGasComposition = new AliDCSSensorArray (UInt_t(timeStart), 
238                                             UInt_t(timeEnd), confTree);
239       }
240    // Pressure values
241      
242        TString runType = GetRunType();
243
244        if( runType == kPhysicsRunType || 
245         runType == kLaserRunType ) {    
246        TString pressureConf = fConfEnv->GetValue("Pressure","ON");
247        pressureConf.ToUpper();
248        if (pressureConf != "OFF" ) { 
249          TClonesArray * array = new TClonesArray("AliDCSSensor",kNumPressureSensors); 
250          for(Int_t j = 0; j < kNumPressureSensors; j++) {
251            AliDCSSensor * sens = new ((*array)[j])AliDCSSensor;
252            sens->SetStringID(kPressureSensorNames[j]);
253          }
254          fPressure = new AliDCSSensorArray(startTime, endTime, array);
255        }
256      }
257 }
258
259 //______________________________________________________________________________________________
260 UInt_t AliTPCPreprocessor::Process(TMap* dcsAliasMap)
261 {
262   // Fills data into TPC calibrations objects
263
264   // Amanda servers provide information directly through dcsAliasMap
265
266   
267   if (!fConfigOK) return 9;
268   UInt_t result = 0;
269   TObjArray *resultArray = new TObjArray();
270   TString errorHandling = fConfEnv->GetValue("ErrorHandling","ON");
271   errorHandling.ToUpper();
272   TObject * status;
273
274   UInt_t dcsResult=0;
275   if (!dcsAliasMap) { 
276      dcsResult=1;
277   } else {
278      if (dcsAliasMap->GetEntries() == 0 ) dcsResult=1;
279   }  
280   status = new TParameter<int>("dcsResult",dcsResult);
281   resultArray->Add(status);
282
283
284   TString runType = GetRunType();
285
286   if ( dcsResult == 0 ) {
287
288   // Temperature sensors are processed by AliTPCCalTemp
289
290     TString tempConf = fConfEnv->GetValue("Temperature","ON");
291     tempConf.ToUpper();
292     if (tempConf != "OFF" ) {
293       UInt_t tempResult = MapTemperature(dcsAliasMap);
294       if ( tempConf != "TRY") result+=tempResult;
295       status = new TParameter<int>("tempResult",tempResult);
296       resultArray->Add(status);
297     }
298
299   // High Voltage recordings
300
301
302     TString hvConf = fConfEnv->GetValue("HighVoltage","ON");
303     hvConf.ToUpper();
304     if (hvConf != "OFF" ) { 
305      UInt_t hvResult = MapHighVoltage(dcsAliasMap);
306      if (hvConf != "TRY") result+=hvResult;
307      status = new TParameter<int>("hvResult",hvResult);
308      resultArray->Add(status);
309     }
310
311     // Goofie values
312
313
314     TString goofieConf = fConfEnv->GetValue("Goofie","ON");
315     goofieConf.ToUpper();
316     if (goofieConf != "OFF" ) { 
317      UInt_t goofieResult = MapGoofie(dcsAliasMap);
318      if (goofieConf != "TRY") result+=goofieResult;
319      status = new TParameter<int>("goofieResult",goofieResult);
320      resultArray->Add(status);
321     }
322
323
324     // Gas composition recordings
325
326
327     TString gasConf = fConfEnv->GetValue("GasComposition","OFF");
328     gasConf.ToUpper();
329     if (gasConf != "OFF" ) { 
330      UInt_t gasResult = MapGasComposition(dcsAliasMap);
331      if (gasConf != "TRY") result+=gasResult;
332      status = new TParameter<int>("gasResult",gasResult);
333      resultArray->Add(status);
334     }
335
336     // Pressure values
337
338     if( runType == kPhysicsRunType || 
339       runType == kLaserRunType ) {    
340
341       TString pressureConf = fConfEnv->GetValue("Pressure","ON");
342       pressureConf.ToUpper();
343       if (pressureConf != "OFF" ) { 
344        UInt_t pressureResult = MapPressure(dcsAliasMap);
345        status = new TParameter<int>("pressureResult",pressureResult);
346        resultArray->Add(status);
347       }
348     }
349   }
350   // Other calibration information will be retrieved through FXS files
351   //  examples:
352   //    TList* fileSourcesDAQ = GetFile(AliShuttleInterface::kDAQ, "pedestals");
353   //    const char* fileNamePed = GetFile(AliShuttleInterface::kDAQ, "pedestals", "LDC1");
354   //
355   //    TList* fileSourcesHLT = GetFile(AliShuttleInterface::kHLT, "calib");
356   //    const char* fileNameHLT = GetFile(AliShuttleInterface::kHLT, "calib", "LDC1");
357
358   // pedestal entries
359
360   if(runType == kPedestalRunType) {
361     Int_t numSources = 1;
362     Int_t pedestalSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
363     TString source = fConfEnv->GetValue("Pedestal","DAQ");
364     source.ToUpper();
365     if (source != "OFF" ) { 
366      if ( source == "HLT") pedestalSource[0] = AliShuttleInterface::kHLT;
367      if (!GetHLTStatus()) pedestalSource[0] = AliShuttleInterface::kDAQ;
368      if (source == "HLTDAQ" ) {
369          numSources=2;
370          pedestalSource[0] = AliShuttleInterface::kHLT;
371          pedestalSource[1] = AliShuttleInterface::kDAQ;
372      }
373      if (source == "DAQHLT" ) numSources=2;
374      UInt_t pedestalResult=0;
375      for (Int_t i=0; i<numSources; i++ ) {      
376        pedestalResult = ExtractPedestals(pedestalSource[i]);
377        if ( pedestalResult == 0 ) break;
378      }
379      result += pedestalResult;
380      status = new TParameter<int>("pedestalResult",pedestalResult);
381      resultArray->Add(status);
382     }
383   }
384
385   // pulser trigger processing
386
387   if(runType == kPulserRunType) {
388     Int_t numSources = 1;
389     Int_t pulserSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
390     TString source = fConfEnv->GetValue("Pulser","DAQ");
391     source.ToUpper();
392     if ( source != "OFF") { 
393      if ( source == "HLT") pulserSource[0] = AliShuttleInterface::kHLT;
394      if (!GetHLTStatus()) pulserSource[0] = AliShuttleInterface::kDAQ;
395      if (source == "HLTDAQ" ) {
396          numSources=2;
397          pulserSource[0] = AliShuttleInterface::kHLT;
398          pulserSource[1] = AliShuttleInterface::kDAQ;
399      }
400      if (source == "DAQHLT" ) numSources=2;
401      if (source == "TRY" ) numSources=2;
402      UInt_t pulserResult=0;
403      for (Int_t i=0; i<numSources; i++ ) {      
404        pulserResult = ExtractPulser(pulserSource[i]);
405        if ( pulserResult == 0 ) break;
406      }
407      if (source != "TRY") result += pulserResult;
408      status = new TParameter<int>("pulserResult",pulserResult);
409      resultArray->Add(status);
410     }
411   }
412
413
414 // raw calibration processing
415
416   if(runType == kPhysicsRunType) {
417     Int_t numSources = 1;
418     Int_t rawSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
419     TString source = fConfEnv->GetValue("Raw","DAQ");
420     source.ToUpper();
421     if ( source != "OFF") { 
422      if ( source == "HLT") rawSource[0] = AliShuttleInterface::kHLT;
423      if (!GetHLTStatus()) rawSource[0] = AliShuttleInterface::kDAQ;
424      if (source == "HLTDAQ" ) {
425          numSources=2;
426          rawSource[0] = AliShuttleInterface::kHLT;
427          rawSource[1] = AliShuttleInterface::kDAQ;
428      }
429      if (source == "DAQHLT" ) numSources=2;
430      if (source == "TRY" ) numSources=2;
431      UInt_t rawResult=0;
432      for (Int_t i=0; i<numSources; i++ ) {      
433        rawResult = ExtractRaw(rawSource[i]);
434        if ( rawResult == 0 ) break;
435      }
436      if (source != "TRY" )result += rawResult;
437      status = new TParameter<int>("rawResult",rawResult);
438      resultArray->Add(status);
439     }
440   }
441
442
443   // Altro configuration
444
445
446   TString altroConf = fConfEnv->GetValue("AltroConf","ON");
447   altroConf.ToUpper();
448   if (altroConf != "OFF" ) { 
449    UInt_t altroResult = ExtractAltro(AliShuttleInterface::kDCS,dcsAliasMap);
450    if (altroConf != "TRY" ) result+=altroResult;
451    status = new TParameter<int>("altroResult",altroResult);
452    resultArray->Add(status);
453  }
454
455
456   // Central Electrode processing
457
458   if( runType == kPhysicsRunType || 
459       runType == kLaserRunType ) {    
460
461     Int_t numSources = 1;
462     Int_t ceSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
463     TString source = fConfEnv->GetValue("CE","DAQ");
464     source.ToUpper();
465     if ( source != "OFF" ) { 
466      if ( source == "HLT") ceSource[0] = AliShuttleInterface::kHLT;
467      if (!GetHLTStatus()) ceSource[0] = AliShuttleInterface::kDAQ;
468      if (source == "HLTDAQ" ) {
469         numSources=2;
470         ceSource[0] = AliShuttleInterface::kHLT;
471         ceSource[1] = AliShuttleInterface::kDAQ;
472      }
473      if (source == "DAQHLT" ) numSources=2;
474      if (source == "TRY" ) numSources=2;
475      UInt_t ceResult=0;
476      for (Int_t i=0; i<numSources; i++ ) {      
477        ceResult = ExtractCE(ceSource[i]);
478        if ( ceResult == 0 ) break;
479      }
480
481    // only flag error if CE result is missing from LASER runs
482    //    -- for PHYSICS run do CE processing if data available
483    
484      if ( runType == kLaserRunType && source != "TRY"  && ceResult<10 ) result += ceResult;
485      status = new TParameter<int>("ceResult",ceResult);
486      resultArray->Add(status);
487
488     numSources = 1;
489     Int_t qaSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
490     source = fConfEnv->GetValue("QA","DAQ");
491     source.ToUpper();
492     if ( source != "OFF" ) { 
493      if ( source == "HLT") qaSource[0] = AliShuttleInterface::kHLT;
494      if (!GetHLTStatus()) qaSource[0] = AliShuttleInterface::kDAQ;
495      if (source == "HLTDAQ" ) {
496         numSources=2;
497         qaSource[0] = AliShuttleInterface::kHLT;
498         qaSource[1] = AliShuttleInterface::kDAQ;
499      }
500      if (source == "DAQHLT" ) numSources=2;
501      if (source == "TRY" ) numSources=2;
502      UInt_t qaResult=0;
503      for (Int_t i=0; i<numSources; i++ ) {      
504        qaResult = ExtractQA(qaSource[i]);
505        if ( qaResult == 0 ) break;
506      }
507 //     result += qaResult;
508      if ( qaResult !=0 ) Log ("ExtractQA failed, no QA entry available.");
509      status = new TParameter<int>("qaResult",qaResult);
510      resultArray->Add(status);
511     }
512    }
513   }
514   
515 // Store component status to OCDB
516
517   AliCDBMetaData metaData;
518   metaData.SetBeamPeriod(0);
519   metaData.SetResponsible("Haavard Helstrup");
520   metaData.SetAliRootVersion(ALIROOT_BRANCH);
521   metaData.SetComment("Preprocessor AliTPC status.");
522   Bool_t storeOK = Store("Calib", "PreprocStatus", resultArray, &metaData, 0,  kFALSE);
523   if (!storeOK) Log ("Unable to store preprocessor status entry");
524  
525   resultArray->Delete();
526   delete resultArray;
527
528   if (errorHandling == "OFF" ) return 0;
529   return result;
530   
531 }
532 //______________________________________________________________________________________________
533 UInt_t AliTPCPreprocessor::MapTemperature(TMap* dcsAliasMap)
534 {
535
536    // extract DCS temperature maps. Perform fits to save space
537
538   UInt_t result=0;
539   TMap *map = fTemp->ExtractDCS(dcsAliasMap);
540   if (map) {
541     fTemp->MakeSplineFit(map);
542     Double_t fitFraction = 1.0*fTemp->NumFits()/fTemp->NumSensors(); 
543     if (fitFraction > kFitFraction ) {
544       AliInfo(Form("Temperature values extracted, fits performed.\n"));
545     } else { 
546       Log ("Too few temperature maps fitted. \n");
547       result = 9;
548     }
549   } else {
550     Log("No temperature map extracted. \n");
551     result=9;
552   }
553   delete map;
554   // Now store the final CDB file
555
556   if ( result == 0 ) {
557         AliCDBMetaData metaData;
558         metaData.SetBeamPeriod(0);
559         metaData.SetResponsible("Haavard Helstrup");
560         metaData.SetAliRootVersion(ALIROOT_BRANCH);
561         metaData.SetComment("Preprocessor AliTPC data base entries.");
562
563         Bool_t storeOK = Store("Calib", "Temperature", fTemp, &metaData, 0, kFALSE);
564         if ( !storeOK )  result=1;
565
566    }
567
568    return result;
569
570 }
571 //______________________________________________________________________________________________
572 UInt_t AliTPCPreprocessor::MapPressure(TMap* dcsAliasMap)
573 {
574
575    // extract DCS pressure maps. Perform fits to save space
576
577   UInt_t result=0;
578   TMap *map = fPressure->ExtractDCS(dcsAliasMap);
579   if (map) {
580     fPressure->MakeSplineFit(map);
581     Double_t fitFraction = 1.0*fPressure->NumFits()/fPressure->NumSensors(); 
582     if (fitFraction > kFitFraction ) {
583       AliInfo(Form("Pressure values extracted, fits performed.\n"));
584     } else { 
585       Log ("Too few pressure maps fitted. \n");
586       result = 9;
587     }
588   } else {
589     Log("No pressure map extracted. \n");
590     result=9;
591   }
592   delete map;
593   return result;
594 }
595
596 //______________________________________________________________________________________________
597 UInt_t AliTPCPreprocessor::MapHighVoltage(TMap* dcsAliasMap)
598 {
599
600    // extract DCS HV maps. Perform fits to save space
601
602   UInt_t result=0;
603   TMap *map = fHighVoltage->ExtractDCS(dcsAliasMap);
604   if (map) {
605     fHighVoltage->ClearFit();
606     fHighVoltage->RemoveGraphDuplicates(kHighVoltageDifference);
607                // don't keep new point if too similar to previous one
608     fHighVoltage->SetGraph(map);
609   } else {
610     Log("No high voltage recordings extracted. \n");
611     result=9;
612   }
613   delete map;
614
615   TString hvStatConf = fConfEnv->GetValue("HighVoltageStat","ON");
616   hvStatConf.ToUpper();
617   if (hvStatConf != "OFF" ) { 
618     TMap *map2 = fHighVoltageStat->ExtractDCS(dcsAliasMap);
619     if (map2) {
620       fHighVoltageStat->ClearFit();
621       fHighVoltageStat->SetGraph(map2);
622     } else {
623        Log("No high voltage status recordings extracted. \n");
624       result=9;
625     }
626     delete map2;
627
628     // add status maps to high voltage sensor array
629
630     fHighVoltage->AddSensors(fHighVoltageStat);
631    }
632   // Now store the final CDB file
633
634   if ( result == 0 ) {
635         AliCDBMetaData metaData;
636         metaData.SetBeamPeriod(0);
637         metaData.SetResponsible("Haavard Helstrup");
638         metaData.SetAliRootVersion(ALIROOT_BRANCH);
639         metaData.SetComment("Preprocessor AliTPC data base entries.");
640
641         Bool_t storeOK = Store("Calib", "HighVoltage", fHighVoltage, &metaData, 0, kFALSE);
642         if ( !storeOK )  result=1;
643
644    }
645
646    return result;
647
648 }
649
650 //______________________________________________________________________________________________
651 UInt_t AliTPCPreprocessor::MapGoofie(TMap* dcsAliasMap)
652 {
653
654    // extract DCS Goofie maps. Do not perform fits (low update rate)
655
656   UInt_t result=0;
657
658   TMap *map = fGoofie->ExtractDCS(dcsAliasMap);
659   if (map) {
660     fGoofie->ClearFit();
661     fGoofie->SetGraph(map);
662   } else {
663     Log("No Goofie recordings extracted. \n");
664     result=9;
665   }
666   delete map;
667
668   // Now store the final CDB file
669
670   if ( result == 0 ) {
671         AliCDBMetaData metaData;
672         metaData.SetBeamPeriod(0);
673         metaData.SetResponsible("Haavard Helstrup");
674         metaData.SetAliRootVersion(ALIROOT_BRANCH);
675         metaData.SetComment("Preprocessor AliTPC data base entries.");
676
677         Bool_t storeOK = Store("Calib", "Goofie", fGoofie, &metaData, 0, kFALSE);
678         if ( !storeOK )  result=1;
679
680    }
681
682    return result;
683
684 }
685
686 //______________________________________________________________________________________________
687 UInt_t AliTPCPreprocessor::MapGasComposition(TMap* dcsAliasMap)
688 {
689
690    // extract DCS HV maps. Perform fits to save space
691
692   UInt_t result=0;
693   TMap *map = fGasComposition->ExtractDCS(dcsAliasMap);
694   if (map) {
695     fHighVoltage->ClearFit();
696     fHighVoltage->RemoveGraphDuplicates(kGasCompositionDifference);
697                // don't keep new point if too similar to previous one
698     fHighVoltage->SetGraph(map);
699   } else {
700     Log("No gas composition recordings extracted. \n");
701     result=9;
702   }
703   delete map;
704
705   // Now store the final CDB file
706
707   if ( result == 0 ) {
708         AliCDBMetaData metaData;
709         metaData.SetBeamPeriod(0);
710         metaData.SetResponsible("Haavard Helstrup");
711         metaData.SetAliRootVersion(ALIROOT_BRANCH);
712         metaData.SetComment("Preprocessor AliTPC data base entries.");
713
714         Bool_t storeOK = Store("Calib", "GasComposition", fGasComposition, &metaData, 0, kFALSE);
715         if ( !storeOK )  result=1;
716
717    }
718
719    return result;
720
721 }
722
723
724 //______________________________________________________________________________________________
725
726 UInt_t AliTPCPreprocessor::ExtractPedestals(Int_t sourceFXS)
727 {
728  //
729  //  Read pedestal file from file exchage server
730  //  Keep original entry from OCDB in case no new pedestals are available
731  //
732  AliTPCCalPad *calPadPedOCDB=0;
733  AliCDBEntry* entry = GetFromOCDB("Calib", "Pedestals");
734  if (entry) calPadPedOCDB = (AliTPCCalPad*)entry->GetObject();
735  if ( calPadPedOCDB==NULL ) {
736      Log("AliTPCPreprocsessor: No previous TPC pedestal entry available.\n");
737  }
738
739  AliTPCCalPad *calPadRMSOCDB=0;
740  entry = GetFromOCDB("Calib", "PadNoise");
741  if (entry) calPadRMSOCDB = (AliTPCCalPad*)entry->GetObject();
742  if ( calPadRMSOCDB==NULL ) {
743      Log("AliTPCPreprocsessor: No previous TPC noise entry available.\n");
744  }
745
746   AliTPCCalPad* calPadPed = new AliTPCCalPad("PedestalsMean","PedestalsMean");
747   AliTPCCalPad* calPadRMS = new AliTPCCalPad("PedestalsRMS","PedestalsRMS");
748
749
750  UInt_t result=0;
751
752  Int_t nSectors = fROC->GetNSectors();
753  TVectorD foundSectorsPed(nSectors);
754  foundSectorsPed=0;
755  TVectorD foundSectorsRMS(nSectors);
756  foundSectorsRMS=0;
757
758  TList* list = GetFileSources(sourceFXS,"pedestals");
759  
760  if (list && list->GetEntries()>0) {
761
762 //  loop through all files from LDCs
763
764     Bool_t changed=false;
765     UInt_t index = 0;
766     while (list->At(index)!=NULL) {
767      TObjString* fileNameEntry = (TObjString*) list->At(index);
768      if (fileNameEntry!=NULL) {
769         TString fileName = GetFile(sourceFXS, "pedestals",
770                                          fileNameEntry->GetString().Data());
771         TFile *f = TFile::Open(fileName);
772         if (!f) {
773           Log ("Error opening pedestal file.");
774           result =2;
775           break;
776         }
777         AliTPCCalibPedestal *calPed;
778         f->GetObject("tpcCalibPedestal",calPed);
779         if ( !calPed ) {
780           Log ("No pedestal calibration object in file.");
781           result = 2;
782           break;
783         }
784
785         //  replace entries for the sectors available in the present file
786
787         changed=true;
788         for (Int_t sector=0; sector<nSectors; sector++) {
789            AliTPCCalROC *rocPed=calPed->GetCalRocPedestal(sector, kFALSE);
790            if ( rocPed ) {
791                AliTPCCalROC* roc=calPadPed->GetCalROC(sector);
792                roc->Add(rocPed,1);
793                foundSectorsPed[sector]++;
794            }
795            AliTPCCalROC *rocRMS=calPed->GetCalRocRMS(sector, kFALSE);
796            if ( rocRMS )  {
797                AliTPCCalROC* roc=calPadRMS->GetCalROC(sector);
798                roc->Add(rocRMS,1);
799                foundSectorsRMS[sector]++;
800            } 
801         }
802         delete calPed; 
803         f->Close();
804       }
805      ++index;
806     }  // while(list)
807
808 //
809 //  Check if calibration is complete -- otherwise take from old OCDB entry
810
811 // inner sectors
812     for (Int_t sector=0; sector<nSectors/2; sector++) {
813        if (foundSectorsPed[sector] < 1 ) {
814           if (calPadPedOCDB) {
815             AliTPCCalROC *rocOCDB=calPadPedOCDB->GetCalROC(sector);
816             calPadPed->SetCalROC(rocOCDB,sector);
817           } else {
818             const int mess_length=100;
819             char message[mess_length];
820             snprintf(message,mess_length,"Missing pedestals for sector %d - also not available from previous OCDB entry.\n",
821                      sector);
822             Log (message);
823             result = 2;
824           }
825        }
826        if (foundSectorsRMS[sector] < 1 ) {
827           if (calPadRMSOCDB) {
828             AliTPCCalROC *rocOCDB=calPadRMSOCDB->GetCalROC(sector);
829             calPadRMS->SetCalROC(rocOCDB,sector);
830           } else {
831             const int mess_length=100;
832             char message[mess_length];
833             snprintf(message,mess_length,"Missing pedestal RMS for sector %d - also not available from previous OCDB entry.\n",
834                      sector);
835             Log (message);
836             result = 2;
837           }       
838        }
839     }
840
841 // outer sectors -- two updates needed
842
843     for (Int_t sector=nSectors/2; sector<nSectors; sector++) {
844        if (foundSectorsPed[sector] < 2 ) {
845           if (calPadPedOCDB) {
846             AliTPCCalROC *rocOCDB=calPadPedOCDB->GetCalROC(sector);
847             calPadPed->SetCalROC(rocOCDB,sector);
848           } else {
849             const int mess_length=100;
850             char message[mess_length];
851             snprintf(message,mess_length,"Missing pedestals for sector %d - also not available from previous OCDB entry.\n",
852                      sector);
853             Log (message);
854             result = 2;
855           }
856
857        }
858        if (foundSectorsRMS[sector] < 2 ) {
859           if (calPadRMSOCDB) {
860             AliTPCCalROC *rocOCDB=calPadRMSOCDB->GetCalROC(sector);
861             calPadRMS->SetCalROC(rocOCDB,sector);
862           } else {
863             const int mess_length=100;
864             char message[mess_length];
865             snprintf(message,mess_length,"Missing pedestal RMS for sector %d - also not available from previous OCDB entry.\n",
866                      sector);
867             Log (message);
868             result = 2;
869           }       
870
871        }       
872     }
873
874     
875 //
876 //  Store updated pedestal entry to OCDB
877 //
878     if (changed) {
879      AliCDBMetaData metaData;
880      metaData.SetBeamPeriod(0);
881      metaData.SetResponsible("Haavard Helstrup");
882      metaData.SetAliRootVersion(ALIROOT_BRANCH);
883      metaData.SetComment("Preprocessor AliTPC data base entries."); 
884  
885      Bool_t storeOK = Store("Calib", "Pedestals", calPadPed, &metaData, 0, kTRUE);
886      if ( !storeOK ) ++result;
887      storeOK = Store("Calib", "PadNoise", calPadRMS, &metaData, 0, kTRUE);
888      if ( !storeOK ) ++result;
889     }
890   } else {
891     Log ("Error: no entries in pedestal file list!");
892     result = 1;
893   }
894
895   delete calPadPed;
896   delete calPadRMS;
897   delete calPadPedOCDB;
898   delete calPadRMSOCDB;
899
900   return result;
901 }
902
903 //______________________________________________________________________________________________
904
905
906 UInt_t AliTPCPreprocessor::ExtractPulser(Int_t sourceFXS)
907 {
908  //
909  //  Read pulser calibration file from file exchage server
910  //  Keep original entry from OCDB in case no new pulser calibration is available
911  //
912  
913  TObjArray *pulserObjects = new TObjArray;
914  TObjArray *pulserObjectsOCDB=0;   
915   
916  AliCDBEntry* entry = GetFromOCDB("Calib", "Pulser");
917  if (entry) pulserObjectsOCDB = (TObjArray*)entry->GetObject();
918  if ( pulserObjectsOCDB==NULL ) {
919      Log("AliTPCPreprocsessor: No previous TPC pulser entry available.\n");
920  }
921
922  AliTPCCalPad *pulserTmean = new AliTPCCalPad("PulserTmean","PulserTmean");
923  pulserObjects->Add(pulserTmean);
924  
925  AliTPCCalPad *pulserTrms = new AliTPCCalPad("PulserTrms","PulserTrms");
926  pulserObjects->Add(pulserTrms);
927  
928  AliTPCCalPad *pulserQmean = new AliTPCCalPad("PulserQmean","PulserQmean");
929  pulserObjects->Add(pulserQmean);
930  
931
932
933  UInt_t result=0;
934
935  Int_t nSectors = fROC->GetNSectors();
936  TVectorD foundTmean(nSectors);
937  foundTmean=0;
938  TVectorD foundTrms(nSectors);
939  foundTrms=0;
940  TVectorD foundQmean(nSectors);
941  foundQmean=0;
942
943
944  TList* list = GetFileSources(sourceFXS,"pulser");
945  
946  if (list && list->GetEntries()>0) {
947
948 //  loop through all files from LDCs
949
950     Bool_t changed=false;
951     UInt_t index = 0;
952     while (list->At(index)!=NULL) {
953      TObjString* fileNameEntry = (TObjString*) list->At(index);
954      if (fileNameEntry!=NULL) {
955         TString fileName = GetFile(sourceFXS, "pulser",
956                                          fileNameEntry->GetString().Data());
957         TFile *f = TFile::Open(fileName);
958         if (!f) {
959           Log ("Error opening pulser file.");
960           result =2;
961           break;
962         }
963         AliTPCCalibPulser *calPulser;
964         f->GetObject("tpcCalibPulser",calPulser);
965         if ( !calPulser ) {
966           Log ("No pulser calibration object in file.");
967           result = 2;
968           break;
969         }
970
971         //  replace entries for the sectors available in the present file
972
973         changed=true;
974         for (Int_t sector=0; sector<nSectors; sector++) {
975            AliTPCCalROC *rocTmean=calPulser->GetCalRocT0(sector);
976            if ( rocTmean ) {
977                AliTPCCalROC* roc=pulserTmean->GetCalROC(sector);
978                roc->Add(rocTmean,1);
979                foundTmean[sector]++;
980            }
981            AliTPCCalROC *rocTrms=calPulser->GetCalRocRMS(sector);
982            if ( rocTrms ) {
983                AliTPCCalROC* roc=pulserTrms->GetCalROC(sector);
984                roc->Add(rocTrms,1);
985                foundTrms[sector]++;
986            }
987            AliTPCCalROC *rocQmean=calPulser->GetCalRocQ(sector);
988            if ( rocQmean ) {
989                AliTPCCalROC* roc=pulserQmean->GetCalROC(sector);
990                roc->Add(rocQmean,1);
991                foundQmean[sector]++;
992            }
993         }
994        delete calPulser;
995        f->Close();
996       }
997      ++index;
998     }  // while(list)
999
1000
1001
1002
1003 if (pulserObjectsOCDB) {
1004   AliTPCCalPad* pulserTmeanOCDB = (AliTPCCalPad*)pulserObjectsOCDB->FindObject("PulserTmean");
1005   AliTPCCalPad* pulserTrmsOCDB  = (AliTPCCalPad*)pulserObjectsOCDB->FindObject("PulserTrms");
1006   AliTPCCalPad* pulserQmeanOCDB = (AliTPCCalPad*)pulserObjectsOCDB->FindObject("PulserQmean");
1007  
1008 //
1009 //  Check if calibration is complete -- otherwise take from old OCDB entry
1010
1011 // inner sectors
1012     for (Int_t sector=0; sector<nSectors/2; sector++) {
1013        if (foundTmean[sector] < 1 ) {
1014          if (pulserTmeanOCDB) {
1015           AliTPCCalROC* rocOCDB = pulserTmeanOCDB->GetCalROC(sector);
1016           if ( rocOCDB ) pulserTmean->SetCalROC(rocOCDB,sector);
1017          }
1018        }
1019        if (foundTrms[sector] < 1 ) {
1020          if (pulserTrmsOCDB) {
1021           AliTPCCalROC* rocOCDB = pulserTrmsOCDB->GetCalROC(sector);
1022           if ( rocOCDB ) pulserTrms->SetCalROC(rocOCDB,sector);
1023          }
1024        }
1025        if (foundQmean[sector] < 1 ) {
1026          if (pulserQmeanOCDB) {
1027           AliTPCCalROC* rocOCDB = pulserQmeanOCDB->GetCalROC(sector);
1028           if ( rocOCDB ) pulserQmean->SetCalROC(rocOCDB,sector);
1029          }
1030        }
1031     }
1032   
1033 // outer sectors -- two updates needed
1034
1035     for (Int_t sector=0; sector<nSectors/2; sector++) {
1036        if (foundTmean[sector] < 2 ) {
1037          if (pulserTmeanOCDB) {
1038           AliTPCCalROC* rocOCDB = pulserTmeanOCDB->GetCalROC(sector);
1039           if ( rocOCDB ) pulserTmean->SetCalROC(rocOCDB,sector);
1040          }
1041        }
1042        if (foundTrms[sector] < 2 ) {
1043          if (pulserTrmsOCDB) {
1044           AliTPCCalROC* rocOCDB = pulserTrmsOCDB->GetCalROC(sector);
1045           if ( rocOCDB ) pulserTrms->SetCalROC(rocOCDB,sector);
1046          }
1047        }
1048        if (foundQmean[sector] < 2 ) {
1049          if (pulserQmeanOCDB) {
1050           AliTPCCalROC* rocOCDB = pulserQmeanOCDB->GetCalROC(sector);
1051           if ( rocOCDB ) pulserQmean->SetCalROC(rocOCDB,sector);
1052          }
1053        }
1054
1055     }
1056 }
1057
1058 //
1059 //  Store updated pedestal entry to OCDB
1060 //
1061     if (changed) {
1062      AliCDBMetaData metaData;
1063      metaData.SetBeamPeriod(0);
1064      metaData.SetResponsible("Haavard Helstrup");
1065      metaData.SetAliRootVersion(ALIROOT_BRANCH);
1066      metaData.SetComment("Preprocessor AliTPC data base entries.");
1067
1068      Bool_t storeOK = Store("Calib", "Pulser", pulserObjects, &metaData, 0, kTRUE);
1069      if ( !storeOK ) ++result;
1070     }  
1071   } else {
1072     Log ("Error: no entries in pulser file list!");
1073     result = 1;
1074   }
1075   pulserObjects->Delete();
1076   delete pulserObjects;
1077   if (pulserObjectsOCDB) pulserObjectsOCDB->Delete();
1078   delete pulserObjectsOCDB;
1079
1080   return result;
1081 }
1082
1083 //______________________________________________________________________________________________
1084
1085
1086 UInt_t AliTPCPreprocessor::ExtractRaw(Int_t sourceFXS)
1087 {
1088  //
1089  //  Read Raw calibration file from file exchage server
1090  //
1091  
1092  UInt_t result=0;
1093  TObjArray* rawArray = new TObjArray;
1094
1095  TList* list = GetFileSources(sourceFXS,"tpcCalibRaw");
1096  
1097  if (list && list->GetEntries()>0) {
1098
1099 //  loop through all files
1100
1101     UInt_t index = 0;
1102     while (list->At(index)!=NULL) {
1103      TObjString* fileNameEntry = (TObjString*) list->At(index);
1104      if (fileNameEntry!=NULL) {
1105         TString fileName = GetFile(sourceFXS, "tpcCalibRaw",
1106                                          fileNameEntry->GetString().Data());
1107         TFile *f = TFile::Open(fileName);
1108         if (!f) {
1109           Log ("Error opening raw file.");
1110           result =2;
1111           break;
1112         }
1113         AliTPCCalibRaw *calRaw;
1114         f->GetObject("tpcCalibRaw",calRaw);
1115         if ( !calRaw ) {
1116           Log ("No raw calibration object in file.");
1117           result = 2;
1118           break;
1119         }
1120        rawArray->Add(calRaw);
1121        f->Close();
1122       }
1123      ++index;
1124     }  // while(list)
1125 //
1126 //  Store updated pedestal entry to OCDB
1127 //
1128      AliCDBMetaData metaData;
1129      metaData.SetBeamPeriod(0);
1130      metaData.SetResponsible("Haavard Helstrup");
1131      metaData.SetAliRootVersion(ALIROOT_BRANCH);
1132      metaData.SetComment("Preprocessor AliTPC data base entries.");
1133
1134      Bool_t storeOK = Store("Calib", "Raw", rawArray, &metaData, 0, kTRUE);
1135      if ( !storeOK ) ++result;
1136   } else {
1137     Log ("Error: no entries in raw file list!");
1138     result = 1;
1139   }
1140   
1141   rawArray->Delete();
1142   delete rawArray;
1143
1144   return result;
1145 }
1146 //______________________________________________________________________________________________
1147
1148 UInt_t AliTPCPreprocessor::ExtractCE(Int_t sourceFXS)
1149 {
1150  //
1151  //  Read Central Electrode file from file exchage server
1152  //
1153  //
1154   AliTPCCalPad *ceTmean=0;
1155   AliTPCCalPad *ceTrms=0;
1156   AliTPCCalPad *ceQmean=0;
1157   TObjArray    *rocTtime=0;
1158   TObjArray    *rocQtime=0;
1159   
1160   TObjArray    *ceObjects= new TObjArray;
1161   
1162   
1163   Int_t nSectors = fROC->GetNSectors();
1164   
1165   ceTmean = new AliTPCCalPad("CETmean","CETmean");
1166   ceObjects->Add(ceTmean);
1167   
1168   ceTrms = new AliTPCCalPad("CETrms","CETrms");
1169   ceObjects->Add(ceTrms);
1170   
1171   ceQmean = new AliTPCCalPad("CEQmean","CEQmean");
1172   ceObjects->Add(ceQmean);
1173   
1174   rocTtime = new TObjArray(nSectors+2);   // also make room for A and C side average
1175   rocTtime->SetName("rocTtime");
1176   ceObjects->Add(rocTtime);
1177   
1178   rocQtime = new TObjArray(nSectors);
1179   rocQtime->SetName("rocQtime");
1180   ceObjects->Add(rocQtime);
1181
1182   //=== new part
1183   TObjArray *arrFitGraphs=new TObjArray;
1184   arrFitGraphs->SetName("ceFitsDrift");
1185   ceObjects->Add(arrFitGraphs);
1186   
1187 // Temperature maps 
1188   
1189   if (fTemp) {
1190     AliTPCSensorTempArray *tempMap = new AliTPCSensorTempArray(*fTemp);
1191     tempMap->SetNameTitle("TempMap","TempMap");
1192     ceObjects->Add(tempMap);
1193   }
1194   
1195 // Pressure maps
1196   
1197   if (fPressure) {
1198     AliDCSSensor *sensor=0, *sensorCopy=0;
1199     for (Int_t isensor=0; isensor<kNumPressureSensors; ++isensor ) {
1200       sensor = fPressure->GetSensor(kPressureSensorNames[isensor]);
1201       if (sensor) {
1202         sensorCopy = new AliDCSSensor(*sensor);
1203         sensorCopy->SetNameTitle(kPressureSensorNames[isensor],kPressureSensorNames[isensor]);
1204         ceObjects->Add(sensorCopy);
1205       }
1206     }
1207   }
1208   
1209   UInt_t result=0;
1210   
1211   TList* list = GetFileSources(sourceFXS,"CE");
1212   
1213   if (list && list->GetEntries()>0) {
1214     
1215 //  loop through all files from LDCs
1216     
1217     UInt_t index = 0;
1218     while (list->At(index)!=NULL) {
1219       TObjString* fileNameEntry = (TObjString*) list->At(index);
1220       if (fileNameEntry!=NULL) {
1221         TString fileName = GetFile(sourceFXS, "CE",
1222                                    fileNameEntry->GetString().Data());
1223         AliTPCCalibCE *calCE=AliTPCCalibCE::ReadFromFile(fileName.Data());        
1224         
1225         if (!calCE) {
1226           Log ("No valid calibCE object.");
1227           result=2;
1228           break;
1229         }
1230         //  replace entries for the sectors available in the present file
1231         
1232         for (Int_t sector=0; sector<nSectors; sector++) {
1233           AliTPCCalROC *rocTmean=calCE->GetCalRocT0(sector);
1234           if ( rocTmean )  ceTmean->SetCalROC(rocTmean,sector);
1235           AliTPCCalROC *rocTrms=calCE->GetCalRocRMS(sector);
1236           if ( rocTrms )  ceTrms->SetCalROC(rocTrms,sector);
1237           AliTPCCalROC *rocQmean=calCE->GetCalRocQ(sector);
1238           if ( rocQmean )  ceQmean->SetCalROC(rocQmean,sector);
1239           TGraph *grT=calCE->MakeGraphTimeCE(sector,0,2); // T time graph
1240           if ( grT ) rocTtime->AddAt(grT,sector);
1241           TGraph *grQ=calCE->MakeGraphTimeCE(sector,0,3); // Q time graph
1242           if ( grQ ) rocQtime->AddAt(grQ,sector);
1243         }
1244         
1245         TGraph *grT=calCE->MakeGraphTimeCE(-1,0,2); // A side average
1246         if ( grT ) {
1247           rocTtime->AddAt(grT,nSectors);
1248         } else {
1249           result=10;
1250         }
1251         grT=calCE->MakeGraphTimeCE(-2,0,2); // C side average
1252         if ( grT ) {
1253           rocTtime->AddAt(grT,nSectors+1);
1254         } else {
1255           result=10;
1256         }
1257
1258         delete calCE;
1259       }
1260       ++index;
1261     }  // while(list)
1262 //
1263 //   Check number of calibrated sectors per side
1264 // 
1265     Int_t aside=0, cside=0;
1266     for (Int_t ind=0; ind<nSectors/4; ind++ ) {
1267         TGraph *grT=(TGraph*)rocTtime->At(ind);
1268         if (grT) aside++;
1269         grT=(TGraph*)rocTtime->At(ind+nSectors/2);
1270         if (grT) aside++;
1271         grT=(TGraph*)rocTtime->At(ind+nSectors/4);
1272         if (grT) cside++;
1273         grT=(TGraph*)rocTtime->At(ind+3*nSectors/4);
1274         if (grT) cside++;
1275      }
1276      if ( (aside<kMinCESectors) && (cside<kMinCESectors) ) {
1277         Log (Form("ExtractCE: Too few fitted sectors: Aside =%d, Cside=%d\n",
1278              aside, cside)) ;
1279         result=10;
1280      }
1281
1282     //
1283     //=== New CE part
1284     //    if it is validated this part needs to be modified again
1285     //    currently its only processed if there is a valid standard CE object
1286     //
1287     list = GetFileSources(sourceFXS,"CEnew");
1288     
1289     if (result==0 && list && list->GetEntries()>0) {
1290       
1291 //  loop through all files from LDCs
1292       
1293       UInt_t index2 = 0;
1294       while (list->At(index2)!=NULL) {
1295         TObjString* fileNameEntry = (TObjString*) list->At(index2);
1296         if (fileNameEntry!=NULL) {
1297           TString fileName = GetFile(sourceFXS, "CEnew",
1298                                      fileNameEntry->GetString().Data());
1299           AliTPCCalibCE *calCE=AliTPCCalibCE::ReadFromFile(fileName.Data());        
1300           
1301           if (!calCE) {
1302             Log ("No valid new calibCE object.");
1303 //             result=2;
1304             break;
1305           }
1306
1307           TIter nextObj(calCE->GetArrFitGraphs());
1308           TObject *obj=0x0;
1309           while ( (obj=nextObj()) ){
1310             arrFitGraphs->Add(obj->Clone());
1311           }
1312           delete calCE;
1313         }
1314         ++index2;
1315       }
1316     }
1317     
1318 //
1319 //  Store updated pedestal entry to OCDB
1320 //
1321     AliCDBMetaData metaData;
1322     metaData.SetBeamPeriod(0);
1323     metaData.SetResponsible("Haavard Helstrup");
1324     metaData.SetAliRootVersion(ALIROOT_BRANCH);
1325     metaData.SetComment("Preprocessor AliTPC data base entries.");
1326     
1327     if ( result == 0 ) {
1328       Bool_t storeOK = Store("Calib", "CE", ceObjects, &metaData, 0, kTRUE);
1329       if ( !storeOK ) ++result;
1330     } else {
1331       Log ("Warning: Average time graphs not available - no OCDB entry written");
1332     }
1333   } else {
1334     Log ("Error: no CE entries available from FXS!");
1335     result = 1;
1336   }
1337   
1338   ceObjects->Delete();
1339   delete ceObjects;
1340   
1341   return result;
1342 }
1343 //______________________________________________________________________________________________
1344
1345 UInt_t AliTPCPreprocessor::ExtractQA(Int_t sourceFXS)
1346 {
1347  //
1348  //  Read Quality Assurance file from file exchage server
1349  //
1350  
1351  UInt_t result=0;
1352
1353  TList* list = GetFileSources(sourceFXS,"QA");
1354  
1355  if (list && list->GetEntries()>0) {
1356
1357 //  only one QA object should be available!
1358
1359     AliTPCdataQA *calQA;
1360
1361     UInt_t nentries = list->GetEntries();  
1362     UInt_t index=0;
1363     if ( nentries > 1) Log ( "More than one QA entry. First one processed");      
1364     TObjString* fileNameEntry = (TObjString*) list->At(index);
1365     if (fileNameEntry!=NULL) {
1366         TString fileName = GetFile(sourceFXS, "QA",
1367                                          fileNameEntry->GetString().Data());
1368         TFile *f = TFile::Open(fileName);
1369         if (!f) {
1370           Log ("Error opening QA file.");
1371           result =2;          
1372         } else {
1373           f->GetObject("tpcCalibQA",calQA);
1374           if ( calQA ) {      
1375 //
1376 //  Store updated pedestal entry to OCDB
1377 //
1378            AliCDBMetaData metaData;
1379            metaData.SetBeamPeriod(0);
1380            metaData.SetResponsible("Haavard Helstrup");
1381            metaData.SetAliRootVersion(ALIROOT_BRANCH);
1382            metaData.SetComment("Preprocessor AliTPC data base entries.");
1383
1384            Bool_t storeOK = Store("Calib", "QA", calQA, &metaData, 0, kFALSE);
1385            if ( !storeOK ) ++result;
1386
1387            delete calQA;
1388           }
1389         }
1390     } else {
1391     Log ("Error: no QA files on FXS!");
1392     result = 2;
1393     }
1394   } else {
1395     Log ("Error: no QA entries in FXS list!");
1396     result = 1;
1397   }
1398   return result;
1399 }
1400
1401 //______________________________________________________________________________________________
1402
1403
1404 UInt_t AliTPCPreprocessor::ExtractAltro(Int_t sourceFXS, TMap* dcsMap)
1405 {
1406  //
1407  //  Read Altro configuration file from file exchage server
1408  //  Keep original entry from OCDB in case no new pulser calibration is available
1409  //
1410  TObjArray    *altroObjects=0;
1411  AliTPCCalPad *acqStart=0;
1412  AliTPCCalPad *zsThr=0;
1413  AliTPCCalPad *acqStop=0;
1414  AliTPCCalPad *FPED=0;
1415  AliTPCCalPad *masked=0;
1416  AliTPCCalPad *k1=0, *k2=0, *k3=0;
1417  AliTPCCalPad *l1=0, *l2=0, *l3=0;
1418  TMap *mapRCUconfig=0;
1419
1420  AliCDBEntry* entry = GetFromOCDB("Calib", "AltroConfig");
1421  if (entry) altroObjects = (TObjArray*)entry->GetObject();
1422  if ( altroObjects==NULL ) {
1423      Log("AliTPCPreprocsessor: No previous TPC altro calibration entry available.\n");
1424      altroObjects = new TObjArray;    
1425  }
1426
1427  acqStart = (AliTPCCalPad*)altroObjects->FindObject("AcqStart");
1428  if ( !acqStart ) {
1429     acqStart = new AliTPCCalPad("AcqStart","AcqStart");
1430     altroObjects->Add(acqStart);
1431  }
1432  zsThr = (AliTPCCalPad*)altroObjects->FindObject("ZsThr");
1433  if ( !zsThr )  { 
1434     zsThr = new AliTPCCalPad("ZsThr","ZsThr");
1435     altroObjects->Add(zsThr);
1436  }
1437  FPED = (AliTPCCalPad*)altroObjects->FindObject("FPED");
1438  if ( !FPED )  { 
1439     FPED = new AliTPCCalPad("FPED","FPED");
1440     altroObjects->Add(FPED);
1441  }
1442  acqStop = (AliTPCCalPad*)altroObjects->FindObject("AcqStop");
1443  if ( !acqStop ) {
1444     acqStop = new AliTPCCalPad("AcqStop","AcqStop");
1445     altroObjects->Add(acqStop);
1446  }
1447  masked = (AliTPCCalPad*)altroObjects->FindObject("Masked");
1448  if ( !masked )  { 
1449     masked = new AliTPCCalPad("Masked","Masked");
1450     altroObjects->Add(masked);
1451  }
1452  k1 = (AliTPCCalPad*)altroObjects->FindObject("K1");
1453  if ( !k1 )  { 
1454     k1 = new AliTPCCalPad("K1","K1");
1455     altroObjects->Add(k1);
1456  }
1457  k2 = (AliTPCCalPad*)altroObjects->FindObject("K2");
1458  if ( !k2 )  { 
1459     k2 = new AliTPCCalPad("K2","K2");
1460     altroObjects->Add(k2);
1461  }
1462  k3 = (AliTPCCalPad*)altroObjects->FindObject("K3");
1463  if ( !k3 )  { 
1464     k3 = new AliTPCCalPad("K3","K3");
1465     altroObjects->Add(k3);
1466  }
1467  l1 = (AliTPCCalPad*)altroObjects->FindObject("L1");
1468  if ( !l1 )  { 
1469     l1 = new AliTPCCalPad("L1","L1");
1470     altroObjects->Add(l1);
1471  }
1472  l2 = (AliTPCCalPad*)altroObjects->FindObject("L2");
1473  if ( !l2 )  { 
1474     l2 = new AliTPCCalPad("L2","L2");
1475     altroObjects->Add(l2);
1476  }
1477  l3 = (AliTPCCalPad*)altroObjects->FindObject("L3");
1478  if ( !l3 )  { 
1479     l3 = new AliTPCCalPad("L3","L3");
1480     altroObjects->Add(l3);
1481  }
1482  mapRCUconfig = (TMap*)altroObjects->FindObject("RCUconfig");
1483  if (!mapRCUconfig) {
1484     mapRCUconfig = new TMap();
1485     mapRCUconfig->SetName("RCUconfig");
1486     altroObjects->Add(mapRCUconfig);
1487  }
1488
1489
1490  UInt_t result=0;
1491  TString idFXS[2]={"AltroConfigA","AltroConfigC"};
1492
1493  Int_t nSectors = fROC->GetNSectors();
1494  Bool_t changed=false;
1495  if (altroObjects == 0 ) altroObjects = new TObjArray;
1496
1497 // extract list of active DDLs
1498
1499  if (dcsMap) {
1500   Bool_t found; 
1501   TString arrDDL(kNumDDL);
1502   arrDDL.Append('x',kNumDDL);
1503   for ( Int_t iDDL = 0; iDDL<kNumDDL; iDDL++ ) {
1504     TString stringID = Form (kAmandaDDL.Data(),iDDL+kFirstDDL);
1505     TPair *pair = (TPair*)dcsMap->FindObject(stringID.Data());
1506     found = false;
1507     if ( pair ) {
1508         TObjArray *valueSet=(TObjArray*)pair->Value();
1509         if ( valueSet) { 
1510           AliDCSValue *val = (AliDCSValue*)valueSet->At(0);
1511           if (val) { 
1512               found = val->GetBool();
1513               if (found){
1514                arrDDL[iDDL] = '1';
1515               } else { 
1516                arrDDL[iDDL] = '0';
1517               }
1518           }    
1519         }
1520     } 
1521   }
1522   TObjString *ddlArray = new TObjString;
1523   ddlArray->SetString(arrDDL);
1524   TMap *activeDDL = new TMap;
1525   activeDDL->SetName("DDLArray");
1526   TObjString *key = new TObjString("DDLArray");
1527   activeDDL->Add(key,ddlArray);
1528   altroObjects->Add(activeDDL);
1529   changed=true;
1530  } else {
1531    Log ("ExtractAltro: No DCS map available. Active DDL list cannot be obtained.");
1532    result = 3;
1533  }        
1534
1535 // extract Altro configuration files
1536
1537  for ( Int_t id=0; id<2; id++) {
1538    TList* list = GetFileSources(sourceFXS,idFXS[id].Data());
1539  
1540    if (list && list->GetEntries()>0) {
1541
1542 //  loop through all files from LDCs
1543
1544     UInt_t index = 0;
1545     while (list->At(index)!=NULL) {
1546      TObjString* fileNameEntry = (TObjString*) list->At(index);
1547      if (fileNameEntry!=NULL) {
1548         TString fileName = GetFile(sourceFXS, idFXS[id].Data(),
1549                                          fileNameEntry->GetString().Data());
1550         TFile *f = TFile::Open(fileName);
1551         if (!f) {
1552           const int mess_length=40;
1553           char message[mess_length];
1554           snprintf(message,mess_length,"Error opening Altro configuration file, id = %d",id);
1555           Log (message);
1556           result =2;
1557           break;
1558         }
1559         TObjArray *altroFXS;
1560         f->GetObject("AltroConfig",altroFXS);
1561         if ( !altroFXS ) {
1562           Log ("No Altro configuration object in file.");
1563           result = 2;
1564           break;
1565         }
1566
1567         //  replace entries for the sectors available in the present file
1568         AliTPCCalPad *acqStartFXS=(AliTPCCalPad*)altroFXS->FindObject("AcqStart");
1569         AliTPCCalPad *zsThrFXS=(AliTPCCalPad*)altroFXS->FindObject("ZsThr");
1570         AliTPCCalPad *acqStopFXS=(AliTPCCalPad*)altroFXS->FindObject("AcqStop");
1571         AliTPCCalPad *FPEDFXS=(AliTPCCalPad*)altroFXS->FindObject("FPED");
1572         AliTPCCalPad *maskedFXS=(AliTPCCalPad*)altroFXS->FindObject("Masked");
1573         AliTPCCalPad *k1FXS=(AliTPCCalPad*)altroFXS->FindObject("K1");
1574         AliTPCCalPad *k2FXS=(AliTPCCalPad*)altroFXS->FindObject("K2");
1575         AliTPCCalPad *k3FXS=(AliTPCCalPad*)altroFXS->FindObject("K3");
1576         AliTPCCalPad *l1FXS=(AliTPCCalPad*)altroFXS->FindObject("L1");
1577         AliTPCCalPad *l2FXS=(AliTPCCalPad*)altroFXS->FindObject("L2");
1578         AliTPCCalPad *l3FXS=(AliTPCCalPad*)altroFXS->FindObject("L3");
1579         TMap *mapRCUconfigFXS = (TMap*)altroFXS->FindObject("RCUconfig");
1580         TIterator *mapFXSiter = mapRCUconfigFXS->MakeIterator();
1581         
1582         changed=true;
1583         for (Int_t sector=0; sector<nSectors; sector++) {
1584             
1585            if (acqStartFXS) {
1586               AliTPCCalROC *rocAcqStart=acqStartFXS->GetCalROC(sector);
1587               if ( rocAcqStart )  acqStart->SetCalROC(rocAcqStart,sector);
1588            }
1589            if (zsThrFXS ) {
1590               AliTPCCalROC *rocZsThr=zsThrFXS->GetCalROC(sector);
1591               if ( rocZsThr )  zsThr->SetCalROC(rocZsThr,sector);
1592            }
1593            if (acqStopFXS) {
1594               AliTPCCalROC *rocAcqStop=acqStopFXS->GetCalROC(sector);
1595               if ( rocAcqStop )  acqStop->SetCalROC(rocAcqStop,sector);
1596            }
1597            if (FPEDFXS ) {
1598               AliTPCCalROC *rocFPED=FPEDFXS->GetCalROC(sector);
1599               if ( rocFPED )  FPED->SetCalROC(rocFPED,sector);
1600            }
1601            if (maskedFXS) {
1602               AliTPCCalROC *rocMasked=maskedFXS->GetCalROC(sector);
1603               if ( rocMasked )  masked->SetCalROC(rocMasked,sector);
1604            }
1605            if (k1FXS) {
1606               AliTPCCalROC *rocK1=k1FXS->GetCalROC(sector);
1607               if ( rocK1 )  k1->SetCalROC(rocK1,sector);
1608            }
1609            if (k2FXS) {
1610               AliTPCCalROC *rocK2=k2FXS->GetCalROC(sector);
1611               if ( rocK2 )  k2->SetCalROC(rocK2,sector);
1612            }
1613            if (k3FXS) {
1614               AliTPCCalROC *rocK3=k3FXS->GetCalROC(sector);
1615               if ( rocK3 )  k3->SetCalROC(rocK3,sector);
1616            }
1617            if (l1FXS) {
1618               AliTPCCalROC *rocL1=l1FXS->GetCalROC(sector);
1619               if ( rocL1 )  l1->SetCalROC(rocL1,sector);
1620            }
1621            if (l2FXS) {
1622               AliTPCCalROC *rocL2=l2FXS->GetCalROC(sector);
1623               if ( rocL2 )  l2->SetCalROC(rocL2,sector);
1624            }
1625            if (l3FXS) {
1626               AliTPCCalROC *rocL3=l3FXS->GetCalROC(sector);
1627               if ( rocL3 )  l3->SetCalROC(rocL3,sector);
1628            }
1629          }
1630          if (mapRCUconfigFXS) {
1631           Int_t mapEntries = mapRCUconfigFXS->GetEntries();
1632           TObjString* keyFXS;
1633           TVectorF* vecFXS;
1634           TVectorF* vec;              // nSectors = 72  (total number of inner/outer sectors)
1635           for (Int_t i=0; i<mapEntries; ++i) {
1636             keyFXS=(TObjString*)mapFXSiter->Next();
1637             vecFXS=(TVectorF*)mapRCUconfigFXS->GetValue(keyFXS);
1638             vec=(TVectorF*)mapRCUconfig->GetValue(keyFXS);
1639             if (!vec) {
1640               vec = new TVectorF(3*nSectors);
1641               *vec = -1;
1642               mapRCUconfig->Add(keyFXS,vec);
1643             }
1644             if (vec->GetNoElements() != 3*nSectors ) {
1645               vec->ResizeTo(3*nSectors);
1646             }
1647             if (id==0) {                        // A side
1648               vec->SetSub(0,vecFXS->GetSub(0,nSectors/2-1));
1649               vec->SetSub(nSectors,vecFXS->GetSub(nSectors,2*nSectors-1));
1650             } else {                             // C side
1651               vec->SetSub(nSectors/2,vecFXS->GetSub(nSectors/2,nSectors-1));
1652               vec->SetSub(2*nSectors,vecFXS->GetSub(2*nSectors,3*nSectors-1));
1653             }
1654           }
1655         }
1656        delete altroFXS;
1657        f->Close();
1658       }
1659      ++index;
1660      }  // while(list)
1661     } else {
1662       Log ("Error: no entries in AltroConfig file list!");
1663       result = 1;
1664     }
1665
1666    }   // for - id
1667 //
1668 //  Store updated pedestal entry to OCDB
1669 //
1670     if (changed) {
1671      AliCDBMetaData metaData;
1672      metaData.SetBeamPeriod(0);
1673      metaData.SetResponsible("Haavard Helstrup");
1674      metaData.SetAliRootVersion(ALIROOT_BRANCH);
1675      metaData.SetComment("Preprocessor AliTPC data base entries.");
1676
1677      Bool_t storeOK = Store("Calib", "AltroConfig", altroObjects, &metaData, 0, kFALSE);
1678      if ( !storeOK ) ++result;
1679     }  
1680
1681   altroObjects->Delete();
1682   delete altroObjects;
1683   
1684   return result;
1685 }