]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDPreprocessor.cxx
PMD DA informations : satya
[u/mrichter/AliRoot.git] / PMD / AliPMDPreprocessor.cxx
1 /***************************************************************************
2  * Copyright(c) 1998-1999, 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  
18  *   PMD Preproccessor Source Code      *
19     
20     
21      0 --> Run Successesfully
22      1 --> No pmd Alias is available
23  
24 **********************************************/
25  
26
27 #include <TFile.h>
28 #include <TTimeStamp.h>
29 #include <TObjString.h>
30 #include <TTree.h>
31
32 #include "AliLog.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"
40
41 ClassImp(AliPMDPreprocessor)
42
43 AliPMDPreprocessor::AliPMDPreprocessor(AliShuttleInterface* shuttle) :
44   AliPreprocessor("PMD", shuttle)
45 {
46   AddRunType("PHYSICS");
47   AddRunType("PEDESTAL");
48 }
49
50 //________________________________________________________ 
51 AliPMDPreprocessor::~AliPMDPreprocessor()
52 {
53   
54 }
55
56 //________________________________________________________ 
57 void AliPMDPreprocessor::Initialize(Int_t run, UInt_t startTime, UInt_t endTime)
58 {
59   AliPreprocessor::Initialize(run, startTime, endTime);
60   
61   AliInfo(Form("\nRun       : %d \nStart Time: %s \nEnd Time  : %s", run,
62                TTimeStamp(startTime).AsString(),
63                TTimeStamp(endTime).AsString()));
64   
65   fRun = run;
66   fStartTime = startTime;
67   fEndTime   = endTime;
68
69 }
70
71 //________________________________________________________ 
72 Bool_t AliPMDPreprocessor::ProcessDAQ()
73 {
74   TString sRunType = GetRunType();
75   Log(Form("RunType %s",sRunType.Data()));
76   if (sRunType !="PHYSICS" || sRunType != "PEDESTAL") 
77     {
78       return kFALSE;
79     }
80   return kTRUE;
81 }
82
83 //________________________________________________________ 
84 Bool_t AliPMDPreprocessor::StorePmdPED()
85 {
86   
87   AliPMDPedestal *pedestal = new AliPMDPedestal();
88   TList* gfsPmdPed = GetFileSources(kDAQ, "PMD_PED.root");
89   
90   if(!gfsPmdPed) 
91     {
92       Log(Form("PMDPED: No Shuttle List for PMD PED "));
93       return kFALSE;
94     }
95   else
96     {
97       AliInfo("PMDPED: Here's the list of sources for PMD_PED.root");
98       gfsPmdPed->Print();
99       
100       TIter iter(gfsPmdPed);
101       TObjString* srcPed;
102       TString nameOfFile;
103       UInt_t resultPmdPed = 0;
104       
105       while((srcPed = dynamic_cast<TObjString*> (iter.Next())))
106         {
107           nameOfFile = GetFile(kDAQ, "PMD_PED.root", srcPed->GetName());
108           if(nameOfFile.Length() == 0) 
109             {
110               Log(Form("PMDPED: Error retrieving file from source %s failed!", srcPed->GetName()));
111               delete gfsPmdPed;
112               return kFALSE;
113             }
114           
115           Log(Form("PMDPED: File with id PMD_PED.root got from %s", srcPed->GetName()));
116           
117           Int_t det, sm, row, col;
118           Float_t mean, rms;
119           
120           TFile *opnFile= new TFile(nameOfFile.Data());
121           if(!opnFile || !opnFile->IsOpen()) 
122             {
123               Log(Form("PMDPED: Error opening file with Id PMD_PED.root from source %s!", srcPed->GetName()));
124               return kFALSE;
125             } 
126           
127           TTree *tree = dynamic_cast<TTree *> (opnFile->Get("ped"));
128           if (!tree) 
129             {
130               Log("PMDPED: Could not find object \"ped\" in PED file!");
131                return kFALSE;
132             }
133           
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);
140           
141           Int_t nEntries = (Int_t) tree->GetEntries();
142           for(Int_t i = 0; i < nEntries; i++)
143             {
144               tree->GetEntry(i);
145               pedestal->SetPedMeanRms(det,sm,row,col,mean,rms);
146             }
147           opnFile->Close();
148           delete opnFile;
149         }
150       
151       AliCDBMetaData mdPED;
152       mdPED.SetBeamPeriod(0);
153       mdPED.SetResponsible("Satyajit Jena");
154       mdPED.SetComment("PMDPED: PMD preprocessor");
155       
156       resultPmdPed = Store("Calib","Ped", pedestal, &mdPED,0,kTRUE);
157       delete pedestal;
158       if(resultPmdPed==0)
159         {
160           Log("PMDPED: Error storing");                        
161           return kFALSE;
162         }
163       else
164         {
165           return kTRUE;
166         }
167
168     } 
169  
170 }
171
172 //________________________________________________________
173 Bool_t AliPMDPreprocessor::StorePmdGAIN()
174 {
175 TList* gfsPmdGain = GetFileSources(kDAQ, "PMDGAINS.root");
176
177 if(!gfsPmdGain)
178     {
179       Log(Form("PMDGAIN: No sources found for PMDGAINS.root!"));
180       return kFALSE;
181     }    
182   else
183     {
184       AliInfo("PMDGAIN: Here's the list of sources for PMDGAINS.root");
185       gfsPmdGain->Print();
186       
187       TIter iter(gfsPmdGain);
188       TObjString* source;
189       TString nameOfFile;
190       UInt_t result = 0;
191  
192       while((source = dynamic_cast<TObjString*> (iter.Next())))
193         {
194           nameOfFile = GetFile(kDAQ, "PMDGAINS.root", source->GetName());
195           if(nameOfFile.Length() == 0) 
196             {
197               Log(Form("PMDGAIN: Error retrieving file from source %s failed!", source->GetName()));
198               delete gfsPmdGain;
199               return kFALSE;
200             }
201           
202           Log(Form("PMDGAIN: File with id PMDGAINS.root got from %s", source->GetName()));
203           
204           Int_t det, sm, row, col;
205           Float_t gain;
206           
207           TFile *opnFile = new TFile(nameOfFile.Data());
208           if (!opnFile || !opnFile->IsOpen()) 
209             {
210               Log(Form("PMDGAIN: Error opening file with Id PMDGAINS.root from source %s!", source->GetName()));
211               return kFALSE;
212             }
213           
214           TTree *tree = dynamic_cast<TTree *> (opnFile->Get("ic"));
215           if (!tree) 
216             {
217               Log("PMDGAIN: Could not find object \"ic\" in DAQ file!");
218               return kFALSE;
219             }
220           
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);
226           
227           Int_t nEntries = (Int_t) tree->GetEntries();
228           AliPMDCalibData *calibda = new AliPMDCalibData();
229           
230           for(Int_t i = 0; i < nEntries; i++)
231             {
232               tree->GetEntry(i);
233               calibda->SetGainFact(det,sm,row,col,gain);
234             }
235           opnFile->Close();
236           delete opnFile;
237           
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);
243           delete calibda;
244         }
245       
246       if (result==0)
247         {
248           Log("PMDGAIN: Error storing");                        
249           return kFALSE;
250         }
251       else
252         {
253           return kTRUE;
254         }
255       
256     }
257 }
258
259
260 //___________________________________________________
261 Bool_t AliPMDPreprocessor::StorePmdHOT()
262 {
263   
264   AliPMDHotData *hotda = new AliPMDHotData();
265   TList* fsPmdHot = GetFileSources(kDAQ, "PMD_HOT.root");
266   
267   if(!fsPmdHot) 
268     {
269       Log(Form("PMDHOT: No sources found for PMD_HOT.root!"));
270       return kFALSE;
271     }
272   else
273     {
274       AliInfo("PMDHOT: Here's the list of sources for PMD_HOT.root");
275       fsPmdHot->Print();
276       
277       TIter iter(fsPmdHot);
278       TObjString* source;
279       UInt_t hotresult = 0;
280       TString nameOfFile;
281
282       while((source = dynamic_cast<TObjString*> (iter.Next())))
283         {
284           nameOfFile = GetFile(kDAQ, "PMD_HOT.root", source->GetName());
285           if(nameOfFile.Length() == 0) 
286             {
287               Log(Form("PMDHOT: Error retrieving file from source %s failed!", source->GetName()));
288               delete fsPmdHot;
289               return kFALSE;
290             }
291           
292           Log(Form("PMDHOT: File with id PMD_HOT.root got from %s", source->GetName()));
293           
294           Int_t det, sm, row, col;
295           Float_t flag;
296           
297           TFile *opnFile = new TFile(nameOfFile.Data());
298           if(!opnFile || !opnFile->IsOpen()) 
299             {
300               Log(Form("PMDHOT: Error opening file with Id PMD_HOT.root from source %s!", source->GetName()));
301               return kFALSE;
302             } 
303
304           TTree *tree = dynamic_cast<TTree *> (opnFile->Get("hot"));
305           if (!tree) 
306             {
307               Log("PMDHOT: Could not find object \"hot\" in DAQ file!");
308               return kFALSE;
309             }
310           
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);
316           
317           Int_t nEntries = (Int_t) tree->GetEntries();
318           for(Int_t j = 0; j < nEntries; j++)
319             {
320               tree->GetEntry(j);
321               hotda->SetHotChannel(det,sm,row,col,flag);
322             }
323           opnFile->Close();
324           delete opnFile;
325         }
326       
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);
332       delete hotda;
333       if(hotresult==0)
334         {
335           Log("PMDHOT: Error storing");                        
336           return kFALSE;
337         }
338       else
339         {
340           return kTRUE;
341         }
342       
343     }
344 }
345
346 //________________________________________________________
347 Bool_t AliPMDPreprocessor::StorePmdMEAN()
348 {
349   AliPMDMeanSm *smmeanda = new AliPMDMeanSm();
350   TList* gfsPmdMean = GetFileSources(kDAQ, "PMD_MEAN_SM.root");
351   
352   if(!gfsPmdMean) 
353     {
354       Log(Form("PMDMEAN: No sources found for PMD_MEAN_SM.root!"));
355       return kFALSE;
356     }
357   else
358     {
359       AliInfo("PMDMEAN: Here's the list of sources for PMD_MEAN_SM.root");
360       gfsPmdMean->Print();
361       
362       TIter iter(gfsPmdMean);
363       TObjString* sourc;
364       UInt_t storeMeanData = 0;
365       TString filenam;
366       
367       while((sourc=dynamic_cast<TObjString*> (iter.Next())))
368         {
369           filenam = GetFile(kDAQ, "PMD_MEAN_SM.root", sourc->GetName());
370           if(filenam.Length() == 0) 
371             {
372               Log(Form("PMDMEAN: Error retrieving file from source %s failed!", sourc->GetName()));
373               delete gfsPmdMean;
374               return kFALSE;
375             }
376           
377           Log(Form("PMDMEAN: File with id PMD_MEAN_SM.root got from %s", sourc->GetName()));
378           
379           Int_t det = 0, sm = 0;
380           Float_t smmean = 0.;
381           
382           TFile *opnFile = new TFile(filenam.Data());
383           if(!opnFile || !opnFile->IsOpen()) 
384             {
385               Log(Form("PMDMEAN: Error opening file with Id PMD_MEAN_SM.root from source %s!", sourc->GetName()));
386               return kFALSE;
387             }
388           
389           TTree *tree = dynamic_cast<TTree *> (opnFile->Get("mean"));
390           if (!tree) 
391             {
392               Log("PMDMEAN: Could not find object \"hot\" in DAQ file!");
393               return kFALSE;
394             }
395           
396           tree->SetBranchAddress("det",  &det);
397           tree->SetBranchAddress("sm",   &sm);
398           tree->SetBranchAddress("smmean", &smmean);
399           
400           Int_t nEntries = (Int_t) tree->GetEntries();
401           for(Int_t j = 0; j < nEntries; j++)
402             {
403               tree->GetEntry(j);
404               smmeanda->SetMeanSm(det,sm,smmean);
405             }
406           opnFile->Close();
407           delete opnFile;
408         }
409       
410       AliCDBMetaData mdMEAN;
411       mdMEAN.SetBeamPeriod(0);
412       mdMEAN.SetResponsible("Satyajit Jena");
413       mdMEAN.SetComment("PMDMEAN: PMD preprocessor");
414       
415       storeMeanData = Store("Calib","SMMEAN", smmeanda, &mdMEAN);
416       delete smmeanda;
417
418       if(storeMeanData==0)
419         {
420           Log("PMDMEAN: Error storing");                        
421           return kFALSE;
422         }
423       else
424         {
425           return kTRUE;
426         }
427     }
428   
429 }
430
431 //_____________________________________________________________
432 Bool_t AliPMDPreprocessor::StorePmdDCS(TMap *sDaqAM)
433 {
434         
435   AliCDBMetaData mdDCS;
436   mdDCS.SetResponsible("Satyajit Jena");
437   mdDCS.SetComment("DCS data for PMD");
438
439   Bool_t resStore = kFALSE;
440   resStore = StoreReferenceData("DCS","Data",sDaqAM,&mdDCS);
441
442   if(resStore == 0)
443     {
444       Log("PMDDP: Error storing");                        
445       return kFALSE;
446     }
447   else
448     {
449       return kTRUE;
450     }
451
452 }
453
454 //_____________________________________________________________________
455
456 UInt_t AliPMDPreprocessor::Process(TMap* pmdDaqAliasMap)
457 {
458   
459   if(!pmdDaqAliasMap)
460     {
461       return 1;
462     }
463   
464   TString runType = GetRunType();
465
466   if(runType == "PEDESTAL")
467     {
468  
469       Log(Form("------------------ PMD Pedestal --------------"));
470       Bool_t pmdPed = StorePmdPED();
471       if(!pmdPed)
472         {
473           Log(Form("ERROR:  Couldn't write PMD pedestal file"));
474           return 1;
475         }
476       else
477         {
478           Log(Form("Storing of PMD Pedestal File is Successful"));
479           return 0;
480         }
481     } 
482   
483   else if (runType == "PHYSICS")
484       {
485         Log(Form("------------------- PMD GAIN----------------"));      
486         Bool_t pmdGAIN = StorePmdGAIN();
487         if (!pmdGAIN)
488           {
489             Log(Form("ERROR:  Couldn't write PMD GAIN file"));
490            
491           }
492         else
493           {
494             Log(Form("Storing of PMD GAIN File is Successful"));
495            
496           }
497         
498         Log(Form("------------------- PMD HOT ----------------"));
499         Bool_t pmdHOT = StorePmdHOT();
500         if (!pmdHOT)
501           {
502             Log(Form("ERROR:  Couldn't write PMD HOT file"));
503            
504           }
505         else
506           {
507             Log(Form("Storing of PMD HOT File is Successful"));
508            
509           }
510         
511         Log(Form("------------------- SM MEAN ----------------"));
512         Bool_t pmdMEAN = StorePmdMEAN();
513         if (!pmdMEAN)
514           {
515             Log(Form("ERROR:  Couldn't write PMD SM MEAN file"));
516            
517           }
518         else
519           {
520             Log(Form("Storing of PMD SM MEAN File is Successful"));
521            
522           }
523         
524         Log(Form("------------------- DCS DP -----------------"));
525         Bool_t pmdDCS = StorePmdDCS(pmdDaqAliasMap);
526         if (!pmdDCS)
527           {
528             Log(Form("ERROR:  Couldn't write PMD DCS DP file"));
529            
530           }
531         else
532           {
533             Log(Form("Storing of PMD DCS dp is Successul"));
534            
535           }
536         
537         if (pmdGAIN || pmdHOT || pmdMEAN || pmdDCS) 
538           return 0;
539         else 
540           return 1;
541         
542       }
543   return 2;
544 }
545
546