]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCPreprocessor.cxx
Add test for sector-by-sector fits from CE calibration. Fail CE entry production
[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) dcsResult=1;
253   if (dcsAliasMap->GetEntries() == 0 ) dcsResult=1;  
254   status = new TParameter<int>("dcsResult",dcsResult);
255   resultArray->Add(status);
256
257
258   TString runType = GetRunType();
259
260   if ( dcsResult == 0 ) {
261
262   // Temperature sensors are processed by AliTPCCalTemp
263
264     TString tempConf = fConfEnv->GetValue("Temperature","ON");
265     tempConf.ToUpper();
266     if (tempConf != "OFF" ) {
267       UInt_t tempResult = MapTemperature(dcsAliasMap);
268       if ( tempConf != "TRY") result+=tempResult;
269       status = new TParameter<int>("tempResult",tempResult);
270       resultArray->Add(status);
271     }
272
273   // High Voltage recordings
274
275
276     TString hvConf = fConfEnv->GetValue("HighVoltage","ON");
277     hvConf.ToUpper();
278     if (hvConf != "OFF" ) { 
279      UInt_t hvResult = MapHighVoltage(dcsAliasMap);
280      if (hvConf != "TRY") result+=hvResult;
281      status = new TParameter<int>("hvResult",hvResult);
282      resultArray->Add(status);
283     }
284
285     // Goofie values
286
287
288     TString goofieConf = fConfEnv->GetValue("Goofie","ON");
289     goofieConf.ToUpper();
290     if (goofieConf != "OFF" ) { 
291      UInt_t goofieResult = MapGoofie(dcsAliasMap);
292      if (goofieConf != "TRY") result+=goofieResult;
293      status = new TParameter<int>("goofieResult",goofieResult);
294      resultArray->Add(status);
295     }
296
297     // Pressure values
298
299     if( runType == kPhysicsRunType || 
300       runType == kLaserRunType ) {    
301
302       TString pressureConf = fConfEnv->GetValue("Pressure","ON");
303       pressureConf.ToUpper();
304       if (pressureConf != "OFF" ) { 
305        UInt_t pressureResult = MapPressure(dcsAliasMap);
306        status = new TParameter<int>("pressureResult",pressureResult);
307        resultArray->Add(status);
308       }
309     }
310   }
311   // Other calibration information will be retrieved through FXS files
312   //  examples:
313   //    TList* fileSourcesDAQ = GetFile(AliShuttleInterface::kDAQ, "pedestals");
314   //    const char* fileNamePed = GetFile(AliShuttleInterface::kDAQ, "pedestals", "LDC1");
315   //
316   //    TList* fileSourcesHLT = GetFile(AliShuttleInterface::kHLT, "calib");
317   //    const char* fileNameHLT = GetFile(AliShuttleInterface::kHLT, "calib", "LDC1");
318
319   // pedestal entries
320
321   if(runType == kPedestalRunType) {
322     Int_t numSources = 1;
323     Int_t pedestalSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
324     TString source = fConfEnv->GetValue("Pedestal","DAQ");
325     source.ToUpper();
326     if (source != "OFF" ) { 
327      if ( source == "HLT") pedestalSource[0] = AliShuttleInterface::kHLT;
328      if (!GetHLTStatus()) pedestalSource[0] = AliShuttleInterface::kDAQ;
329      if (source == "HLTDAQ" ) {
330          numSources=2;
331          pedestalSource[0] = AliShuttleInterface::kHLT;
332          pedestalSource[1] = AliShuttleInterface::kDAQ;
333      }
334      if (source == "DAQHLT" ) numSources=2;
335      UInt_t pedestalResult=0;
336      for (Int_t i=0; i<numSources; i++ ) {      
337        pedestalResult = ExtractPedestals(pedestalSource[i]);
338        if ( pedestalResult == 0 ) break;
339      }
340      result += pedestalResult;
341      status = new TParameter<int>("pedestalResult",pedestalResult);
342      resultArray->Add(status);
343     }
344   }
345
346   // pulser trigger processing
347
348   if(runType == kPulserRunType) {
349     Int_t numSources = 1;
350     Int_t pulserSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
351     TString source = fConfEnv->GetValue("Pulser","DAQ");
352     source.ToUpper();
353     if ( source != "OFF") { 
354      if ( source == "HLT") pulserSource[0] = AliShuttleInterface::kHLT;
355      if (!GetHLTStatus()) pulserSource[0] = AliShuttleInterface::kDAQ;
356      if (source == "HLTDAQ" ) {
357          numSources=2;
358          pulserSource[0] = AliShuttleInterface::kHLT;
359          pulserSource[1] = AliShuttleInterface::kDAQ;
360      }
361      if (source == "DAQHLT" ) numSources=2;
362      if (source == "TRY" ) numSources=2;
363      UInt_t pulserResult=0;
364      for (Int_t i=0; i<numSources; i++ ) {      
365        pulserResult = ExtractPulser(pulserSource[i]);
366        if ( pulserResult == 0 ) break;
367      }
368      if (source != "TRY") result += pulserResult;
369      status = new TParameter<int>("pulserResult",pulserResult);
370      resultArray->Add(status);
371     }
372   }
373
374
375 // raw calibration processing
376
377   if(runType == kPhysicsRunType) {
378     Int_t numSources = 1;
379     Int_t rawSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
380     TString source = fConfEnv->GetValue("Raw","DAQ");
381     source.ToUpper();
382     if ( source != "OFF") { 
383      if ( source == "HLT") rawSource[0] = AliShuttleInterface::kHLT;
384      if (!GetHLTStatus()) rawSource[0] = AliShuttleInterface::kDAQ;
385      if (source == "HLTDAQ" ) {
386          numSources=2;
387          rawSource[0] = AliShuttleInterface::kHLT;
388          rawSource[1] = AliShuttleInterface::kDAQ;
389      }
390      if (source == "DAQHLT" ) numSources=2;
391      if (source == "TRY" ) numSources=2;
392      UInt_t rawResult=0;
393      for (Int_t i=0; i<numSources; i++ ) {      
394        rawResult = ExtractRaw(rawSource[i]);
395        if ( rawResult == 0 ) break;
396      }
397      if (source != "TRY" )result += rawResult;
398      status = new TParameter<int>("rawResult",rawResult);
399      resultArray->Add(status);
400     }
401   }
402
403
404   // Altro configuration
405
406
407   TString altroConf = fConfEnv->GetValue("AltroConf","ON");
408   altroConf.ToUpper();
409   if (altroConf != "OFF" ) { 
410    UInt_t altroResult = ExtractAltro(AliShuttleInterface::kDCS,dcsAliasMap);
411    if (altroConf != "TRY" ) result+=altroResult;
412    status = new TParameter<int>("altroResult",altroResult);
413    resultArray->Add(status);
414  }
415
416
417   // Central Electrode processing
418
419   if( runType == kPhysicsRunType || 
420       runType == kLaserRunType ) {    
421
422     Int_t numSources = 1;
423     Int_t ceSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
424     TString source = fConfEnv->GetValue("CE","DAQ");
425     source.ToUpper();
426     if ( source != "OFF" ) { 
427      if ( source == "HLT") ceSource[0] = AliShuttleInterface::kHLT;
428      if (!GetHLTStatus()) ceSource[0] = AliShuttleInterface::kDAQ;
429      if (source == "HLTDAQ" ) {
430         numSources=2;
431         ceSource[0] = AliShuttleInterface::kHLT;
432         ceSource[1] = AliShuttleInterface::kDAQ;
433      }
434      if (source == "DAQHLT" ) numSources=2;
435      if (source == "TRY" ) numSources=2;
436      UInt_t ceResult=0;
437      for (Int_t i=0; i<numSources; i++ ) {      
438        ceResult = ExtractCE(ceSource[i]);
439        if ( ceResult == 0 ) break;
440      }
441
442    // only flag error if CE result is missing from LASER runs
443    //    -- for PHYSICS run do CE processing if data available
444    
445      if ( runType == kLaserRunType && source != "TRY"  && ceResult<10 ) result += ceResult;
446      status = new TParameter<int>("ceResult",ceResult);
447      resultArray->Add(status);
448
449     numSources = 1;
450     Int_t qaSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
451     source = fConfEnv->GetValue("QA","DAQ");
452     source.ToUpper();
453     if ( source != "OFF" ) { 
454      if ( source == "HLT") qaSource[0] = AliShuttleInterface::kHLT;
455      if (!GetHLTStatus()) qaSource[0] = AliShuttleInterface::kDAQ;
456      if (source == "HLTDAQ" ) {
457         numSources=2;
458         qaSource[0] = AliShuttleInterface::kHLT;
459         qaSource[1] = AliShuttleInterface::kDAQ;
460      }
461      if (source == "DAQHLT" ) numSources=2;
462      if (source == "TRY" ) numSources=2;
463      UInt_t qaResult=0;
464      for (Int_t i=0; i<numSources; i++ ) {      
465        qaResult = ExtractQA(qaSource[i]);
466        if ( qaResult == 0 ) break;
467      }
468 //     result += qaResult;
469      if ( qaResult !=0 ) Log ("ExtractQA failed, no QA entry available.");
470      status = new TParameter<int>("qaResult",qaResult);
471      resultArray->Add(status);
472     }
473    }
474   }
475   
476 // Store component status to OCDB
477
478   AliCDBMetaData metaData;
479   metaData.SetBeamPeriod(0);
480   metaData.SetResponsible("Haavard Helstrup");
481   metaData.SetAliRootVersion(ALIROOT_SVN_BRANCH);
482   metaData.SetComment("Preprocessor AliTPC status.");
483   Store("Calib", "PreprocStatus", resultArray, &metaData, 0, kFALSE);
484   resultArray->Delete();
485   delete resultArray;
486
487   if (errorHandling == "OFF" ) return 0;
488   return result;
489   
490 }
491 //______________________________________________________________________________________________
492 UInt_t AliTPCPreprocessor::MapTemperature(TMap* dcsAliasMap)
493 {
494
495    // extract DCS temperature maps. Perform fits to save space
496
497   UInt_t result=0;
498   TMap *map = fTemp->ExtractDCS(dcsAliasMap);
499   if (map) {
500     fTemp->MakeSplineFit(map);
501     Double_t fitFraction = 1.0*fTemp->NumFits()/fTemp->NumSensors(); 
502     if (fitFraction > kFitFraction ) {
503       AliInfo(Form("Temperature values extracted, fits performed.\n"));
504     } else { 
505       Log ("Too few temperature maps fitted. \n");
506       result = 9;
507     }
508   } else {
509     Log("No temperature map extracted. \n");
510     result=9;
511   }
512   delete map;
513   // Now store the final CDB file
514
515   if ( result == 0 ) {
516         AliCDBMetaData metaData;
517         metaData.SetBeamPeriod(0);
518         metaData.SetResponsible("Haavard Helstrup");
519         metaData.SetAliRootVersion(ALIROOT_SVN_BRANCH);
520         metaData.SetComment("Preprocessor AliTPC data base entries.");
521
522         Bool_t storeOK = Store("Calib", "Temperature", fTemp, &metaData, 0, kFALSE);
523         if ( !storeOK )  result=1;
524
525    }
526
527    return result;
528
529 }
530 //______________________________________________________________________________________________
531 UInt_t AliTPCPreprocessor::MapPressure(TMap* dcsAliasMap)
532 {
533
534    // extract DCS pressure maps. Perform fits to save space
535
536   UInt_t result=0;
537   TMap *map = fPressure->ExtractDCS(dcsAliasMap);
538   if (map) {
539     fPressure->MakeSplineFit(map);
540     Double_t fitFraction = 1.0*fPressure->NumFits()/fPressure->NumSensors(); 
541     if (fitFraction > kFitFraction ) {
542       AliInfo(Form("Pressure values extracted, fits performed.\n"));
543     } else { 
544       Log ("Too few pressure maps fitted. \n");
545       result = 9;
546     }
547   } else {
548     Log("No pressure map extracted. \n");
549     result=9;
550   }
551   delete map;
552   return result;
553 }
554
555 //______________________________________________________________________________________________
556 UInt_t AliTPCPreprocessor::MapHighVoltage(TMap* dcsAliasMap)
557 {
558
559    // extract DCS HV maps. Perform fits to save space
560
561   UInt_t result=0;
562   TMap *map = fHighVoltage->ExtractDCS(dcsAliasMap);
563   if (map) {
564     fHighVoltage->ClearFit();
565     fHighVoltage->RemoveGraphDuplicates(kHighVoltageDifference);
566                // don't keep new point if too similar to previous one
567     fHighVoltage->SetGraph(map);
568   } else {
569     Log("No high voltage recordings extracted. \n");
570     result=9;
571   }
572   delete map;
573
574   TString hvStatConf = fConfEnv->GetValue("HighVoltageStat","ON");
575   hvStatConf.ToUpper();
576   if (hvStatConf != "OFF" ) { 
577     TMap *map2 = fHighVoltageStat->ExtractDCS(dcsAliasMap);
578     if (map2) {
579       fHighVoltageStat->ClearFit();
580       fHighVoltageStat->SetGraph(map2);
581     } else {
582        Log("No high voltage status recordings extracted. \n");
583       result=9;
584     }
585     delete map2;
586
587     // add status maps to high voltage sensor array
588
589     fHighVoltage->AddSensors(fHighVoltageStat);
590    }
591   // Now store the final CDB file
592
593   if ( result == 0 ) {
594         AliCDBMetaData metaData;
595         metaData.SetBeamPeriod(0);
596         metaData.SetResponsible("Haavard Helstrup");
597         metaData.SetAliRootVersion(ALIROOT_SVN_BRANCH);
598         metaData.SetComment("Preprocessor AliTPC data base entries.");
599
600         Bool_t storeOK = Store("Calib", "HighVoltage", fHighVoltage, &metaData, 0, kFALSE);
601         if ( !storeOK )  result=1;
602
603    }
604
605    return result;
606
607 }
608
609 //______________________________________________________________________________________________
610 UInt_t AliTPCPreprocessor::MapGoofie(TMap* dcsAliasMap)
611 {
612
613    // extract DCS Goofie maps. Do not perform fits (low update rate)
614
615   UInt_t result=0;
616
617   TMap *map = fGoofie->ExtractDCS(dcsAliasMap);
618   if (map) {
619     fGoofie->ClearFit();
620     fGoofie->SetGraph(map);
621   } else {
622     Log("No Goofie recordings extracted. \n");
623     result=9;
624   }
625   delete map;
626
627   // Now store the final CDB file
628
629   if ( result == 0 ) {
630         AliCDBMetaData metaData;
631         metaData.SetBeamPeriod(0);
632         metaData.SetResponsible("Haavard Helstrup");
633         metaData.SetAliRootVersion(ALIROOT_SVN_BRANCH);
634         metaData.SetComment("Preprocessor AliTPC data base entries.");
635
636         Bool_t storeOK = Store("Calib", "Goofie", fGoofie, &metaData, 0, kFALSE);
637         if ( !storeOK )  result=1;
638
639    }
640
641    return result;
642
643 }
644
645
646 //______________________________________________________________________________________________
647
648 UInt_t AliTPCPreprocessor::ExtractPedestals(Int_t sourceFXS)
649 {
650  //
651  //  Read pedestal file from file exchage server
652  //  Keep original entry from OCDB in case no new pedestals are available
653  //
654  AliTPCCalPad *calPadPed=0;
655  AliCDBEntry* entry = GetFromOCDB("Calib", "Pedestals");
656  if (entry) calPadPed = (AliTPCCalPad*)entry->GetObject();
657  if ( calPadPed==NULL ) {
658      Log("AliTPCPreprocsessor: No previous TPC pedestal entry available.\n");
659      calPadPed = new AliTPCCalPad("PedestalsMean","PedestalsMean");
660  }
661
662  AliTPCCalPad *calPadRMS=0;
663  entry = GetFromOCDB("Calib", "PadNoise");
664  if (entry) calPadRMS = (AliTPCCalPad*)entry->GetObject();
665  if ( calPadRMS==NULL ) {
666      Log("AliTPCPreprocsessor: No previous TPC noise entry available.\n");
667      calPadRMS = new AliTPCCalPad("PedestalsRMS","PedestalsRMS");
668  }
669
670
671  UInt_t result=0;
672
673  Int_t nSectors = fROC->GetNSectors();
674  TList* list = GetFileSources(sourceFXS,"pedestals");
675  
676  if (list && list->GetEntries()>0) {
677
678 //  loop through all files from LDCs
679
680     Bool_t changed=false;
681     UInt_t index = 0;
682     while (list->At(index)!=NULL) {
683      TObjString* fileNameEntry = (TObjString*) list->At(index);
684      if (fileNameEntry!=NULL) {
685         TString fileName = GetFile(sourceFXS, "pedestals",
686                                          fileNameEntry->GetString().Data());
687         TFile *f = TFile::Open(fileName);
688         if (!f) {
689           Log ("Error opening pedestal file.");
690           result =2;
691           break;
692         }
693         AliTPCCalibPedestal *calPed;
694         f->GetObject("tpcCalibPedestal",calPed);
695         if ( !calPed ) {
696           Log ("No pedestal calibration object in file.");
697           result = 2;
698           break;
699         }
700
701         //  replace entries for the sectors available in the present file
702
703         changed=true;
704         for (Int_t sector=0; sector<nSectors; sector++) {
705            AliTPCCalROC *rocPed=calPed->GetCalRocPedestal(sector, kFALSE);
706            if ( rocPed )  calPadPed->SetCalROC(rocPed,sector);
707            AliTPCCalROC *rocRMS=calPed->GetCalRocRMS(sector, kFALSE);
708            if ( rocRMS )  calPadRMS->SetCalROC(rocRMS,sector);
709         }
710         delete calPed; 
711         f->Close();
712       }
713      ++index;
714     }  // while(list)
715 //
716 //  Store updated pedestal entry to OCDB
717 //
718     if (changed) {
719      AliCDBMetaData metaData;
720      metaData.SetBeamPeriod(0);
721      metaData.SetResponsible("Haavard Helstrup");
722      metaData.SetAliRootVersion(ALIROOT_SVN_BRANCH);
723      metaData.SetComment("Preprocessor AliTPC data base entries."); 
724  
725      Bool_t storeOK = Store("Calib", "Pedestals", calPadPed, &metaData, 0, kTRUE);
726      if ( !storeOK ) ++result;
727      storeOK = Store("Calib", "PadNoise", calPadRMS, &metaData, 0, kTRUE);
728      if ( !storeOK ) ++result;
729     }
730   } else {
731     Log ("Error: no entries in pedestal file list!");
732     result = 1;
733   }
734
735   delete calPadPed;
736   delete calPadRMS;
737
738   return result;
739 }
740
741 //______________________________________________________________________________________________
742
743
744 UInt_t AliTPCPreprocessor::ExtractPulser(Int_t sourceFXS)
745 {
746  //
747  //  Read pulser calibration file from file exchage server
748  //  Keep original entry from OCDB in case no new pulser calibration is available
749  //
750  TObjArray    *pulserObjects=0;
751  AliTPCCalPad *pulserTmean=0;
752  AliTPCCalPad *pulserTrms=0;
753  AliTPCCalPad *pulserQmean=0;
754  AliCDBEntry* entry = GetFromOCDB("Calib", "Pulser");
755  if (entry) pulserObjects = (TObjArray*)entry->GetObject();
756  if ( pulserObjects==NULL ) {
757      Log("AliTPCPreprocsessor: No previous TPC pulser entry available.\n");
758      pulserObjects = new TObjArray;    
759  }
760
761  pulserTmean = (AliTPCCalPad*)pulserObjects->FindObject("PulserTmean");
762  if ( !pulserTmean ) {
763     pulserTmean = new AliTPCCalPad("PulserTmean","PulserTmean");
764     pulserObjects->Add(pulserTmean);
765  }
766  pulserTrms = (AliTPCCalPad*)pulserObjects->FindObject("PulserTrms");
767  if ( !pulserTrms )  { 
768     pulserTrms = new AliTPCCalPad("PulserTrms","PulserTrms");
769     pulserObjects->Add(pulserTrms);
770  }
771  pulserQmean = (AliTPCCalPad*)pulserObjects->FindObject("PulserQmean");
772  if ( !pulserQmean )  { 
773     pulserQmean = new AliTPCCalPad("PulserQmean","PulserQmean");
774     pulserObjects->Add(pulserQmean);
775  }
776
777
778  UInt_t result=0;
779
780  Int_t nSectors = fROC->GetNSectors();
781  TList* list = GetFileSources(sourceFXS,"pulser");
782  
783  if (list && list->GetEntries()>0) {
784
785 //  loop through all files from LDCs
786
787     Bool_t changed=false;
788     UInt_t index = 0;
789     while (list->At(index)!=NULL) {
790      TObjString* fileNameEntry = (TObjString*) list->At(index);
791      if (fileNameEntry!=NULL) {
792         TString fileName = GetFile(sourceFXS, "pulser",
793                                          fileNameEntry->GetString().Data());
794         TFile *f = TFile::Open(fileName);
795         if (!f) {
796           Log ("Error opening pulser file.");
797           result =2;
798           break;
799         }
800         AliTPCCalibPulser *calPulser;
801         f->GetObject("tpcCalibPulser",calPulser);
802         if ( !calPulser ) {
803           Log ("No pulser calibration object in file.");
804           result = 2;
805           break;
806         }
807
808         //  replace entries for the sectors available in the present file
809
810         changed=true;
811         for (Int_t sector=0; sector<nSectors; sector++) {
812            AliTPCCalROC *rocTmean=calPulser->GetCalRocT0(sector);
813            if ( rocTmean )  pulserTmean->SetCalROC(rocTmean,sector);
814            AliTPCCalROC *rocTrms=calPulser->GetCalRocRMS(sector);
815            if ( rocTrms )  pulserTrms->SetCalROC(rocTrms,sector);
816            AliTPCCalROC *rocQmean=calPulser->GetCalRocQ(sector);
817            if ( rocQmean )  pulserQmean->SetCalROC(rocQmean,sector);
818         }
819        delete calPulser;
820        f->Close();
821       }
822      ++index;
823     }  // while(list)
824 //
825 //  Store updated pedestal entry to OCDB
826 //
827     if (changed) {
828      AliCDBMetaData metaData;
829      metaData.SetBeamPeriod(0);
830      metaData.SetResponsible("Haavard Helstrup");
831      metaData.SetAliRootVersion(ALIROOT_SVN_BRANCH);
832      metaData.SetComment("Preprocessor AliTPC data base entries.");
833
834      Bool_t storeOK = Store("Calib", "Pulser", pulserObjects, &metaData, 0, kTRUE);
835      if ( !storeOK ) ++result;
836     }  
837   } else {
838     Log ("Error: no entries in pulser file list!");
839     result = 1;
840   }
841   pulserObjects->Delete();
842   delete pulserObjects;
843
844   return result;
845 }
846
847 //______________________________________________________________________________________________
848
849
850 UInt_t AliTPCPreprocessor::ExtractRaw(Int_t sourceFXS)
851 {
852  //
853  //  Read Raw calibration file from file exchage server
854  //
855  
856  UInt_t result=0;
857  TObjArray* rawArray = new TObjArray;
858
859  TList* list = GetFileSources(sourceFXS,"tpcCalibRaw");
860  
861  if (list && list->GetEntries()>0) {
862
863 //  loop through all files
864
865     UInt_t index = 0;
866     while (list->At(index)!=NULL) {
867      TObjString* fileNameEntry = (TObjString*) list->At(index);
868      if (fileNameEntry!=NULL) {
869         TString fileName = GetFile(sourceFXS, "tpcCalibRaw",
870                                          fileNameEntry->GetString().Data());
871         TFile *f = TFile::Open(fileName);
872         if (!f) {
873           Log ("Error opening raw file.");
874           result =2;
875           break;
876         }
877         AliTPCCalibRaw *calRaw;
878         f->GetObject("tpcCalibRaw",calRaw);
879         if ( !calRaw ) {
880           Log ("No raw calibration object in file.");
881           result = 2;
882           break;
883         }
884        rawArray->Add(calRaw);
885        f->Close();
886       }
887      ++index;
888     }  // while(list)
889 //
890 //  Store updated pedestal entry to OCDB
891 //
892      AliCDBMetaData metaData;
893      metaData.SetBeamPeriod(0);
894      metaData.SetResponsible("Haavard Helstrup");
895      metaData.SetAliRootVersion(ALIROOT_SVN_BRANCH);
896      metaData.SetComment("Preprocessor AliTPC data base entries.");
897
898      Bool_t storeOK = Store("Calib", "Raw", rawArray, &metaData, 0, kTRUE);
899      if ( !storeOK ) ++result;
900   } else {
901     Log ("Error: no entries in raw file list!");
902     result = 1;
903   }
904   
905   rawArray->Delete();
906   delete rawArray;
907
908   return result;
909 }
910 //______________________________________________________________________________________________
911
912 UInt_t AliTPCPreprocessor::ExtractCE(Int_t sourceFXS)
913 {
914  //
915  //  Read Central Electrode file from file exchage server
916  //
917  //
918   AliTPCCalPad *ceTmean=0;
919   AliTPCCalPad *ceTrms=0;
920   AliTPCCalPad *ceQmean=0;
921   TObjArray    *rocTtime=0;
922   TObjArray    *rocQtime=0;
923   
924   TObjArray    *ceObjects= new TObjArray;
925   
926   
927   Int_t nSectors = fROC->GetNSectors();
928   
929   ceTmean = new AliTPCCalPad("CETmean","CETmean");
930   ceObjects->Add(ceTmean);
931   
932   ceTrms = new AliTPCCalPad("CETrms","CETrms");
933   ceObjects->Add(ceTrms);
934   
935   ceQmean = new AliTPCCalPad("CEQmean","CEQmean");
936   ceObjects->Add(ceQmean);
937   
938   rocTtime = new TObjArray(nSectors+2);   // also make room for A and C side average
939   rocTtime->SetName("rocTtime");
940   ceObjects->Add(rocTtime);
941   
942   rocQtime = new TObjArray(nSectors);
943   rocQtime->SetName("rocQtime");
944   ceObjects->Add(rocQtime);
945
946   //=== new part
947   TObjArray *arrFitGraphs=new TObjArray;
948   arrFitGraphs->SetName("ceFitsDrift");
949   ceObjects->Add(arrFitGraphs);
950   
951 // Temperature maps 
952   
953   if (fTemp) {
954     AliTPCSensorTempArray *tempMap = new AliTPCSensorTempArray(*fTemp);
955     tempMap->SetNameTitle("TempMap","TempMap");
956     ceObjects->Add(tempMap);
957   }
958   
959 // Pressure maps
960   
961   if (fPressure) {
962     AliDCSSensor *sensor=0, *sensorCopy=0;
963     for (Int_t isensor=0; isensor<kNumPressureSensors; ++isensor ) {
964       sensor = fPressure->GetSensor(kPressureSensorNames[isensor]);
965       if (sensor) {
966         sensorCopy = new AliDCSSensor(*sensor);
967         sensorCopy->SetNameTitle(kPressureSensorNames[isensor],kPressureSensorNames[isensor]);
968         ceObjects->Add(sensorCopy);
969       }
970     }
971   }
972   
973   UInt_t result=0;
974   
975   TList* list = GetFileSources(sourceFXS,"CE");
976   
977   if (list && list->GetEntries()>0) {
978     
979 //  loop through all files from LDCs
980     
981     UInt_t index = 0;
982     while (list->At(index)!=NULL) {
983       TObjString* fileNameEntry = (TObjString*) list->At(index);
984       if (fileNameEntry!=NULL) {
985         TString fileName = GetFile(sourceFXS, "CE",
986                                    fileNameEntry->GetString().Data());
987         TFile *f = TFile::Open(fileName);
988         if (!f) {
989           Log ("Error opening central electrode file.");
990           result =2;
991           break;
992         }
993         AliTPCCalibCE *calCE;
994         f->GetObject("tpcCalibCE",calCE);
995         
996         if (!calCE) {
997           Log ("No valid calibCE object.");
998           result=2;
999           break;
1000         }
1001         //  replace entries for the sectors available in the present file
1002         
1003         for (Int_t sector=0; sector<nSectors; sector++) {
1004           AliTPCCalROC *rocTmean=calCE->GetCalRocT0(sector);
1005           if ( rocTmean )  ceTmean->SetCalROC(rocTmean,sector);
1006           AliTPCCalROC *rocTrms=calCE->GetCalRocRMS(sector);
1007           if ( rocTrms )  ceTrms->SetCalROC(rocTrms,sector);
1008           AliTPCCalROC *rocQmean=calCE->GetCalRocQ(sector);
1009           if ( rocQmean )  ceQmean->SetCalROC(rocQmean,sector);
1010           TGraph *grT=calCE->MakeGraphTimeCE(sector,0,2); // T time graph
1011           if ( grT ) rocTtime->AddAt(grT,sector);
1012           TGraph *grQ=calCE->MakeGraphTimeCE(sector,0,3); // Q time graph
1013           if ( grQ ) rocQtime->AddAt(grQ,sector);
1014         }
1015         
1016         TGraph *grT=calCE->MakeGraphTimeCE(-1,0,2); // A side average
1017         if ( grT ) {
1018           rocTtime->AddAt(grT,nSectors);
1019         } else {
1020           result=10;
1021         }
1022         grT=calCE->MakeGraphTimeCE(-2,0,2); // C side average
1023         if ( grT ) {
1024           rocTtime->AddAt(grT,nSectors+1);
1025         } else {
1026           result=10;
1027         }
1028
1029         delete calCE;
1030         f->Close();
1031       }
1032       ++index;
1033     }  // while(list)
1034 //
1035 //   Check number of calibrated sectors per side
1036 // 
1037     Int_t aside=0, cside=0;
1038     for (Int_t ind=0; ind<nSectors/4; ind++ ) {
1039         TGraph *grT=(TGraph*)rocTtime->At(ind);
1040         if (grT) aside++;
1041         grT=(TGraph*)rocTtime->At(ind+nSectors/2);
1042         if (grT) aside++;
1043         grT=(TGraph*)rocTtime->At(ind+nSectors/4);
1044         if (grT) cside++;
1045         grT=(TGraph*)rocTtime->At(ind+3*nSectors/4);
1046         if (grT) cside++;
1047      }
1048      if ( (aside<kMinCESectors) && (cside<kMinCESectors) ) result=10;
1049
1050     //
1051     //=== New CE part
1052     //    if it is validated this part needs to be modified again
1053     //    currently its only processed if there is a valid standard CE object
1054     //
1055     list = GetFileSources(sourceFXS,"CEnew");
1056     
1057     if (result==0 && list && list->GetEntries()>0) {
1058       
1059 //  loop through all files from LDCs
1060       
1061       UInt_t index2 = 0;
1062       while (list->At(index2)!=NULL) {
1063         TObjString* fileNameEntry = (TObjString*) list->At(index2);
1064         if (fileNameEntry!=NULL) {
1065           TString fileName = GetFile(sourceFXS, "CEnew",
1066                                      fileNameEntry->GetString().Data());
1067           TFile *f = TFile::Open(fileName);
1068           if (!f) {
1069             Log ("Error opening new central electrode file.");
1070 //             result =2;
1071             break;
1072           }
1073           AliTPCCalibCE *calCE;
1074           f->GetObject("tpcCalibCE",calCE);
1075           
1076           if (!calCE) {
1077             Log ("No valid new calibCE object.");
1078 //             result=2;
1079             break;
1080           }
1081
1082           TIter nextObj(calCE->GetArrFitGraphs());
1083           TObject *obj=0x0;
1084           while ( (obj=nextObj()) ){
1085             arrFitGraphs->Add(obj->Clone());
1086           }
1087           delete calCE;
1088           f->Close();          
1089         }
1090         ++index2;
1091       }
1092     }
1093     
1094 //
1095 //  Store updated pedestal entry to OCDB
1096 //
1097     AliCDBMetaData metaData;
1098     metaData.SetBeamPeriod(0);
1099     metaData.SetResponsible("Haavard Helstrup");
1100     metaData.SetAliRootVersion(ALIROOT_SVN_BRANCH);
1101     metaData.SetComment("Preprocessor AliTPC data base entries.");
1102     
1103     if ( result == 0 ) {
1104       Bool_t storeOK = Store("Calib", "CE", ceObjects, &metaData, 0, kTRUE);
1105       if ( !storeOK ) ++result;
1106     } else {
1107       Log ("Warning: Average time graphs not available - no OCDB entry written");
1108     }
1109   } else {
1110     Log ("Error: no CE entries available from FXS!");
1111     result = 1;
1112   }
1113   
1114   ceObjects->Delete();
1115   delete ceObjects;
1116   
1117   return result;
1118 }
1119 //______________________________________________________________________________________________
1120
1121 UInt_t AliTPCPreprocessor::ExtractQA(Int_t sourceFXS)
1122 {
1123  //
1124  //  Read Quality Assurance file from file exchage server
1125  //
1126  
1127  UInt_t result=0;
1128
1129  TList* list = GetFileSources(sourceFXS,"QA");
1130  
1131  if (list && list->GetEntries()>0) {
1132
1133 //  only one QA objetc should be available!
1134
1135     AliTPCdataQA *calQA;
1136
1137     UInt_t nentries = list->GetEntries();  
1138     UInt_t index=0;
1139     if ( nentries > 1) Log ( "More than one QA entry. First one processed");      
1140     TObjString* fileNameEntry = (TObjString*) list->At(index);
1141     if (fileNameEntry!=NULL) {
1142         TString fileName = GetFile(sourceFXS, "QA",
1143                                          fileNameEntry->GetString().Data());
1144         TFile *f = TFile::Open(fileName);
1145         if (!f) {
1146           Log ("Error opening QA file.");
1147           result =2;          
1148         } else {
1149           f->GetObject("tpcCalibQA",calQA);
1150           if ( calQA ) {      
1151 //
1152 //  Store updated pedestal entry to OCDB
1153 //
1154            AliCDBMetaData metaData;
1155            metaData.SetBeamPeriod(0);
1156            metaData.SetResponsible("Haavard Helstrup");
1157            metaData.SetAliRootVersion(ALIROOT_SVN_BRANCH);
1158            metaData.SetComment("Preprocessor AliTPC data base entries.");
1159
1160            Bool_t storeOK = Store("Calib", "QA", calQA, &metaData, 0, kFALSE);
1161            if ( !storeOK ) ++result;
1162
1163            delete calQA;
1164           }
1165         }
1166     } else {
1167     Log ("Error: no QA files on FXS!");
1168     result = 2;
1169     }
1170   } else {
1171     Log ("Error: no QA entries in FXS list!");
1172     result = 1;
1173   }
1174   return result;
1175 }
1176
1177 //______________________________________________________________________________________________
1178
1179
1180 UInt_t AliTPCPreprocessor::ExtractAltro(Int_t sourceFXS, TMap* dcsMap)
1181 {
1182  //
1183  //  Read Altro configuration file from file exchage server
1184  //  Keep original entry from OCDB in case no new pulser calibration is available
1185  //
1186  TObjArray    *altroObjects=0;
1187  AliTPCCalPad *acqStart=0;
1188  AliTPCCalPad *zsThr=0;
1189  AliTPCCalPad *acqStop=0;
1190  AliTPCCalPad *FPED=0;
1191  AliTPCCalPad *masked=0;
1192  AliTPCCalPad *k1=0, *k2=0, *k3=0;
1193  AliTPCCalPad *l1=0, *l2=0, *l3=0;
1194  TMap *mapRCUconfig=0;
1195
1196  AliCDBEntry* entry = GetFromOCDB("Calib", "AltroConfig");
1197  if (entry) altroObjects = (TObjArray*)entry->GetObject();
1198  if ( altroObjects==NULL ) {
1199      Log("AliTPCPreprocsessor: No previous TPC altro calibration entry available.\n");
1200      altroObjects = new TObjArray;    
1201  }
1202
1203  acqStart = (AliTPCCalPad*)altroObjects->FindObject("AcqStart");
1204  if ( !acqStart ) {
1205     acqStart = new AliTPCCalPad("AcqStart","AcqStart");
1206     altroObjects->Add(acqStart);
1207  }
1208  zsThr = (AliTPCCalPad*)altroObjects->FindObject("ZsThr");
1209  if ( !zsThr )  { 
1210     zsThr = new AliTPCCalPad("ZsThr","ZsThr");
1211     altroObjects->Add(zsThr);
1212  }
1213  FPED = (AliTPCCalPad*)altroObjects->FindObject("FPED");
1214  if ( !FPED )  { 
1215     FPED = new AliTPCCalPad("FPED","FPED");
1216     altroObjects->Add(FPED);
1217  }
1218  acqStop = (AliTPCCalPad*)altroObjects->FindObject("AcqStop");
1219  if ( !acqStop ) {
1220     acqStop = new AliTPCCalPad("AcqStop","AcqStop");
1221     altroObjects->Add(acqStop);
1222  }
1223  masked = (AliTPCCalPad*)altroObjects->FindObject("Masked");
1224  if ( !masked )  { 
1225     masked = new AliTPCCalPad("Masked","Masked");
1226     altroObjects->Add(masked);
1227  }
1228  k1 = (AliTPCCalPad*)altroObjects->FindObject("K1");
1229  if ( !k1 )  { 
1230     k1 = new AliTPCCalPad("K1","K1");
1231     altroObjects->Add(k1);
1232  }
1233  k2 = (AliTPCCalPad*)altroObjects->FindObject("K2");
1234  if ( !k2 )  { 
1235     k2 = new AliTPCCalPad("K2","K2");
1236     altroObjects->Add(k2);
1237  }
1238  k3 = (AliTPCCalPad*)altroObjects->FindObject("K3");
1239  if ( !k3 )  { 
1240     k3 = new AliTPCCalPad("K3","K3");
1241     altroObjects->Add(k3);
1242  }
1243  l1 = (AliTPCCalPad*)altroObjects->FindObject("L1");
1244  if ( !l1 )  { 
1245     l1 = new AliTPCCalPad("L1","L1");
1246     altroObjects->Add(l1);
1247  }
1248  l2 = (AliTPCCalPad*)altroObjects->FindObject("L2");
1249  if ( !l2 )  { 
1250     l2 = new AliTPCCalPad("L2","L2");
1251     altroObjects->Add(l2);
1252  }
1253  l3 = (AliTPCCalPad*)altroObjects->FindObject("L3");
1254  if ( !l3 )  { 
1255     l3 = new AliTPCCalPad("L3","L3");
1256     altroObjects->Add(l3);
1257  }
1258  mapRCUconfig = (TMap*)altroObjects->FindObject("RCUconfig");
1259  if (!mapRCUconfig) {
1260     mapRCUconfig = new TMap();
1261     mapRCUconfig->SetName("RCUconfig");
1262     altroObjects->Add(mapRCUconfig);
1263  }
1264
1265
1266  UInt_t result=0;
1267  TString idFXS[2]={"AltroConfigA","AltroConfigC"};
1268
1269  Int_t nSectors = fROC->GetNSectors();
1270  Bool_t changed=false;
1271  if (altroObjects == 0 ) altroObjects = new TObjArray;
1272
1273 // extract list of active DDLs
1274
1275   Bool_t found; 
1276   TString arrDDL(kNumDDL);
1277   arrDDL.Append('x',kNumDDL);
1278   for ( Int_t iDDL = 0; iDDL<kNumDDL; iDDL++ ) {
1279     TString stringID = Form (kAmandaDDL.Data(),iDDL+kFirstDDL);
1280     TPair *pair = (TPair*)dcsMap->FindObject(stringID.Data());
1281     found = false;
1282     if ( pair ) {
1283         TObjArray *valueSet=(TObjArray*)pair->Value();
1284         if ( valueSet) { 
1285           AliDCSValue *val = (AliDCSValue*)valueSet->At(0);
1286           if (val) { 
1287               found = val->GetBool();
1288               if (found){
1289                arrDDL[iDDL] = '1';
1290               } else { 
1291                arrDDL[iDDL] = '0';
1292               }
1293           }    
1294         }
1295     } 
1296   }
1297   TObjString *ddlArray = new TObjString;
1298   ddlArray->SetString(arrDDL);
1299   TMap *activeDDL = new TMap;
1300   activeDDL->SetName("DDLArray");
1301   TObjString *key = new TObjString("DDLArray");
1302   activeDDL->Add(key,ddlArray);
1303   altroObjects->Add(activeDDL);
1304   changed=true;
1305   
1306
1307 // extract Altro configuration files
1308
1309  for ( Int_t id=0; id<2; id++) {
1310    TList* list = GetFileSources(sourceFXS,idFXS[id].Data());
1311  
1312    if (list && list->GetEntries()>0) {
1313
1314 //  loop through all files from LDCs
1315
1316     UInt_t index = 0;
1317     while (list->At(index)!=NULL) {
1318      TObjString* fileNameEntry = (TObjString*) list->At(index);
1319      if (fileNameEntry!=NULL) {
1320         TString fileName = GetFile(sourceFXS, idFXS[id].Data(),
1321                                          fileNameEntry->GetString().Data());
1322         TFile *f = TFile::Open(fileName);
1323         if (!f) {
1324           char message[40];
1325           sprintf(message,"Error opening Altro configuration file, id = %d",id);
1326           Log (message);
1327           result =2;
1328           break;
1329         }
1330         TObjArray *altroFXS;
1331         f->GetObject("AltroConfig",altroFXS);
1332         if ( !altroFXS ) {
1333           Log ("No Altro configuration object in file.");
1334           result = 2;
1335           break;
1336         }
1337
1338         //  replace entries for the sectors available in the present file
1339         AliTPCCalPad *acqStartFXS=(AliTPCCalPad*)altroFXS->FindObject("AcqStart");
1340         AliTPCCalPad *zsThrFXS=(AliTPCCalPad*)altroFXS->FindObject("ZsThr");
1341         AliTPCCalPad *acqStopFXS=(AliTPCCalPad*)altroFXS->FindObject("AcqStop");
1342         AliTPCCalPad *FPEDFXS=(AliTPCCalPad*)altroFXS->FindObject("FPED");
1343         AliTPCCalPad *maskedFXS=(AliTPCCalPad*)altroFXS->FindObject("Masked");
1344         AliTPCCalPad *k1FXS=(AliTPCCalPad*)altroFXS->FindObject("K1");
1345         AliTPCCalPad *k2FXS=(AliTPCCalPad*)altroFXS->FindObject("K2");
1346         AliTPCCalPad *k3FXS=(AliTPCCalPad*)altroFXS->FindObject("K3");
1347         AliTPCCalPad *l1FXS=(AliTPCCalPad*)altroFXS->FindObject("L1");
1348         AliTPCCalPad *l2FXS=(AliTPCCalPad*)altroFXS->FindObject("L2");
1349         AliTPCCalPad *l3FXS=(AliTPCCalPad*)altroFXS->FindObject("L3");
1350         TMap *mapRCUconfigFXS = (TMap*)altroFXS->FindObject("RCUconfig");
1351         TIterator *mapFXSiter = mapRCUconfigFXS->MakeIterator();
1352         
1353         changed=true;
1354         for (Int_t sector=0; sector<nSectors; sector++) {
1355             
1356            if (acqStartFXS) {
1357               AliTPCCalROC *rocAcqStart=acqStartFXS->GetCalROC(sector);
1358               if ( rocAcqStart )  acqStart->SetCalROC(rocAcqStart,sector);
1359            }
1360            if (zsThrFXS ) {
1361               AliTPCCalROC *rocZsThr=zsThrFXS->GetCalROC(sector);
1362               if ( rocZsThr )  zsThr->SetCalROC(rocZsThr,sector);
1363            }
1364            if (acqStopFXS) {
1365               AliTPCCalROC *rocAcqStop=acqStopFXS->GetCalROC(sector);
1366               if ( rocAcqStop )  acqStop->SetCalROC(rocAcqStop,sector);
1367            }
1368            if (FPEDFXS ) {
1369               AliTPCCalROC *rocFPED=FPEDFXS->GetCalROC(sector);
1370               if ( rocFPED )  FPED->SetCalROC(rocFPED,sector);
1371            }
1372            if (maskedFXS) {
1373               AliTPCCalROC *rocMasked=maskedFXS->GetCalROC(sector);
1374               if ( rocMasked )  masked->SetCalROC(rocMasked,sector);
1375            }
1376            if (k1FXS) {
1377               AliTPCCalROC *rocK1=k1FXS->GetCalROC(sector);
1378               if ( rocK1 )  k1->SetCalROC(rocK1,sector);
1379            }
1380            if (k2FXS) {
1381               AliTPCCalROC *rocK2=k2FXS->GetCalROC(sector);
1382               if ( rocK2 )  k2->SetCalROC(rocK2,sector);
1383            }
1384            if (k3FXS) {
1385               AliTPCCalROC *rocK3=k3FXS->GetCalROC(sector);
1386               if ( rocK3 )  k3->SetCalROC(rocK3,sector);
1387            }
1388            if (l1FXS) {
1389               AliTPCCalROC *rocL1=l1FXS->GetCalROC(sector);
1390               if ( rocL1 )  l1->SetCalROC(rocL1,sector);
1391            }
1392            if (l2FXS) {
1393               AliTPCCalROC *rocL2=l2FXS->GetCalROC(sector);
1394               if ( rocL2 )  l2->SetCalROC(rocL2,sector);
1395            }
1396            if (l3FXS) {
1397               AliTPCCalROC *rocL3=l3FXS->GetCalROC(sector);
1398               if ( rocL3 )  l3->SetCalROC(rocL3,sector);
1399            }
1400          }
1401          if (mapRCUconfigFXS) {
1402           Int_t mapEntries = mapRCUconfigFXS->GetEntries();
1403           TObjString* keyFXS;
1404           TVectorF* vecFXS;
1405           TVectorF* vec;              // nSectors = 72  (total number of inner/outer sectors)
1406           for (Int_t i=0; i<mapEntries; ++i) {
1407             keyFXS=(TObjString*)mapFXSiter->Next();
1408             vecFXS=(TVectorF*)mapRCUconfigFXS->GetValue(keyFXS);
1409             vec=(TVectorF*)mapRCUconfig->GetValue(keyFXS);
1410             if (!vec) {
1411               vec = new TVectorF(3*nSectors);
1412               *vec = -1;
1413               mapRCUconfig->Add(keyFXS,vec);
1414             }
1415             if (vec->GetNoElements() != 3*nSectors ) {
1416               vec->ResizeTo(3*nSectors);
1417             }
1418             if (id==0) {                        // A side
1419               vec->SetSub(0,vecFXS->GetSub(0,nSectors/2-1));
1420               vec->SetSub(nSectors,vecFXS->GetSub(nSectors,2*nSectors-1));
1421             } else {                             // C side
1422               vec->SetSub(nSectors/2,vecFXS->GetSub(nSectors/2,nSectors-1));
1423               vec->SetSub(2*nSectors,vecFXS->GetSub(2*nSectors,3*nSectors-1));
1424             }
1425           }
1426         }
1427        delete altroFXS;
1428        f->Close();
1429       }
1430      ++index;
1431      }  // while(list)
1432     } else {
1433       Log ("Error: no entries in AltroConfig file list!");
1434       result = 1;
1435     }
1436
1437    }   // for - id
1438 //
1439 //  Store updated pedestal entry to OCDB
1440 //
1441     if (changed) {
1442      AliCDBMetaData metaData;
1443      metaData.SetBeamPeriod(0);
1444      metaData.SetResponsible("Haavard Helstrup");
1445      metaData.SetAliRootVersion(ALIROOT_SVN_BRANCH);
1446      metaData.SetComment("Preprocessor AliTPC data base entries.");
1447
1448      Bool_t storeOK = Store("Calib", "AltroConfig", altroObjects, &metaData, 0, kFALSE);
1449      if ( !storeOK ) ++result;
1450     }  
1451
1452   altroObjects->Delete();
1453   delete altroObjects;
1454   
1455   return result;
1456 }