]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALCalibHistoProducer.cxx
Typo corrected.
[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),fIsOldRCUFormat(kFALSE), 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),fIsOldRCUFormat(kFALSE), 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),fIsOldRCUFormat(copy.fIsOldRCUFormat),
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         fIsOldRCUFormat = copy.fIsOldRCUFormat;
126         fNSuperModules = copy.fNSuperModules;
127         fNCellsEta = copy.fNCellsEta;
128         fNCellsPhi = copy.fNCellsPhi;
129         fNCellsPhiHalfSM = copy.fNCellsPhiHalfSM;
130         
131         fRawReader  = (AliRawReader*)copy. fRawReader->Clone();
132         fHistoFile      = (TFile*)copy.fHistoFile->Clone();
133
134         for(Int_t ism=0; ism<fNSuperModules; ism++) {
135           fAmpProf[ism] = copy. fAmpProf[ism];
136           fSMInstalled[ism]= copy.fSMInstalled[ism];
137           for(Int_t icol=0; icol<fNCellsEta; icol++) 
138             for(Int_t irow=0; irow<fNCellsPhi; irow++) 
139               fAmpHisto[ism][icol][irow]= copy.fAmpHisto[ism][icol][irow];
140         }
141
142         return (*this);
143 }
144 //-----------------------------------------------------------------------------
145 void AliEMCALCalibHistoProducer::Init()
146 {
147   // initializes input data stream supplied by rawReader
148   // Checks existence of histograms which might have been left
149   // from the previous runs to continue their filling
150   fHistoFile =  new TFile(fHistoFileName,"update");
151   char hname[128];
152   Int_t nRow =  fNCellsPhi ;
153
154   for(Int_t supermodule=0; supermodule<fNSuperModules; supermodule++) {
155     //Check installed supermodules
156     if(fSMInstalled[supermodule]==kFALSE) continue;
157     //Check created profiles
158     sprintf(hname,"mod%d",supermodule);
159     TProfile* prof = (TProfile*)fHistoFile->Get(hname);
160     if(prof)
161       fAmpProf[supermodule]=prof;
162     
163     //Check created histograms
164     if(supermodule > 10) nRow = fNCellsPhiHalfSM ; //Supermodules 11 and 12 are half supermodules
165     for(Int_t column=0; column<fNCellsEta; column++) {
166       for(Int_t row=0; row<nRow; row++) {
167         sprintf(hname,"mod%dcol%drow%d",supermodule,column,row);
168         TH1F* hist = (TH1F*)fHistoFile->Get(hname);
169         if(hist) 
170           fAmpHisto[supermodule][column][row]=hist;
171       }
172     }
173   }
174   
175 }
176 //-----------------------------------------------------------------------------
177 void AliEMCALCalibHistoProducer::Run()
178 {
179   // Reads raw data stream and fills amplitude histograms
180   // The histograms are written to file every fUpdatingRate events
181   //Also fills profiles to study the stability of supermodules during runs.
182
183   Init();
184   
185 //   TH1F* gHighGain = 0;
186 //   TH1F* gLowGain = 0;
187   Int_t iBin = 0;
188   Int_t iEvent = 0;
189   Int_t runNum = 0;
190   Int_t nProfFreq = 1000; //Number of events with which a bin of the TProfile if filled
191   Int_t nEvtBins = 1000; //Total number of the profile survey bins.
192
193   AliCaloRawStream in(fRawReader,"EMCAL");
194   if(fIsOldRCUFormat)
195     in.SetOldRCUFormat(kTRUE);
196
197   // Read raw data event by event
198
199   while (fRawReader->NextEvent()) {
200     runNum = fRawReader->GetRunNumber();
201     Float_t energy = 0;
202      while ( in.Next() ) { 
203
204       if(fSMInstalled[in.GetModule()]==kFALSE) continue;
205        
206       if (in.GetSignal() > energy) {
207         energy = (Double_t) in.GetSignal();
208       }    
209
210       iBin++;
211
212       if(iBin==in.GetTimeLength()) {
213         iBin=0;
214             
215         Int_t mod = in.GetModule();
216         Int_t col = in.GetColumn();
217         Int_t row = in.GetRow();
218         Int_t evtbin = iEvent/nProfFreq;
219         char hname[128];
220
221         //Check if histogram/profile already exist, if not create it.
222         if(!fAmpHisto[mod][col][row]) {
223           sprintf(hname,"mod%dcol%drow%d",mod,col,row);
224           fAmpHisto[mod][col][row] = new TH1F(hname,hname,1024,-0.5,1023.);
225         }
226         if(!fAmpProf[mod]) {
227           sprintf(hname,"mod%d",mod);
228           fAmpProf[mod] = new TProfile(hname,hname,nEvtBins,0.,nEvtBins);
229         }
230         //Fill histogram/profile 
231         Bool_t lowGainFlag = in.IsLowGain();
232         if(!lowGainFlag) {
233           fAmpHisto[mod][col][row]->Fill(energy);
234           fAmpProf[mod]->Fill(evtbin, energy);
235         }
236       }
237     }
238
239     // update histograms in local file every 100th event
240     if(iEvent%fUpdatingRate == 0) {
241       AliInfo(Form("Updating histo file, event %d, run %d\n",iEvent,runNum));
242       UpdateHistoFile();
243     } 
244     iEvent++;
245   }
246
247   UpdateHistoFile(); 
248   AliInfo(Form("%d events of run %d processed.",iEvent,runNum));
249 }
250
251 //-----------------------------------------------------------------------------
252 void AliEMCALCalibHistoProducer::UpdateHistoFile()
253 {
254   // Write histograms to file
255
256   if(!fHistoFile) return;
257   if(!fHistoFile->IsOpen()) return;
258
259   TH1F* hist=0;
260   TProfile* prof =0;
261  
262   Int_t nRow =  fNCellsPhi ;
263   for(Int_t supermodule=0; supermodule<fNSuperModules; supermodule++) {
264     
265     prof = fAmpProf[supermodule]; 
266     if(prof) prof->Write(prof->GetName(),TObject::kWriteDelete);
267     
268     if(supermodule > 10)  nRow = fNCellsPhiHalfSM ; //Supermodules 11 and 12 are half supermodules
269     for(Int_t column=0; column<fNCellsEta; column++) {
270       for(Int_t row=0; row<nRow; row++) {
271         hist = fAmpHisto[supermodule][column][row]; 
272         if(hist) hist->Write(hist->GetName(),TObject::kWriteDelete);
273       }
274     }
275   }
276   
277 }