]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALCalibHistoProducer.cxx
Reference calibration histogram file, one histogram per channel with accumulated...
[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  *
21 */
22 ///////////////////////////////////////////////////////////////////////////////
23 // Class AliEMCALCalibHistoProducer accumulating histograms
24 // with amplitudes per EMCAL channel
25 // It is intended to run at DAQ computers (LDC, GDC, HLT or MOOD)
26 // and it fills the histograms with amplitudes per channel.
27 // Usage example see in EMCAL/macros/Shuttle/AliEMCALCalibHistoProducer.C
28 //
29 // Author: Boris Polichtchouk, 4 October 2006
30 // Adapted for EMCAL by Gustavo Conesa Balbastre, October 2006
31 ///////////////////////////////////////////////////////////////////////////////
32
33 #include "AliLog.h"
34 #include "AliEMCALCalibHistoProducer.h"
35 #include "TH1.h"
36 #include "TFile.h"
37 #include "TProfile.h"
38 #include "AliRawReader.h"
39 #include "AliCaloRawStream.h"
40
41 ClassImp(AliEMCALCalibHistoProducer)
42
43 //-----------------------------------------------------------------------------
44 AliEMCALCalibHistoProducer::AliEMCALCalibHistoProducer(AliRawReader* rawReader) : 
45   fRawReader(rawReader),fHistoFile(0),fHistoFileName("calibHisto.root"),
46   fUpdatingRate(100),fIsOldRCUFormat(kFALSE), fNSuperModules(12),  fNCellsEta (48),   
47   fNCellsPhi(24),  fNCellsPhiHalfSM(12)
48 {
49   // Constructor
50
51   for(Int_t ism=0; ism<fNSuperModules; ism++) {
52     fAmpProf[ism] = 0;
53     fSMInstalled[ism]=kTRUE;
54     for(Int_t icol=0; icol<fNCellsEta; icol++) 
55       for(Int_t irow=0; irow<fNCellsPhi; irow++) 
56           fAmpHisto[ism][icol][irow]=0;
57   }
58
59 }
60
61 //-----------------------------------------------------------------------------
62 AliEMCALCalibHistoProducer::AliEMCALCalibHistoProducer() :
63   fRawReader(0),fHistoFile(0),fUpdatingRate(0),fIsOldRCUFormat(kFALSE),
64   fNSuperModules(12), fNCellsEta (48), fNCellsPhi(24), fNCellsPhiHalfSM(12)
65 {
66   // Default constructor
67 }
68
69 //-----------------------------------------------------------------------------
70 AliEMCALCalibHistoProducer::~AliEMCALCalibHistoProducer()
71 {
72   // Destructor
73   if(fHistoFile) {
74     fHistoFile->Close();
75     delete fHistoFile;
76   }
77 }
78 //-----------------------------------------------------------------------------
79 void AliEMCALCalibHistoProducer::Init()
80 {
81   // initializes input data stream supplied by rawReader
82   // Checks existence of histograms which might have been left
83   // from the previous runs to continue their filling
84   fHistoFile =  new TFile(fHistoFileName,"update");
85   char hname[128];
86   Int_t nRow =  fNCellsPhi ;
87
88   for(Int_t supermodule=0; supermodule<fNSuperModules; supermodule++) {
89     //Check installed supermodules
90     if(fSMInstalled[supermodule]==kFALSE) continue;
91     //Check created profiles
92     sprintf(hname,"mod%d",supermodule);
93     TProfile* prof = (TProfile*)fHistoFile->Get(hname);
94     if(prof)
95       fAmpProf[supermodule]=prof;
96     
97     //Check created histograms
98     if(supermodule > 10) nRow = fNCellsPhiHalfSM ; //Supermodules 11 and 12 are half supermodules
99     for(Int_t column=0; column<fNCellsEta; column++) {
100       for(Int_t row=0; row<nRow; row++) {
101         sprintf(hname,"mod%dcol%drow%d",supermodule,column,row);
102         TH1F* hist = (TH1F*)fHistoFile->Get(hname);
103         if(hist) 
104           fAmpHisto[supermodule][column][row]=hist;
105       }
106     }
107   }
108   
109 }
110 //-----------------------------------------------------------------------------
111 void AliEMCALCalibHistoProducer::Run()
112 {
113   // Reads raw data stream and fills amplitude histograms
114   // The histograms are written to file every fUpdatingRate events
115   //Also fills profiles to study the stability of supermodules during runs.
116
117   Init();
118   
119 //   TH1F* gHighGain = 0;
120 //   TH1F* gLowGain = 0;
121   Int_t iBin = 0;
122   Int_t iEvent = 0;
123   Int_t runNum = 0;
124   Int_t nProfFreq = 1000; //Number of events with which a bin of the TProfile if filled
125   Int_t nEvtBins = 1000; //Total number of the profile survey bins.
126
127   AliCaloRawStream in(fRawReader,"EMCAL");
128   if(fIsOldRCUFormat)
129     in.SetOldRCUFormat(kTRUE);
130
131   // Read raw data event by event
132
133   while (fRawReader->NextEvent()) {
134     runNum = fRawReader->GetRunNumber();
135     Float_t energy = 0;
136      while ( in.Next() ) { 
137
138       if(fSMInstalled[in.GetModule()]==kFALSE) continue;
139        
140       if (in.GetSignal() > energy) {
141         energy = (Double_t) in.GetSignal();
142       }    
143
144       iBin++;
145
146       if(iBin==in.GetTimeLength()) {
147         iBin=0;
148             
149         Int_t mod = in.GetModule();
150         Int_t col = in.GetColumn();
151         Int_t row = in.GetRow();
152         Int_t evtbin = iEvent/nProfFreq;
153         char hname[128];
154
155         //Check if histogram/profile already exist, if not create it.
156         if(!fAmpHisto[mod][col][row]) {
157           sprintf(hname,"mod%dcol%drow%d",mod,col,row);
158           fAmpHisto[mod][col][row] = new TH1F(hname,hname,1024,-0.5,1023.);
159         }
160         if(!fAmpProf[mod]) {
161           sprintf(hname,"mod%d",mod);
162           fAmpProf[mod] = new TProfile(hname,hname,nEvtBins,0.,nEvtBins);
163         }
164         //Fill histogram/profile 
165         Bool_t lowGainFlag = in.IsLowGain();
166         if(!lowGainFlag) {
167           fAmpHisto[mod][col][row]->Fill(energy);
168           fAmpProf[mod]->Fill(evtbin, energy);
169         }
170       }
171     }
172
173     // update histograms in local file every 100th event
174     if(iEvent%fUpdatingRate == 0) {
175       AliInfo(Form("Updating histo file, event %d, run %d\n",iEvent,runNum));
176       UpdateHistoFile();
177     } 
178     iEvent++;
179   }
180
181   UpdateHistoFile(); 
182   AliInfo(Form("%d events of run %d processed.",iEvent,runNum));
183 }
184
185 //-----------------------------------------------------------------------------
186 void AliEMCALCalibHistoProducer::UpdateHistoFile()
187 {
188   // Write histograms to file
189
190   if(!fHistoFile) return;
191   if(!fHistoFile->IsOpen()) return;
192
193   TH1F* hist=0;
194   TProfile* prof =0;
195  
196   Int_t nRow =  fNCellsPhi ;
197   for(Int_t supermodule=0; supermodule<fNSuperModules; supermodule++) {
198     
199     prof = fAmpProf[supermodule]; 
200     if(prof) prof->Write(prof->GetName(),TObject::kWriteDelete);
201     
202     if(supermodule > 10)  nRow = fNCellsPhiHalfSM ; //Supermodules 11 and 12 are half supermodules
203     for(Int_t column=0; column<fNCellsEta; column++) {
204       for(Int_t row=0; row<nRow; row++) {
205         hist = fAmpHisto[supermodule][column][row]; 
206         if(hist) hist->Write(hist->GetName(),TObject::kWriteDelete);
207       }
208     }
209   }
210   
211 }