]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EMCAL/AliEMCALCalibHistoProducer.cxx
move minuit initialization to the unfolding method since it leaks and it is not used...
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALCalibHistoProducer.cxx
CommitLineData
9e788b10 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 *
ec955173 19 * $Log$
20 * Revision 1.1 2006/12/07 16:32:16 gustavo
21 * First shuttle code, online calibration histograms producer, EMCAL preprocessor
22 *
9e788b10 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
9e788b10 36#include "TH1.h"
37#include "TFile.h"
38#include "TProfile.h"
ec955173 39
40
41#include "AliLog.h"
9e788b10 42#include "AliRawReader.h"
32cd4c24 43#include "AliCaloRawStreamV3.h"
ec955173 44#include "AliEMCALCalibHistoProducer.h"
9e788b10 45
46ClassImp(AliEMCALCalibHistoProducer)
47
48//-----------------------------------------------------------------------------
49AliEMCALCalibHistoProducer::AliEMCALCalibHistoProducer(AliRawReader* rawReader) :
ec955173 50 TObject(),fRawReader(rawReader),fHistoFile(0),fHistoFileName("calibEmcHisto.root"),
9a090ccd 51 fUpdatingRate(100), fNSuperModules(12), fNCellsEta (48),
9e788b10 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
ec955173 64}
65//-----------------------------------------------------------------------------
66AliEMCALCalibHistoProducer::AliEMCALCalibHistoProducer() :
67 fRawReader(0x0),fHistoFile(0),fHistoFileName(""),
9a090ccd 68 fUpdatingRate(0), fNSuperModules(12), fNCellsEta (48),
ec955173 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
9e788b10 81}
82
83//-----------------------------------------------------------------------------
ec955173 84AliEMCALCalibHistoProducer::AliEMCALCalibHistoProducer(const AliEMCALCalibHistoProducer & copy) :
85 TObject(copy),fRawReader((AliRawReader*)copy. fRawReader->Clone()),
86 fHistoFile((TFile*)copy.fHistoFile->Clone()),fHistoFileName(copy.fHistoFileName),
9a090ccd 87 fUpdatingRate(copy.fUpdatingRate),
ec955173 88 fNSuperModules(copy.fNSuperModules), fNCellsEta (copy.fNCellsEta),
89 fNCellsPhi(copy.fNCellsPhi), fNCellsPhiHalfSM(copy.fNCellsPhiHalfSM)
9e788b10 90{
ec955173 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
9e788b10 101}
102
103//-----------------------------------------------------------------------------
104AliEMCALCalibHistoProducer::~AliEMCALCalibHistoProducer()
105{
106 // Destructor
107 if(fHistoFile) {
108 fHistoFile->Close();
109 delete fHistoFile;
110 }
111}
ec955173 112
113//------------------------------------------------------------------------------
114//
115AliEMCALCalibHistoProducer& 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;
ec955173 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}
9e788b10 143//-----------------------------------------------------------------------------
144void 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//-----------------------------------------------------------------------------
176void 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;
32cd4c24 186 char hname[128];
9e788b10 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
32cd4c24 192 AliCaloRawStreamV3 in(fRawReader,"EMCAL");
9e788b10 193
194 // Read raw data event by event
195
196 while (fRawReader->NextEvent()) {
197 runNum = fRawReader->GetRunNumber();
9e788b10 198
32cd4c24 199 while (in.NextDDL()) {
200 while (in.NextChannel()) {
9e788b10 201
32cd4c24 202 if(fSMInstalled[in.GetModule()]==kFALSE) continue;
9e788b10 203
56613d93 204 // loop over samples
205 int nsamples = 0;
206 Int_t maxSample = 0;
207 while (in.NextBunch()) {
208 const UShort_t *sig = in.GetSignals();
209 nsamples += in.GetBunchLength();
210 for (Int_t i = 0; i < in.GetBunchLength(); i++) {
211 if (sig[i] > maxSample) maxSample = sig[i];
212 }
213 } // bunches
214
215 if (nsamples > 0) { // this check is needed for when we have zero-supp. on, but not sparse readout
216
32cd4c24 217 // indices
9e788b10 218 Int_t mod = in.GetModule();
219 Int_t col = in.GetColumn();
220 Int_t row = in.GetRow();
221 Int_t evtbin = iEvent/nProfFreq;
32cd4c24 222 Bool_t HighGainFlag = in.IsHighGain();
223
9e788b10 224 //Check if histogram/profile already exist, if not create it.
225 if(!fAmpHisto[mod][col][row]) {
226 sprintf(hname,"mod%dcol%drow%d",mod,col,row);
227 fAmpHisto[mod][col][row] = new TH1F(hname,hname,1024,-0.5,1023.);
228 }
229 if(!fAmpProf[mod]) {
230 sprintf(hname,"mod%d",mod);
231 fAmpProf[mod] = new TProfile(hname,hname,nEvtBins,0.,nEvtBins);
232 }
56613d93 233
9e788b10 234 //Fill histogram/profile
32cd4c24 235 if(HighGainFlag) {
236 fAmpHisto[mod][col][row]->Fill(maxSample);
237 fAmpProf[mod]->Fill(evtbin, maxSample);
9e788b10 238 }
56613d93 239
240 } // nsamples>0 check, some data found for this channel; not only trailer/header
241
32cd4c24 242 } // channels
243 } // DDL's
9e788b10 244
245 // update histograms in local file every 100th event
246 if(iEvent%fUpdatingRate == 0) {
247 AliInfo(Form("Updating histo file, event %d, run %d\n",iEvent,runNum));
248 UpdateHistoFile();
249 }
250 iEvent++;
251 }
252
253 UpdateHistoFile();
254 AliInfo(Form("%d events of run %d processed.",iEvent,runNum));
255}
256
257//-----------------------------------------------------------------------------
258void AliEMCALCalibHistoProducer::UpdateHistoFile()
259{
260 // Write histograms to file
261
262 if(!fHistoFile) return;
263 if(!fHistoFile->IsOpen()) return;
264
265 TH1F* hist=0;
266 TProfile* prof =0;
267
268 Int_t nRow = fNCellsPhi ;
269 for(Int_t supermodule=0; supermodule<fNSuperModules; supermodule++) {
270
271 prof = fAmpProf[supermodule];
272 if(prof) prof->Write(prof->GetName(),TObject::kWriteDelete);
273
274 if(supermodule > 10) nRow = fNCellsPhiHalfSM ; //Supermodules 11 and 12 are half supermodules
275 for(Int_t column=0; column<fNCellsEta; column++) {
276 for(Int_t row=0; row<nRow; row++) {
277 hist = fAmpHisto[supermodule][column][row];
278 if(hist) hist->Write(hist->GetName(),TObject::kWriteDelete);
279 }
280 }
281 }
282
283}