]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDPreprocessor.cxx
Several changes:
[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 //  Source File : AliPMDPreprocessor.cxx               //
18 //                                                     //
19 //                                                     //
20 //-----------------------------------------------------//
21
22 // --- ROOT system
23 #include <TFile.h>
24 #include <TTimeStamp.h>
25 #include <TObjString.h>
26 #include <TTree.h>
27
28 #include "AliLog.h"
29 #include "AliShuttleInterface.h"
30 #include "AliCDBMetaData.h"
31 #include "AliPMDCalibData.h"
32 #include "AliPMDPedestal.h"
33 #include "AliPMDPreprocessor.h"
34
35
36 ClassImp(AliPMDPreprocessor)
37   
38 //______________________________________________________________________________________________
39 AliPMDPreprocessor::AliPMDPreprocessor(AliShuttleInterface* shuttle) :
40   AliPreprocessor("PMD", shuttle)
41 {
42   // constructor
43   AddRunType("PHYSICS");
44   AddRunType("PEDESTAL");
45 }
46
47 //______________________________________________________________________________________________
48 AliPMDPreprocessor::~AliPMDPreprocessor()
49 {
50   // destructor
51 }
52
53 //______________________________________________________________________________________________
54 void AliPMDPreprocessor::Initialize(Int_t run, UInt_t startTime,
55         UInt_t endTime)
56 {
57   // Creates AliPMDDataDAQ object
58
59     AliPreprocessor::Initialize(run, startTime, endTime);
60
61     AliInfo(Form("\n\tRun %d \n\tStartTime %s \n\tEndTime %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 RunType = GetRunType();
75     Log(Form("RunType %s",RunType.Data()));
76     if (RunType !="PHYSICS" || RunType != "PEDESTAL") {
77         return kFALSE;
78     }
79
80     return kTRUE;
81 }
82
83
84 //______________________________________________________________________________________________
85 UInt_t AliPMDPreprocessor::Process(TMap* pdaqAliasMap)
86 {
87     
88     if(!pdaqAliasMap) return 1;
89     TString runType = GetRunType();
90     if(runType == "PEDESTAL"){
91         AliPMDPedestal *pedestal = new AliPMDPedestal();
92         
93         TList* filesources = GetFileSources(kDAQ, "PMD_PED.root");
94         
95         if(!filesources) {
96             Log(Form("No sources found for PMD_PED.root!"));
97             return 1;
98         }
99         
100         AliInfo("Here's the list of sources for PMD_PED.root");
101         filesources->Print();
102         
103         TIter iter(filesources);
104         TObjString* source;
105         UInt_t result = 0;
106         TString filename;
107         while((source=dynamic_cast<TObjString*> (iter.Next()))){
108             filename = GetFile(kDAQ, "PMD_PED.root", source->GetName());
109             if(filename.Length() == 0) {
110                 Log(Form("Error retrieving file from source %s failed!", source->GetName()));
111                 delete filesources;
112                 return 1;
113             }
114             
115             Log(Form("File with id PMD_PED.root got from %s", source->GetName()));
116
117             Int_t det, sm, row, col;
118             Float_t mean, rms;
119
120             TFile *f= new TFile(filename.Data());
121             if(!f || !f->IsOpen()) 
122             {
123                 Log(Form("Error opening file with Id PMD_PED.root from source %s!", source->GetName()));
124                 return 1;
125             } 
126             TTree *tree = dynamic_cast<TTree *> (f->Get("ped"));
127             if (!tree) 
128             {
129                 Log("Could not find object \"ped\" in PED file!");
130                 return 1;
131             }
132             
133             tree->SetBranchAddress("det",  &det);
134             tree->SetBranchAddress("sm",   &sm);
135             tree->SetBranchAddress("row",  &row);
136             tree->SetBranchAddress("col",  &col);
137             tree->SetBranchAddress("mean", &mean);
138             tree->SetBranchAddress("rms",  &rms);
139
140             Int_t nEntries = (Int_t) tree->GetEntries();
141             for(Int_t i = 0; i < nEntries; i++)
142             {
143                 tree->GetEntry(i);
144                 pedestal->SetPedMeanRms(det,sm,row,col,mean,rms);
145             }
146             f->Close();
147             delete f;
148         }
149         AliCDBMetaData metaData;
150         metaData.SetBeamPeriod(0);
151         metaData.SetComment("test PMD preprocessor");
152         
153         result = Store("Calib","Ped", pedestal, &metaData,0,kTRUE);
154         delete pedestal;
155         if(result==0)
156         {
157             Log("Error storing");                        
158             return 1;
159         }
160         else
161         {
162             return 0;
163         }
164         
165     }else if (runType == "PHYSICS"){
166         
167         AliPMDCalibData *calibda = new AliPMDCalibData();
168         
169         TList* filesources = GetFileSources(kDAQ, "PMDGAINS.root");
170         
171         if(!filesources) {
172             Log(Form("No sources found for PMDGAINS.root!"));
173             return 1;
174         }
175         
176         AliInfo("Here's the list of sources for PMDGAINS.root");
177         filesources->Print();
178         
179         TIter iter(filesources);
180         TObjString* source;
181         UInt_t result = 0;
182         TString filename;
183         while((source=dynamic_cast<TObjString*> (iter.Next()))){
184             filename = GetFile(kDAQ, "PMDGAINS.root", source->GetName());
185             if(filename.Length() == 0) {
186                 Log(Form("Error retrieving file from source %s failed!", source->GetName()));
187                 delete filesources;
188                 return 1;
189             }
190             
191             Log(Form("File with id PMDGAINS.root got from %s", source->GetName()));
192
193             Int_t det, sm, row, col;
194             Float_t gain;
195
196             TFile *f1= new TFile(filename.Data());
197             if(!f1 || !f1->IsOpen()) 
198             {
199                 Log(Form("Error opening file with Id PMDGAINS.root from source %s!", source->GetName()));
200                 return 1;
201             } 
202             TTree *tree = dynamic_cast<TTree *> (f1->Get("ic"));
203             if (!tree) 
204             {
205                 Log("Could not find object \"ic\" in DAQ file!");
206                 return 1;
207             }
208             
209             tree->SetBranchAddress("det",  &det);
210             tree->SetBranchAddress("sm",   &sm);
211             tree->SetBranchAddress("row",  &row);
212             tree->SetBranchAddress("col",  &col);
213             tree->SetBranchAddress("gain", &gain);
214
215
216
217             Int_t nEntries = (Int_t) tree->GetEntries();
218             for(Int_t i = 0; i < nEntries; i++)
219             {
220                 tree->GetEntry(i);
221
222                 //if(DET>1 || SM>23 || ROW>95 || COL>95) {
223                 //  printf("Error! gain[%d,%d,%d,%d] = %f\n",
224                 //   DET,SM,ROW,COL,GAIN);
225                 //  continue;
226                 //}
227
228                 calibda->SetGainFact(det,sm,row,col,gain);
229             }
230             f1->Close();
231             delete f1;
232         }
233
234         AliCDBMetaData metaData;
235         metaData.SetBeamPeriod(0);
236         metaData.SetComment("test PMD preprocessor");
237         result = Store("Calib","Gain", calibda, &metaData);
238         delete calibda;
239         if(result==0)
240         {
241             Log("Error storing");                        
242             return 1;
243         }
244         else
245         {
246             return 0;
247         }
248         
249   // Store DCS data for reference
250   AliCDBMetaData metadata;
251   metadata.SetComment("DCS data for PMD");
252   Bool_t resStore = kFALSE;
253   resStore = StoreReferenceData("DCS","Data",pdaqAliasMap,&metadata);
254         if(resStore==0)
255         {
256             Log("Error storing");                        
257             return 1;
258         }
259         else
260         {
261             return 0;
262         }
263
264     }
265     
266     return 2;
267 }
268
269