]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSCalibHistoProducer.cxx
Forward declaration
[u/mrichter/AliRoot.git] / PHOS / AliPHOSCalibHistoProducer.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
18 ///////////////////////////////////////////////////////////////////////////////
19 // Class AliPHOSCalibHistoProducer accumulating histograms
20 // with amplitudes per PHOS channel
21 // It is intended to run at DAQ computers (LDC, GDC, HLT or MOOD)
22 // and it fills the histograms with amplitudes per channel.
23 // Usage example see in PHOS/macros/Shuttle/AliPHOSCalibHistoProducer.C
24 //
25 // Author: Boris Polichtchouk, 4 October 2006
26 ///////////////////////////////////////////////////////////////////////////////
27
28 #include "AliLog.h"
29 #include "AliPHOSCalibHistoProducer.h"
30 #include "TH1.h"
31 #include "TFile.h"
32 #include "AliRawReader.h"
33 #include "AliCaloRawStream.h"
34
35 ClassImp(AliPHOSCalibHistoProducer)
36
37 //-----------------------------------------------------------------------------
38 AliPHOSCalibHistoProducer::AliPHOSCalibHistoProducer(AliRawReader* rawReader) : 
39   fRawReader(rawReader),fHistoFile(0),fUpdatingRate(100),fIsOldRCUFormat(kFALSE)
40 {
41   // Constructor: initializes input data stream supplied by rawReader
42   // Checks existence of histograms which might have been left
43   // from the previous runs to continues their filling
44
45   fHistoFile =  new TFile("calibHisto.root","update");
46
47   for(Int_t module=0; module<5; module++) {
48     for(Int_t column=0; column<56; column++) {
49       for(Int_t row=0; row<64; row++) {
50         char hname[128];
51         sprintf(hname,"mod%dcol%drow%d",module,column,row);
52         TH1F* hist = (TH1F*)fHistoFile->Get(hname);
53         if(hist) 
54           fAmpHisto[module][column][row]=hist;
55         else
56           fAmpHisto[module][column][row]=0;
57       }
58     }
59   }
60 }
61
62 //-----------------------------------------------------------------------------
63 AliPHOSCalibHistoProducer::AliPHOSCalibHistoProducer() :
64   fRawReader(0),fHistoFile(0),fUpdatingRate(0),fIsOldRCUFormat(kFALSE)
65 {
66   // Default constructor
67 }
68
69 //-----------------------------------------------------------------------------
70 AliPHOSCalibHistoProducer::~AliPHOSCalibHistoProducer()
71 {
72   // Destructor
73   if(fHistoFile) delete fHistoFile;
74 }
75
76 //-----------------------------------------------------------------------------
77 void AliPHOSCalibHistoProducer::Run()
78 {
79   // Reads raw data stream and fills amplitude histograms
80   // The histograms are written to file every fUpdatingRate events
81
82   TH1F* gHighGain = 0;
83   TH1F* gLowGain = 0;
84   Int_t iBin = 0;
85   Int_t iEvent = 0;
86   Int_t runNum = 0;
87
88   AliCaloRawStream in(fRawReader,"PHOS");
89   if(fIsOldRCUFormat)
90     in.SetOldRCUFormat(kTRUE);
91
92   // Read raw data event by event
93
94   while (fRawReader->NextEvent()) {
95     runNum = fRawReader->GetRunNumber();
96
97     while ( in.Next() ) { 
98        
99       if(!gHighGain) gHighGain = new TH1F("gHighGain","High gain",
100                                           in.GetTimeLength(),0,in.GetTimeLength());
101       else
102         if(gHighGain->GetNbinsX() != in.GetTimeLength()) {
103           delete gHighGain;
104           gHighGain = new TH1F("gHighGain","High gain",in.GetTimeLength(),0,in.GetTimeLength());
105         }
106
107       if(!gLowGain)  gLowGain = new TH1F("gLowGain","Low gain",
108                                          in.GetTimeLength(),0,in.GetTimeLength());
109       else
110         if(gLowGain->GetNbinsX() != in.GetTimeLength()) {
111           delete gLowGain;
112           gLowGain = new TH1F("gLowGain","Low gain",in.GetTimeLength(),0,in.GetTimeLength());
113         }
114
115       Bool_t lowGainFlag = in.IsLowGain();
116       
117       if(lowGainFlag) 
118         gLowGain->SetBinContent(in.GetTimeLength()-iBin-1,in.GetSignal());
119       else {
120         gHighGain->SetBinContent(in.GetTimeLength()-iBin-1,in.GetSignal());
121       }
122
123       iBin++;
124
125       if(iBin==in.GetTimeLength()) {
126         iBin=0;
127
128         Double_t energy;
129
130         if(!lowGainFlag) {
131           energy = gHighGain->GetMaximum(); // no pedestal subtraction!
132         }
133         else {
134           energy = gLowGain->GetMaximum(); // no pedestal subtraction!
135         }
136             
137         Int_t mod = in.GetModule();
138         Int_t col = in.GetColumn();
139         Int_t row = in.GetRow();
140
141         if(fAmpHisto[mod][col][row]) {
142           if(!lowGainFlag) {
143             fAmpHisto[mod][col][row]->Fill(energy);
144           }
145         }
146         else {
147           char hname[128];
148           sprintf(hname,"mod%dcol%drow%d",mod,col,row);
149           fAmpHisto[mod][col][row] = new TH1F(hname,hname,100,0.,1000.);
150         }
151
152
153       }
154     }
155
156     // update histograms in local file every 100th event
157     if(iEvent%fUpdatingRate == 0) {
158       AliInfo(Form("Updating histo file, event %d, run %d\n",iEvent,runNum));
159       UpdateHistoFile();
160     } 
161     iEvent++;
162   }
163
164   UpdateHistoFile(); 
165   AliInfo(Form("%d events of run %d processed.",iEvent,runNum));
166 }
167
168 //-----------------------------------------------------------------------------
169 void AliPHOSCalibHistoProducer::UpdateHistoFile()
170 {
171   // Write histograms to file
172
173   if(!fHistoFile) return;
174   if(!fHistoFile->IsOpen()) return;
175
176   TH1F* hist=0;
177
178   for(Int_t module=0; module<5; module++) {
179     for(Int_t column=0; column<56; column++) {
180       for(Int_t row=0; row<64; row++) {
181         hist = fAmpHisto[module][column][row]; 
182         if(hist) hist->Write(hist->GetName(),TObject::kWriteDelete);
183       }
184     }
185   }
186
187 }