1 /***************************************************************************
2 * Copyright(c) 1998-1999, 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 **************************************************************************/
16 /*********************************************
18 * PMD Preproccessor Source Code *
21 0 --> Run Successesfully
22 1 --> No pmd Alias is available
24 **********************************************/
28 #include <TTimeStamp.h>
29 #include <TObjString.h>
33 #include "AliShuttleInterface.h"
34 #include "AliCDBMetaData.h"
35 #include "AliPMDCalibData.h"
36 #include "AliPMDHotData.h"
37 #include "AliPMDMeanSm.h"
38 #include "AliPMDPedestal.h"
39 #include "AliPMDPreprocessor.h"
41 ClassImp(AliPMDPreprocessor)
43 AliPMDPreprocessor::AliPMDPreprocessor(AliShuttleInterface* shuttle) :
44 AliPreprocessor("PMD", shuttle)
46 AddRunType("PHYSICS");
47 AddRunType("PEDESTAL");
50 //________________________________________________________
51 AliPMDPreprocessor::~AliPMDPreprocessor()
56 //________________________________________________________
57 void AliPMDPreprocessor::Initialize(Int_t run, UInt_t startTime, UInt_t endTime)
59 AliPreprocessor::Initialize(run, startTime, endTime);
61 AliInfo(Form("\nRun : %d \nStart Time: %s \nEnd Time : %s", run,
62 TTimeStamp(startTime).AsString(),
63 TTimeStamp(endTime).AsString()));
66 fStartTime = startTime;
71 //________________________________________________________
72 Bool_t AliPMDPreprocessor::ProcessDAQ()
74 TString sRunType = GetRunType();
75 Log(Form("RunType %s",sRunType.Data()));
76 if (sRunType !="PHYSICS" || sRunType != "PEDESTAL")
83 //________________________________________________________
84 Bool_t AliPMDPreprocessor::StorePmdPED()
87 AliPMDPedestal *pedestal = new AliPMDPedestal();
88 TList* gfsPmdPed = GetFileSources(kDAQ, "PMD_PED.root");
92 Log(Form("PMDPED: No Shuttle List for PMD PED "));
97 AliInfo("PMDPED: Here's the list of sources for PMD_PED.root");
100 TIter iter(gfsPmdPed);
103 UInt_t resultPmdPed = 0;
105 while((srcPed = dynamic_cast<TObjString*> (iter.Next())))
107 nameOfFile = GetFile(kDAQ, "PMD_PED.root", srcPed->GetName());
108 if(nameOfFile.Length() == 0)
110 Log(Form("PMDPED: Error retrieving file from source %s failed!", srcPed->GetName()));
115 Log(Form("PMDPED: File with id PMD_PED.root got from %s", srcPed->GetName()));
117 Int_t det, sm, row, col;
120 TFile *opnFile= new TFile(nameOfFile.Data());
121 if(!opnFile || !opnFile->IsOpen())
123 Log(Form("PMDPED: Error opening file with Id PMD_PED.root from source %s!", srcPed->GetName()));
127 TTree *tree = dynamic_cast<TTree *> (opnFile->Get("ped"));
130 Log("PMDPED: Could not find object \"ped\" in PED file!");
134 tree->SetBranchAddress("det", &det);
135 tree->SetBranchAddress("sm", &sm);
136 tree->SetBranchAddress("row", &row);
137 tree->SetBranchAddress("col", &col);
138 tree->SetBranchAddress("mean", &mean);
139 tree->SetBranchAddress("rms", &rms);
141 Int_t nEntries = (Int_t) tree->GetEntries();
142 for(Int_t i = 0; i < nEntries; i++)
145 pedestal->SetPedMeanRms(det,sm,row,col,mean,rms);
151 AliCDBMetaData mdPED;
152 mdPED.SetBeamPeriod(0);
153 mdPED.SetResponsible("Satyajit Jena");
154 mdPED.SetComment("PMDPED: PMD preprocessor");
156 resultPmdPed = Store("Calib","Ped", pedestal, &mdPED,0,kTRUE);
160 Log("PMDPED: Error storing");
172 //________________________________________________________
173 Bool_t AliPMDPreprocessor::StorePmdGAIN()
175 TList* gfsPmdGain = GetFileSources(kDAQ, "PMDGAINS.root");
179 Log(Form("PMDGAIN: No sources found for PMDGAINS.root!"));
184 AliInfo("PMDGAIN: Here's the list of sources for PMDGAINS.root");
187 TIter iter(gfsPmdGain);
192 while((source = dynamic_cast<TObjString*> (iter.Next())))
194 nameOfFile = GetFile(kDAQ, "PMDGAINS.root", source->GetName());
195 if(nameOfFile.Length() == 0)
197 Log(Form("PMDGAIN: Error retrieving file from source %s failed!", source->GetName()));
202 Log(Form("PMDGAIN: File with id PMDGAINS.root got from %s", source->GetName()));
204 Int_t det, sm, row, col;
207 TFile *opnFile = new TFile(nameOfFile.Data());
208 if (!opnFile || !opnFile->IsOpen())
210 Log(Form("PMDGAIN: Error opening file with Id PMDGAINS.root from source %s!", source->GetName()));
214 TTree *tree = dynamic_cast<TTree *> (opnFile->Get("ic"));
217 Log("PMDGAIN: Could not find object \"ic\" in DAQ file!");
221 tree->SetBranchAddress("det", &det);
222 tree->SetBranchAddress("sm", &sm);
223 tree->SetBranchAddress("row", &row);
224 tree->SetBranchAddress("col", &col);
225 tree->SetBranchAddress("gain", &gain);
227 Int_t nEntries = (Int_t) tree->GetEntries();
228 AliPMDCalibData *calibda = new AliPMDCalibData();
230 for(Int_t i = 0; i < nEntries; i++)
233 calibda->SetGainFact(det,sm,row,col,gain);
238 AliCDBMetaData mdGAIN;
239 mdGAIN.SetBeamPeriod(0);
240 mdGAIN.SetResponsible("Satyajit Jena");
241 mdGAIN.SetComment("PMDGAIN: PMD GAIN Data");
242 result = Store("Calib","Gain", calibda, &mdGAIN);
248 Log("PMDGAIN: Error storing");
260 //___________________________________________________
261 Bool_t AliPMDPreprocessor::StorePmdHOT()
264 AliPMDHotData *hotda = new AliPMDHotData();
265 TList* fsPmdHot = GetFileSources(kDAQ, "PMD_HOT.root");
269 Log(Form("PMDHOT: No sources found for PMD_HOT.root!"));
274 AliInfo("PMDHOT: Here's the list of sources for PMD_HOT.root");
277 TIter iter(fsPmdHot);
279 UInt_t hotresult = 0;
282 while((source = dynamic_cast<TObjString*> (iter.Next())))
284 nameOfFile = GetFile(kDAQ, "PMD_HOT.root", source->GetName());
285 if(nameOfFile.Length() == 0)
287 Log(Form("PMDHOT: Error retrieving file from source %s failed!", source->GetName()));
292 Log(Form("PMDHOT: File with id PMD_HOT.root got from %s", source->GetName()));
294 Int_t det, sm, row, col;
297 TFile *opnFile = new TFile(nameOfFile.Data());
298 if(!opnFile || !opnFile->IsOpen())
300 Log(Form("PMDHOT: Error opening file with Id PMD_HOT.root from source %s!", source->GetName()));
304 TTree *tree = dynamic_cast<TTree *> (opnFile->Get("hot"));
307 Log("PMDHOT: Could not find object \"hot\" in DAQ file!");
311 tree->SetBranchAddress("det", &det);
312 tree->SetBranchAddress("sm", &sm);
313 tree->SetBranchAddress("row", &row);
314 tree->SetBranchAddress("col", &col);
315 tree->SetBranchAddress("flag", &flag);
317 Int_t nEntries = (Int_t) tree->GetEntries();
318 for(Int_t j = 0; j < nEntries; j++)
321 hotda->SetHotChannel(det,sm,row,col,flag);
327 AliCDBMetaData metaData;
328 metaData.SetBeamPeriod(0);
329 metaData.SetResponsible("Satyajit Jena");
330 metaData.SetComment("PMDHOT: PMD preprocessor");
331 hotresult = Store("Calib","Hot", hotda, &metaData);
335 Log("PMDHOT: Error storing");
346 //________________________________________________________
347 Bool_t AliPMDPreprocessor::StorePmdMEAN()
349 AliPMDMeanSm *smmeanda = new AliPMDMeanSm();
350 TList* gfsPmdMean = GetFileSources(kDAQ, "PMD_MEAN_SM.root");
354 Log(Form("PMDMEAN: No sources found for PMD_MEAN_SM.root!"));
359 AliInfo("PMDMEAN: Here's the list of sources for PMD_MEAN_SM.root");
362 TIter iter(gfsPmdMean);
364 UInt_t storeMeanData = 0;
367 while((sourc=dynamic_cast<TObjString*> (iter.Next())))
369 filenam = GetFile(kDAQ, "PMD_MEAN_SM.root", sourc->GetName());
370 if(filenam.Length() == 0)
372 Log(Form("PMDMEAN: Error retrieving file from source %s failed!", sourc->GetName()));
377 Log(Form("PMDMEAN: File with id PMD_MEAN_SM.root got from %s", sourc->GetName()));
379 Int_t det = 0, sm = 0;
382 TFile *opnFile = new TFile(filenam.Data());
383 if(!opnFile || !opnFile->IsOpen())
385 Log(Form("PMDMEAN: Error opening file with Id PMD_MEAN_SM.root from source %s!", sourc->GetName()));
389 TTree *tree = dynamic_cast<TTree *> (opnFile->Get("mean"));
392 Log("PMDMEAN: Could not find object \"hot\" in DAQ file!");
396 tree->SetBranchAddress("det", &det);
397 tree->SetBranchAddress("sm", &sm);
398 tree->SetBranchAddress("smmean", &smmean);
400 Int_t nEntries = (Int_t) tree->GetEntries();
401 for(Int_t j = 0; j < nEntries; j++)
404 smmeanda->SetMeanSm(det,sm,smmean);
410 AliCDBMetaData mdMEAN;
411 mdMEAN.SetBeamPeriod(0);
412 mdMEAN.SetResponsible("Satyajit Jena");
413 mdMEAN.SetComment("PMDMEAN: PMD preprocessor");
415 storeMeanData = Store("Calib","SMMEAN", smmeanda, &mdMEAN);
420 Log("PMDMEAN: Error storing");
431 //_____________________________________________________________
432 Bool_t AliPMDPreprocessor::StorePmdDCS(TMap *sDaqAM)
435 AliCDBMetaData mdDCS;
436 mdDCS.SetResponsible("Satyajit Jena");
437 mdDCS.SetComment("DCS data for PMD");
439 Bool_t resStore = kFALSE;
440 resStore = StoreReferenceData("DCS","Data",sDaqAM,&mdDCS);
444 Log("PMDDP: Error storing");
454 //_____________________________________________________________________
456 UInt_t AliPMDPreprocessor::Process(TMap* pmdDaqAliasMap)
464 TString runType = GetRunType();
466 if(runType == "PEDESTAL")
469 Log(Form("------------------ PMD Pedestal --------------"));
470 Bool_t pmdPed = StorePmdPED();
473 Log(Form("ERROR: Couldn't write PMD pedestal file"));
478 Log(Form("Storing of PMD Pedestal File is Successful"));
483 else if (runType == "PHYSICS")
485 Log(Form("------------------- PMD GAIN----------------"));
486 Bool_t pmdGAIN = StorePmdGAIN();
489 Log(Form("ERROR: Couldn't write PMD GAIN file"));
494 Log(Form("Storing of PMD GAIN File is Successful"));
498 Log(Form("------------------- PMD HOT ----------------"));
499 Bool_t pmdHOT = StorePmdHOT();
502 Log(Form("ERROR: Couldn't write PMD HOT file"));
507 Log(Form("Storing of PMD HOT File is Successful"));
511 Log(Form("------------------- SM MEAN ----------------"));
512 Bool_t pmdMEAN = StorePmdMEAN();
515 Log(Form("ERROR: Couldn't write PMD SM MEAN file"));
520 Log(Form("Storing of PMD SM MEAN File is Successful"));
524 Log(Form("------------------- DCS DP -----------------"));
525 Bool_t pmdDCS = StorePmdDCS(pmdDaqAliasMap);
528 Log(Form("ERROR: Couldn't write PMD DCS DP file"));
533 Log(Form("Storing of PMD DCS dp is Successul"));
537 if (pmdGAIN || pmdHOT || pmdMEAN || pmdDCS)