1 /**************************************************************************
2 * Copyright(c) 2007, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 #include "AliTPCPreprocessor.h"
18 #include "AliShuttleInterface.h"
20 #include "AliCDBMetaData.h"
21 #include "AliDCSValue.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"
34 #include <TTimeStamp.h>
36 const Int_t kValCutTemp = 100; // discard temperatures > 100 degrees
37 const Int_t kDiffCutTemp = 5; // discard temperature differences > 5 degrees
38 const TString kPedestalRunType = "PEDESTAL_RUN"; // pedestal run identifier
39 const TString kPulserRunType = "PULSER_RUN"; // pulser run identifier
40 const TString kPhysicsRunType = "PHYSICS"; // physics run identifier
41 const TString kStandAloneRunType = "STANDALONE"; // standalone run identifier
42 const TString kDaqRunType = "DAQ"; // DAQ run identifier
43 const TString kAmandaTemp = "tpc_PT_%d.Temperature"; // Amanda string for temperature entries
44 //const Double_t kFitFraction = 0.7; // Fraction of DCS sensor fits required
45 const Double_t kFitFraction = -1.0; // Don't require minimum number of fits in commissioning run
48 // This class is the SHUTTLE preprocessor for the TPC detector.
51 ClassImp(AliTPCPreprocessor)
53 //______________________________________________________________________________________________
54 AliTPCPreprocessor::AliTPCPreprocessor(AliShuttleInterface* shuttle) :
55 AliPreprocessor("TPC",shuttle),
56 fConfEnv(0), fTemp(0), fHighVoltage(0), fConfigOK(kTRUE), fROC(0)
59 fROC = AliTPCROC::Instance();
61 //______________________________________________________________________________________________
62 // AliTPCPreprocessor::AliTPCPreprocessor(const AliTPCPreprocessor& org) :
63 // AliPreprocessor(org),
64 // fConfEnv(0), fTemp(0), fHighVoltage(0), fConfigOK(kTRUE)
66 // // copy constructor not implemented
67 // // -- missing underlying copy constructor in AliPreprocessor
69 // Fatal("AliTPCPreprocessor", "copy constructor not implemented");
71 // // fTemp = new AliTPCSensorTempArray(*(org.fTemp));
74 //______________________________________________________________________________________________
75 AliTPCPreprocessor::~AliTPCPreprocessor()
82 //______________________________________________________________________________________________
83 AliTPCPreprocessor& AliTPCPreprocessor::operator = (const AliTPCPreprocessor& )
85 Fatal("operator =", "assignment operator not implemented");
90 //______________________________________________________________________________________________
91 void AliTPCPreprocessor::Initialize(Int_t run, UInt_t startTime,
94 // Creates AliTestDataDCS object
96 AliPreprocessor::Initialize(run, startTime, endTime);
98 AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %s", run,
99 TTimeStamp(startTime).AsString(),
100 TTimeStamp(endTime).AsString()));
102 // Preprocessor configuration
104 AliCDBEntry* entry = GetFromOCDB("Config", "Preprocessor");
105 if (entry) fConfEnv = (TEnv*) entry->GetObject();
107 Log("AliTPCPreprocsessor: Preprocessor Config OCDB entry missing.\n");
112 // Temperature sensors
115 entry = GetFromOCDB("Config", "Temperature");
116 if (entry) confTree = (TTree*) entry->GetObject();
118 Log("AliTPCPreprocsessor: Temperature Config OCDB entry missing.\n");
122 fTemp = new AliTPCSensorTempArray(fStartTime, fEndTime, confTree, kAmandaTemp);
123 fTemp->SetValCut(kValCutTemp);
124 fTemp->SetDiffCut(kDiffCutTemp);
127 // High voltage measurements
129 if (false) { // high voltage maps not yet implemented...
132 entry = GetFromOCDB("Config", "HighVoltage");
133 if (entry) confTree = (TTree*) entry->GetObject();
135 Log("AliTPCPreprocsessor: High Voltage Config OCDB entry missing.\n");
139 fHighVoltage = new AliDCSSensorArray(fStartTime, fEndTime, confTree);
143 //______________________________________________________________________________________________
144 UInt_t AliTPCPreprocessor::Process(TMap* dcsAliasMap)
146 // Fills data into TPC calibrations objects
148 // Amanda servers provide information directly through dcsAliasMap
150 if (!dcsAliasMap) return 9;
151 if (dcsAliasMap->GetEntries() == 0 ) return 9;
152 if (!fConfigOK) return 9;
154 TString runType = GetRunType();
156 // Temperature sensors are processed by AliTPCCalTemp
159 UInt_t tempResult = MapTemperature(dcsAliasMap);
160 UInt_t result=tempResult;
162 // High Voltage recordings
165 if (false) { // High Voltage maps not yet implemented..
167 UInt_t hvResult = MapHighVoltage(dcsAliasMap);
172 // Other calibration information will be retrieved through FXS files
174 // TList* fileSourcesDAQ = GetFile(AliShuttleInterface::kDAQ, "pedestals");
175 // const char* fileNamePed = GetFile(AliShuttleInterface::kDAQ, "pedestals", "LDC1");
177 // TList* fileSourcesHLT = GetFile(AliShuttleInterface::kHLT, "calib");
178 // const char* fileNameHLT = GetFile(AliShuttleInterface::kHLT, "calib", "LDC1");
182 if(runType == kPedestalRunType) {
183 Int_t numSources = 1;
184 Int_t pedestalSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
185 TString source = fConfEnv->GetValue("Pedestal","DAQ");
187 if ( source == "HLT") pedestalSource[0] = AliShuttleInterface::kHLT;
188 if (!GetHLTStatus()) pedestalSource[0] = AliShuttleInterface::kDAQ;
189 if (source == "HLTDAQ" ) {
191 pedestalSource[0] = AliShuttleInterface::kHLT;
192 pedestalSource[1] = AliShuttleInterface::kDAQ;
194 if (source == "DAQHLT" ) numSources=2;
195 UInt_t pedestalResult=0;
196 for (Int_t i=0; i<numSources; i++ ) {
197 UInt_t pedestalResult = ExtractPedestals(pedestalSource[i]);
198 if ( pedestalResult == 0 ) break;
200 result += pedestalResult;
204 // pulser trigger processing
206 if(runType == kPulserRunType) {
207 Int_t numSources = 1;
208 Int_t pulserSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
209 TString source = fConfEnv->GetValue("Pulser","DAQ");
211 if ( source == "HLT") pulserSource[0] = AliShuttleInterface::kHLT;
212 if (!GetHLTStatus()) pulserSource[0] = AliShuttleInterface::kDAQ;
213 if (source == "HLTDAQ" ) {
215 pulserSource[0] = AliShuttleInterface::kHLT;
216 pulserSource[1] = AliShuttleInterface::kDAQ;
218 if (source == "DAQHLT" ) numSources=2;
219 UInt_t pulserResult=0;
220 for (Int_t i=0; i<numSources; i++ ) {
221 pulserResult = ExtractPulser(pulserSource[i]);
222 if ( pulserResult == 0 ) break;
224 result += pulserResult;
230 // Central Electrode processing
232 if( runType == kPhysicsRunType || runType == kStandAloneRunType ||
233 runType == kDaqRunType ) {
235 Int_t numSources = 1;
236 Int_t ceSource[2] = {AliShuttleInterface::kDAQ,AliShuttleInterface::kHLT} ;
237 TString source = fConfEnv->GetValue("CE","DAQ");
239 if ( source == "HLT") ceSource[0] = AliShuttleInterface::kHLT;
240 if (!GetHLTStatus()) ceSource[0] = AliShuttleInterface::kDAQ;
241 if (source == "HLTDAQ" ) {
243 ceSource[0] = AliShuttleInterface::kHLT;
244 ceSource[1] = AliShuttleInterface::kDAQ;
246 if (source == "DAQHLT" ) numSources=2;
248 for (Int_t i=0; i<numSources; i++ ) {
249 ceResult = ExtractCE(ceSource[i]);
250 if ( ceResult == 0 ) break;
258 //______________________________________________________________________________________________
259 UInt_t AliTPCPreprocessor::MapTemperature(TMap* dcsAliasMap)
262 // extract DCS temperature maps. Perform fits to save space
265 TMap *map = fTemp->ExtractDCS(dcsAliasMap);
267 fTemp->MakeSplineFit(map);
268 Double_t fitFraction = 1.0*fTemp->NumFits()/fTemp->NumSensors();
269 if (fitFraction > kFitFraction ) {
270 AliInfo(Form("Temperature values extracted, fits performed.\n"));
272 Log ("Too few temperature maps fitted. \n");
276 Log("No temperature map extracted. \n");
280 // Now store the final CDB file
283 AliCDBMetaData metaData;
284 metaData.SetBeamPeriod(0);
285 metaData.SetResponsible("Haavard Helstrup");
286 metaData.SetComment("Preprocessor AliTPC data base entries.");
288 Bool_t storeOK = Store("Calib", "Temperature", fTemp, &metaData, 0, kFALSE);
289 if ( !storeOK ) result=1;
297 //______________________________________________________________________________________________
298 UInt_t AliTPCPreprocessor::MapHighVoltage(TMap* dcsAliasMap)
301 // extract DCS HV maps. Perform fits to save space
304 TMap *map = fHighVoltage->ExtractDCS(dcsAliasMap);
306 fHighVoltage->MakeSplineFit(map);
307 Double_t fitFraction = 1.0*fHighVoltage->NumFits()/fHighVoltage->NumSensors();
308 if (fitFraction > kFitFraction ) {
309 AliInfo(Form("High volatge recordings extracted, fits performed.\n"));
311 Log ("Too few high voltage recordings fitted. \n");
315 Log("No high voltage recordings extracted. \n");
319 // Now store the final CDB file
322 AliCDBMetaData metaData;
323 metaData.SetBeamPeriod(0);
324 metaData.SetResponsible("Haavard Helstrup");
325 metaData.SetComment("Preprocessor AliTPC data base entries.");
327 Bool_t storeOK = Store("Calib", "HighVoltage", fHighVoltage, &metaData, 0, kFALSE);
328 if ( !storeOK ) result=1;
337 //______________________________________________________________________________________________
339 UInt_t AliTPCPreprocessor::ExtractPedestals(Int_t sourceFXS)
342 // Read pedestal file from file exchage server
343 // Keep original entry from OCDB in case no new pedestals are available
345 AliTPCCalPad *calPadPed=0;
346 AliCDBEntry* entry = GetFromOCDB("Calib", "Pedestals");
347 if (entry) calPadPed = (AliTPCCalPad*)entry->GetObject();
348 if ( calPadPed==NULL ) {
349 Log("AliTPCPreprocsessor: No previous TPC pedestal entry available.\n");
350 calPadPed = new AliTPCCalPad("PedestalsMean","PedestalsMean");
353 AliTPCCalPad *calPadRMS=0;
354 entry = GetFromOCDB("Calib", "Noise");
355 if (entry) calPadRMS = (AliTPCCalPad*)entry->GetObject();
356 if ( calPadRMS==NULL ) {
357 Log("AliTPCPreprocsessor: No previous TPC noise entry available.\n");
358 calPadRMS = new AliTPCCalPad("PedestalsRMS","PedestalsRMS");
364 Int_t nSectors = fROC->GetNSectors();
365 TList* list = GetFileSources(sourceFXS,"pedestals");
367 if (list && list->GetEntries()>0) {
369 // loop through all files from LDCs
371 Bool_t changed=false;
373 while (list->At(index)!=NULL) {
374 TObjString* fileNameEntry = (TObjString*) list->At(index);
375 if (fileNameEntry!=NULL) {
376 TString fileName = GetFile(sourceFXS, "pedestals",
377 fileNameEntry->GetString().Data());
378 TFile *f = TFile::Open(fileName);
380 Log ("Error opening pedestal file.");
384 AliTPCCalibPedestal *calPed;
385 f->GetObject("AliTPCCalibPedestal",calPed);
387 Log ("No pedestal calibration object in file.");
392 // replace entries for the sectors available in the present file
395 for (Int_t sector=0; sector<nSectors; sector++) {
396 AliTPCCalROC *rocPed=calPed->GetCalRocPedestal(sector, kFALSE);
397 if ( rocPed ) calPadPed->SetCalROC(rocPed,sector);
398 AliTPCCalROC *rocRMS=calPed->GetCalRocRMS(sector, kFALSE);
399 if ( rocRMS ) calPadRMS->SetCalROC(rocRMS,sector);
405 // Store updated pedestal entry to OCDB
408 AliCDBMetaData metaData;
409 metaData.SetBeamPeriod(0);
410 metaData.SetResponsible("Haavard Helstrup");
411 metaData.SetComment("Preprocessor AliTPC data base entries.");
413 Bool_t storeOK = Store("Calib", "Pedestals", calPadPed, &metaData, 0, kTRUE);
414 if ( !storeOK ) ++result;
415 storeOK = Store("Calib", "PadNoise", calPadRMS, &metaData, 0, kTRUE);
416 if ( !storeOK ) ++result;
419 Log ("Error: no entries in input file list!");
426 //______________________________________________________________________________________________
429 UInt_t AliTPCPreprocessor::ExtractPulser(Int_t sourceFXS)
432 // Read pulser calibration file from file exchage server
433 // Keep original entry from OCDB in case no new pulser calibration is available
435 TObjArray *pulserObjects=0;
436 AliTPCCalPad *pulserTmean=0;
437 AliTPCCalPad *pulserTrms=0;
438 AliTPCCalPad *pulserQmean=0;
439 AliCDBEntry* entry = GetFromOCDB("Calib", "Pulser");
440 if (entry) pulserObjects = (TObjArray*)entry->GetObject();
441 if ( pulserObjects==NULL ) {
442 Log("AliTPCPreprocsessor: No previous TPC pulser entry available.\n");
443 pulserObjects = new TObjArray;
446 pulserTmean = (AliTPCCalPad*)pulserObjects->FindObject("PulserTmean");
447 if ( !pulserTmean ) {
448 pulserTmean = new AliTPCCalPad("PulserTmean","PulserTmean");
449 pulserObjects->Add(pulserTmean);
451 pulserTrms = (AliTPCCalPad*)pulserObjects->FindObject("PulserTrms");
453 pulserTrms = new AliTPCCalPad("PulserTrms","PulserTrms");
454 pulserObjects->Add(pulserTrms);
456 pulserQmean = (AliTPCCalPad*)pulserObjects->FindObject("PulserQmean");
457 if ( !pulserQmean ) {
458 pulserQmean = new AliTPCCalPad("PulserQmean","PulserQmean");
459 pulserObjects->Add(pulserQmean);
465 Int_t nSectors = fROC->GetNSectors();
466 TList* list = GetFileSources(sourceFXS,"pulser");
468 if (list && list->GetEntries()>0) {
470 // loop through all files from LDCs
472 Bool_t changed=false;
474 while (list->At(index)!=NULL) {
475 TObjString* fileNameEntry = (TObjString*) list->At(index);
476 if (fileNameEntry!=NULL) {
477 TString fileName = GetFile(sourceFXS, "pulser",
478 fileNameEntry->GetString().Data());
479 TFile *f = TFile::Open(fileName);
481 Log ("Error opening pulser file.");
485 AliTPCCalibPulser *calPulser;
486 f->GetObject("AliTPCCalibPulser",calPulser);
488 Log ("No pulser calibration object in file.");
493 // replace entries for the sectors available in the present file
496 for (Int_t sector=0; sector<nSectors; sector++) {
497 AliTPCCalROC *rocTmean=calPulser->GetCalRocT0(sector);
498 if ( rocTmean ) pulserTmean->SetCalROC(rocTmean,sector);
499 AliTPCCalROC *rocTrms=calPulser->GetCalRocRMS(sector);
500 if ( rocTrms ) pulserTrms->SetCalROC(rocTrms,sector);
501 AliTPCCalROC *rocQmean=calPulser->GetCalRocQ(sector);
502 if ( rocQmean ) pulserQmean->SetCalROC(rocQmean,sector);
508 // Store updated pedestal entry to OCDB
511 AliCDBMetaData metaData;
512 metaData.SetBeamPeriod(0);
513 metaData.SetResponsible("Haavard Helstrup");
514 metaData.SetComment("Preprocessor AliTPC data base entries.");
516 Bool_t storeOK = Store("Calib", "Pulser", pulserObjects, &metaData, 0, kTRUE);
517 if ( !storeOK ) ++result;
520 Log ("Error: no entries in input file list!");
527 UInt_t AliTPCPreprocessor::ExtractCE(Int_t sourceFXS)
530 // Read Central Electrode file from file exchage server
531 // Keep original entry from OCDB in case no new CE calibration is available
533 TObjArray *ceObjects=0;
534 AliTPCCalPad *ceTmean=0;
535 AliTPCCalPad *ceTrms=0;
536 AliTPCCalPad *ceQmean=0;
537 AliCDBEntry* entry = GetFromOCDB("Calib", "CE");
538 if (entry) ceObjects = (TObjArray*)entry->GetObject();
539 if ( ceObjects==NULL ) {
540 Log("AliTPCPreprocsessor: No previous TPC central electrode entry available.\n");
541 ceObjects = new TObjArray;
544 ceTmean = (AliTPCCalPad*)ceObjects->FindObject("CETmean");
546 ceTmean = new AliTPCCalPad("CETmean","CETmean");
547 ceObjects->Add(ceTmean);
549 ceTrms = (AliTPCCalPad*)ceObjects->FindObject("CETrms");
551 ceTrms = new AliTPCCalPad("CETrms","CETrms");
552 ceObjects->Add(ceTrms);
554 ceQmean = (AliTPCCalPad*)ceObjects->FindObject("CEQmean");
556 ceQmean = new AliTPCCalPad("CEQmean","CEQmean");
557 ceObjects->Add(ceQmean);
563 Int_t nSectors = fROC->GetNSectors();
564 TList* list = GetFileSources(sourceFXS,"CE");
566 if (list && list->GetEntries()>0) {
568 // loop through all files from LDCs
571 while (list->At(index)!=NULL) {
572 TObjString* fileNameEntry = (TObjString*) list->At(index);
573 if (fileNameEntry!=NULL) {
574 TString fileName = GetFile(sourceFXS, "CE",
575 fileNameEntry->GetString().Data());
576 TFile *f = TFile::Open(fileName);
578 Log ("Error opening central electrode file.");
582 AliTPCCalibCE *calCE;
583 f->GetObject("AliTPCCalibCE",calCE);
585 // replace entries for the sectors available in the present file
587 for (Int_t sector=0; sector<nSectors; sector++) {
588 AliTPCCalROC *rocTmean=calCE->GetCalRocT0(sector);
589 if ( rocTmean ) ceTmean->SetCalROC(rocTmean,sector);
590 AliTPCCalROC *rocTrms=calCE->GetCalRocRMS(sector);
591 if ( rocTrms ) ceTrms->SetCalROC(rocTrms,sector);
592 AliTPCCalROC *rocQmean=calCE->GetCalRocQ(sector);
593 if ( rocQmean ) ceQmean->SetCalROC(rocQmean,sector);
599 // Store updated pedestal entry to OCDB
601 AliCDBMetaData metaData;
602 metaData.SetBeamPeriod(0);
603 metaData.SetResponsible("Haavard Helstrup");
604 metaData.SetComment("Preprocessor AliTPC data base entries.");
606 Bool_t storeOK = Store("Calib", "CE", ceObjects, &metaData, 0, kTRUE);
607 if ( !storeOK ) ++result;
610 Log ("Error: no entries!");