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