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