]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALCalibHistoProducer.cxx
using non-uniform field option for tracking (Jacek)
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALCalibHistoProducer.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 /* $Id$ */
17 /* History of cvs commits:
18  *
19  * $Log$
20  * Revision 1.1  2006/12/07 16:32:16  gustavo
21  * First shuttle code, online calibration histograms producer, EMCAL preprocessor
22  * 
23  *
24 */
25 ///////////////////////////////////////////////////////////////////////////////
26 // Class AliEMCALCalibHistoProducer accumulating histograms
27 // with amplitudes per EMCAL channel
28 // It is intended to run at DAQ computers (LDC, GDC, HLT or MOOD)
29 // and it fills the histograms with amplitudes per channel.
30 // Usage example see in EMCAL/macros/Shuttle/AliEMCALCalibHistoProducer.C
31 //
32 // Author: Boris Polichtchouk, 4 October 2006
33 // Adapted for EMCAL by Gustavo Conesa Balbastre, October 2006
34 ///////////////////////////////////////////////////////////////////////////////
35
36 #include "TH1.h"
37 #include "TFile.h"
38 #include "TProfile.h"
39
40
41 #include "AliLog.h"
42 #include "AliRawReader.h"
43 #include "AliCaloRawStream.h"
44 #include "AliEMCALCalibHistoProducer.h"
45
46 ClassImp(AliEMCALCalibHistoProducer)
47
48 //-----------------------------------------------------------------------------
49 AliEMCALCalibHistoProducer::AliEMCALCalibHistoProducer(AliRawReader* rawReader) : 
50   TObject(),fRawReader(rawReader),fHistoFile(0),fHistoFileName("calibEmcHisto.root"),
51   fUpdatingRate(100), fNSuperModules(12),  fNCellsEta (48),   
52   fNCellsPhi(24),  fNCellsPhiHalfSM(12)
53 {
54   // Constructor
55
56   for(Int_t ism=0; ism<fNSuperModules; ism++) {
57     fAmpProf[ism] = 0;
58     fSMInstalled[ism]=kTRUE;
59     for(Int_t icol=0; icol<fNCellsEta; icol++) 
60       for(Int_t irow=0; irow<fNCellsPhi; irow++) 
61           fAmpHisto[ism][icol][irow]=0;
62   }
63
64 }
65 //-----------------------------------------------------------------------------
66 AliEMCALCalibHistoProducer::AliEMCALCalibHistoProducer() : 
67   fRawReader(0x0),fHistoFile(0),fHistoFileName(""),
68   fUpdatingRate(0), fNSuperModules(12),  fNCellsEta (48),   
69   fNCellsPhi(24),  fNCellsPhiHalfSM(12)
70 {
71   // default Constructor
72
73   for(Int_t ism=0; ism<fNSuperModules; ism++) {
74     fAmpProf[ism] = 0;
75     fSMInstalled[ism]=kTRUE;
76     for(Int_t icol=0; icol<fNCellsEta; icol++) 
77       for(Int_t irow=0; irow<fNCellsPhi; irow++) 
78           fAmpHisto[ism][icol][irow]=0;
79   }
80
81 }
82
83 //-----------------------------------------------------------------------------
84 AliEMCALCalibHistoProducer::AliEMCALCalibHistoProducer(const AliEMCALCalibHistoProducer & copy) :
85   TObject(copy),fRawReader((AliRawReader*)copy. fRawReader->Clone()),
86   fHistoFile((TFile*)copy.fHistoFile->Clone()),fHistoFileName(copy.fHistoFileName),
87   fUpdatingRate(copy.fUpdatingRate),
88   fNSuperModules(copy.fNSuperModules), fNCellsEta (copy.fNCellsEta), 
89   fNCellsPhi(copy.fNCellsPhi), fNCellsPhiHalfSM(copy.fNCellsPhiHalfSM)
90 {
91   //copy constructor
92
93  for(Int_t ism=0; ism<fNSuperModules; ism++) {
94     fAmpProf[ism] = copy. fAmpProf[ism];
95     fSMInstalled[ism]= copy.fSMInstalled[ism];
96     for(Int_t icol=0; icol<fNCellsEta; icol++) 
97       for(Int_t irow=0; irow<fNCellsPhi; irow++) 
98           fAmpHisto[ism][icol][irow]= copy.fAmpHisto[ism][icol][irow];
99   }
100
101 }
102
103 //-----------------------------------------------------------------------------
104 AliEMCALCalibHistoProducer::~AliEMCALCalibHistoProducer()
105 {
106   // Destructor
107   if(fHistoFile) {
108     fHistoFile->Close();
109     delete fHistoFile;
110   }
111 }
112
113 //------------------------------------------------------------------------------
114 //
115 AliEMCALCalibHistoProducer& AliEMCALCalibHistoProducer::operator=(const AliEMCALCalibHistoProducer& copy)
116 {
117         //
118         // Assignment operator.
119         // Besides copying all parameters, duplicates all collections.  
120         //
121                 if (&copy == this) return *this;
122         TObject::operator=(copy);
123         fHistoFileName = copy.fHistoFileName;
124         fUpdatingRate = copy.fUpdatingRate;
125         fNSuperModules = copy.fNSuperModules;
126         fNCellsEta = copy.fNCellsEta;
127         fNCellsPhi = copy.fNCellsPhi;
128         fNCellsPhiHalfSM = copy.fNCellsPhiHalfSM;
129         
130         fRawReader  = (AliRawReader*)copy. fRawReader->Clone();
131         fHistoFile      = (TFile*)copy.fHistoFile->Clone();
132
133         for(Int_t ism=0; ism<fNSuperModules; ism++) {
134           fAmpProf[ism] = copy. fAmpProf[ism];
135           fSMInstalled[ism]= copy.fSMInstalled[ism];
136           for(Int_t icol=0; icol<fNCellsEta; icol++) 
137             for(Int_t irow=0; irow<fNCellsPhi; irow++) 
138               fAmpHisto[ism][icol][irow]= copy.fAmpHisto[ism][icol][irow];
139         }
140
141         return (*this);
142 }
143 //-----------------------------------------------------------------------------
144 void AliEMCALCalibHistoProducer::Init()
145 {
146   // initializes input data stream supplied by rawReader
147   // Checks existence of histograms which might have been left
148   // from the previous runs to continue their filling
149   fHistoFile =  new TFile(fHistoFileName,"update");
150   char hname[128];
151   Int_t nRow =  fNCellsPhi ;
152
153   for(Int_t supermodule=0; supermodule<fNSuperModules; supermodule++) {
154     //Check installed supermodules
155     if(fSMInstalled[supermodule]==kFALSE) continue;
156     //Check created profiles
157     sprintf(hname,"mod%d",supermodule);
158     TProfile* prof = (TProfile*)fHistoFile->Get(hname);
159     if(prof)
160       fAmpProf[supermodule]=prof;
161     
162     //Check created histograms
163     if(supermodule > 10) nRow = fNCellsPhiHalfSM ; //Supermodules 11 and 12 are half supermodules
164     for(Int_t column=0; column<fNCellsEta; column++) {
165       for(Int_t row=0; row<nRow; row++) {
166         sprintf(hname,"mod%dcol%drow%d",supermodule,column,row);
167         TH1F* hist = (TH1F*)fHistoFile->Get(hname);
168         if(hist) 
169           fAmpHisto[supermodule][column][row]=hist;
170       }
171     }
172   }
173   
174 }
175 //-----------------------------------------------------------------------------
176 void AliEMCALCalibHistoProducer::Run()
177 {
178   // Reads raw data stream and fills amplitude histograms
179   // The histograms are written to file every fUpdatingRate events
180   //Also fills profiles to study the stability of supermodules during runs.
181
182   Init();
183   
184 //   TH1F* gHighGain = 0;
185 //   TH1F* gLowGain = 0;
186   Int_t iBin = 0;
187   Int_t iEvent = 0;
188   Int_t runNum = 0;
189   Int_t nProfFreq = 1000; //Number of events with which a bin of the TProfile if filled
190   Int_t nEvtBins = 1000; //Total number of the profile survey bins.
191
192   AliCaloRawStream in(fRawReader,"EMCAL");
193
194   // Read raw data event by event
195
196   while (fRawReader->NextEvent()) {
197     runNum = fRawReader->GetRunNumber();
198     Float_t energy = 0;
199      while ( in.Next() ) { 
200
201       if(fSMInstalled[in.GetModule()]==kFALSE) continue;
202        
203       if (in.GetSignal() > energy) {
204         energy = (Double_t) in.GetSignal();
205       }    
206
207       iBin++;
208
209       if(iBin==in.GetTimeLength()) {
210         iBin=0;
211             
212         Int_t mod = in.GetModule();
213         Int_t col = in.GetColumn();
214         Int_t row = in.GetRow();
215         Int_t evtbin = iEvent/nProfFreq;
216         char hname[128];
217
218         //Check if histogram/profile already exist, if not create it.
219         if(!fAmpHisto[mod][col][row]) {
220           sprintf(hname,"mod%dcol%drow%d",mod,col,row);
221           fAmpHisto[mod][col][row] = new TH1F(hname,hname,1024,-0.5,1023.);
222         }
223         if(!fAmpProf[mod]) {
224           sprintf(hname,"mod%d",mod);
225           fAmpProf[mod] = new TProfile(hname,hname,nEvtBins,0.,nEvtBins);
226         }
227         //Fill histogram/profile 
228         Bool_t lowGainFlag = in.IsLowGain();
229         if(!lowGainFlag) {
230           fAmpHisto[mod][col][row]->Fill(energy);
231           fAmpProf[mod]->Fill(evtbin, energy);
232         }
233       }
234     }
235
236     // update histograms in local file every 100th event
237     if(iEvent%fUpdatingRate == 0) {
238       AliInfo(Form("Updating histo file, event %d, run %d\n",iEvent,runNum));
239       UpdateHistoFile();
240     } 
241     iEvent++;
242   }
243
244   UpdateHistoFile(); 
245   AliInfo(Form("%d events of run %d processed.",iEvent,runNum));
246 }
247
248 //-----------------------------------------------------------------------------
249 void AliEMCALCalibHistoProducer::UpdateHistoFile()
250 {
251   // Write histograms to file
252
253   if(!fHistoFile) return;
254   if(!fHistoFile->IsOpen()) return;
255
256   TH1F* hist=0;
257   TProfile* prof =0;
258  
259   Int_t nRow =  fNCellsPhi ;
260   for(Int_t supermodule=0; supermodule<fNSuperModules; supermodule++) {
261     
262     prof = fAmpProf[supermodule]; 
263     if(prof) prof->Write(prof->GetName(),TObject::kWriteDelete);
264     
265     if(supermodule > 10)  nRow = fNCellsPhiHalfSM ; //Supermodules 11 and 12 are half supermodules
266     for(Int_t column=0; column<fNCellsEta; column++) {
267       for(Int_t row=0; row<nRow; row++) {
268         hist = fAmpHisto[supermodule][column][row]; 
269         if(hist) hist->Write(hist->GetName(),TObject::kWriteDelete);
270       }
271     }
272   }
273   
274 }