]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSCalibHistoProducer.cxx
Updated zero-misalignment object
[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 "AliPHOSDigit.h"
34 #include "AliPHOSRawDecoder.h"
35 #include "AliPHOSRawDigiProducer.h"
36 #include "AliPHOSGeometry.h"
37 #include "TClonesArray.h"
38
39 ClassImp(AliPHOSCalibHistoProducer)
40
41 //-----------------------------------------------------------------------------
42 AliPHOSCalibHistoProducer::AliPHOSCalibHistoProducer() : 
43   fRawReader(0),fHistoFile(0),fUpdatingRate(100),fIsOldRCUFormat(kFALSE),
44   fDigits(0),fEvents(0)
45 {
46   // Constructor: initializes data members
47   // Checks existence of histograms which might have been left
48   // from the previous runs to continues their filling
49
50   fHistoFile =  new TFile("calibHisto.root","update");
51   fDigits = new TClonesArray("AliPHOSDigit",100);
52
53   for(Int_t module=0; module<5; module++) {
54     for(Int_t column=0; column<56; column++) {
55       for(Int_t row=0; row<64; row++) {
56         char hname[128];
57         sprintf(hname,"mod%dcol%drow%d",module,column,row);
58         TH1F* hist = (TH1F*)fHistoFile->Get(hname);
59         if(hist) 
60           fAmpHisto[module][column][row]=hist;
61         else
62           fAmpHisto[module][column][row]=0;
63       }
64     }
65   }
66 }
67
68 //-----------------------------------------------------------------------------
69 AliPHOSCalibHistoProducer::~AliPHOSCalibHistoProducer()
70 {
71   // Destructor
72   
73   UpdateHistoFile();
74
75   if(fHistoFile) delete fHistoFile;
76   if (fDigits) { fDigits->Delete(); delete fDigits; }
77
78 }
79
80 //-----------------------------------------------------------------------------
81 AliPHOSCalibHistoProducer::AliPHOSCalibHistoProducer(const AliPHOSCalibHistoProducer &histoproducer) :
82   TObject(histoproducer),fRawReader(histoproducer.fRawReader),fHistoFile(histoproducer.fHistoFile),
83   fUpdatingRate(histoproducer.fUpdatingRate),fIsOldRCUFormat(histoproducer.fIsOldRCUFormat),
84   fDigits(histoproducer.fDigits),fEvents(histoproducer.fEvents)
85 {
86   //Copy constructor.
87
88   for(Int_t module=0; module<5; module++) {
89     for(Int_t column=0; column<56; column++) {
90       for(Int_t row=0; row<64; row++) {
91         char hname[128];
92         sprintf(hname,"mod%dcol%drow%d",module,column,row);
93         TH1F* hist = (TH1F*)histoproducer.fHistoFile->Get(hname);
94         if(hist) 
95           fAmpHisto[module][column][row]= new TH1F(*hist);
96         else
97           fAmpHisto[module][column][row]=0;
98       }
99     }
100   }
101 }
102
103 //-----------------------------------------------------------------------------
104 AliPHOSCalibHistoProducer& AliPHOSCalibHistoProducer::operator= 
105 (const AliPHOSCalibHistoProducer &histoproducer)
106 {
107   //Assignment operator.
108
109   if(this != &histoproducer) {
110
111     fRawReader = histoproducer.fRawReader;
112     fHistoFile = histoproducer.fHistoFile;
113     fUpdatingRate = histoproducer.fUpdatingRate;
114     fIsOldRCUFormat = histoproducer.fIsOldRCUFormat;
115     fDigits = histoproducer.fDigits;
116     fEvents = histoproducer.fEvents;
117
118     for(Int_t module=0; module<5; module++) {
119       for(Int_t column=0; column<56; column++) {
120         for(Int_t row=0; row<64; row++) {
121           if(fAmpHisto[module][column][row]){
122             delete fAmpHisto[module][column][row];
123             fAmpHisto[module][column][row] = histoproducer.fAmpHisto[module][column][row];
124           }
125           else
126           fAmpHisto[module][column][row] = histoproducer.fAmpHisto[module][column][row];
127         }
128       }
129     }
130
131
132   }
133
134   return *this;
135 }
136 //-----------------------------------------------------------------------------
137 void AliPHOSCalibHistoProducer::Run()
138 {
139   // Reads raw data of current event and fills amplitude histograms
140   // The histograms are written to file every fUpdatingRate events
141
142   Int_t relId[4];
143   AliPHOSDigit* digit;
144   Double_t energy;
145   
146   AliPHOSGeometry* geo = AliPHOSGeometry::GetInstance();
147   if(!geo) geo = AliPHOSGeometry::GetInstance("IHEP");
148   if(!geo) AliFatal(Form("Cannot get PHOS geometry"));
149
150   AliPHOSRawDecoder dc(fRawReader);
151   if(fIsOldRCUFormat)
152     dc.SetOldRCUFormat(kTRUE);
153   
154   AliPHOSRawDigiProducer producer;  
155   producer.MakeDigits(fDigits,&dc);
156     
157 //   printf("% digits.\n",fDigits->GetEntries());
158   for(Int_t iDigit=0; iDigit<fDigits->GetEntries(); iDigit++) {
159     digit = (AliPHOSDigit*)fDigits->At(iDigit);
160     energy = digit->GetEnergy(); // no pedestal subtraction!
161 //     printf("[absId,E]: [%d,%.3f]\n",digit->GetId(),digit->GetEnergy());
162     geo->AbsToRelNumbering(digit->GetId(),relId);           
163     Int_t mod = relId[0]-1;
164     Int_t col = relId[3]-1;
165     Int_t row = relId[2]-1;
166 //     printf("\t(mod,col,row): (%d,%d,%d)\n",mod,col,row);
167
168     if(fAmpHisto[mod][col][row]) {
169       fAmpHisto[mod][col][row]->Fill(energy);
170     }
171     else {
172       char hname[128];
173       sprintf(hname,"mod%dcol%drow%d",mod,col,row);
174       fAmpHisto[mod][col][row] = new TH1F(hname,hname,100,0.,1000.);
175     }
176     
177     // update histograms in local file every 100th event
178     if(fEvents%fUpdatingRate == 0) {
179       AliInfo(Form("Updating histo file, event %d, run %d\n",fEvents,fRawReader->GetRunNumber()));
180       UpdateHistoFile();
181     } 
182   }
183   
184   fDigits->Delete();
185   fEvents++;
186
187 //   }
188
189 //   UpdateHistoFile();
190 //   AliInfo(Form("%d events of run %d processed.",iEvent,runNum));
191 }
192
193 //-----------------------------------------------------------------------------
194 void AliPHOSCalibHistoProducer::UpdateHistoFile()
195 {
196   // Write histograms to file
197
198   if(!fHistoFile) return;
199   if(!fHistoFile->IsOpen()) return;
200
201   TH1F* hist=0;
202
203   for(Int_t module=0; module<5; module++) {
204     for(Int_t column=0; column<56; column++) {
205       for(Int_t row=0; row<64; row++) {
206         hist = fAmpHisto[module][column][row]; 
207         if(hist) hist->Write(hist->GetName(),TObject::kWriteDelete);
208       }
209     }
210   }
211
212 }