]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSCalibHistoProducer.cxx
Remove declaration of unimplemented method
[u/mrichter/AliRoot.git] / PHOS / AliPHOSCalibHistoProducer.cxx
CommitLineData
1ab07e55 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 "AliPHOSRawStream.h"
34
35ClassImp(AliPHOSCalibHistoProducer)
36
37//-----------------------------------------------------------------------------
38AliPHOSCalibHistoProducer::AliPHOSCalibHistoProducer(AliRawReader* rawReader) :
39 fRawReader(rawReader),fHistoFile(0),fUpdatingRate(100)
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//-----------------------------------------------------------------------------
63AliPHOSCalibHistoProducer::AliPHOSCalibHistoProducer() :
64 fRawReader(0),fHistoFile(0),fUpdatingRate(0)
65{
66 // Default constructor
67}
68
69//-----------------------------------------------------------------------------
70AliPHOSCalibHistoProducer::~AliPHOSCalibHistoProducer()
71{
72 // Destructor
73 if(fHistoFile) delete fHistoFile;
74}
75
76//-----------------------------------------------------------------------------
77void 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 AliPHOSRawStream in(fRawReader);
89 in.SetOldRCUFormat(kTRUE);
90
91 // Read raw data event by event
92
93 while (fRawReader->NextEvent()) {
94 runNum = fRawReader->GetRunNumber();
95
96 while ( in.Next() ) {
97
98 if(!gHighGain) gHighGain = new TH1F("gHighGain","High gain",
99 in.GetTimeLength(),0,in.GetTimeLength());
100 if(!gLowGain) gLowGain = new TH1F("gLowGain","Low gain",
101 in.GetTimeLength(),0,in.GetTimeLength());
102
103 Bool_t lowGainFlag = in.IsLowGain();
104
105 if(lowGainFlag)
106 gLowGain->SetBinContent(in.GetTimeLength()-iBin-1,in.GetSignal());
107 else {
108 gHighGain->SetBinContent(in.GetTimeLength()-iBin-1,in.GetSignal());
109 }
110
111 iBin++;
112
113 if(iBin==in.GetTimeLength()) {
114 iBin=0;
115
116 Double_t energy;
117
118 if(!lowGainFlag) {
119 energy = gHighGain->GetMaximum(); // no pedestal subtraction!
120 }
121 else {
122 energy = gLowGain->GetMaximum(); // no pedestal subtraction!
123 }
124
125 Int_t mod = in.GetModule();
126 Int_t col = in.GetColumn();
127 Int_t row = in.GetRow();
128
129 if(fAmpHisto[mod][col][row]) {
130 if(!lowGainFlag) {
131 fAmpHisto[mod][col][row]->Fill(energy);
132 }
133 }
134 else {
135 char hname[128];
136 sprintf(hname,"mod%dcol%drow%d",mod,col,row);
137 fAmpHisto[mod][col][row] = new TH1F(hname,hname,100,0.,1000.);
138 }
139
140
141 }
142 }
143
144 // update histograms in local file every 100th event
145 if(iEvent%fUpdatingRate == 0) {
146 AliInfo(Form("Updating histo file, event %d, run %d\n",iEvent,runNum));
147 UpdateHistoFile();
148 }
149 iEvent++;
150 }
151
152 UpdateHistoFile();
153 AliInfo(Form("%d events of run %d processed.",iEvent,runNum));
154}
155
156//-----------------------------------------------------------------------------
157void AliPHOSCalibHistoProducer::UpdateHistoFile()
158{
159 // Write histograms to file
160
161 if(!fHistoFile) return;
162 if(!fHistoFile->IsOpen()) return;
163
164 TH1F* hist=0;
165
166 for(Int_t module=0; module<5; module++) {
167 for(Int_t column=0; column<56; column++) {
168 for(Int_t row=0; row<64; row++) {
169 hist = fAmpHisto[module][column][row];
170 if(hist) hist->Write(hist->GetName(),TObject::kWriteDelete);
171 }
172 }
173 }
174
175}