]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSCalibHistoProducer.cxx
Runloader is updated when moving to next file (quick fix).
[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(AliRawReader* rawReader) : 
43   fRawReader(rawReader),fHistoFile(0),fUpdatingRate(100),fIsOldRCUFormat(kFALSE)
44 {
45   // Constructor: initializes input data stream supplied by rawReader
46   // Checks existence of histograms which might have been left
47   // from the previous runs to continues their filling
48
49   fHistoFile =  new TFile("calibHisto.root","update");
50
51   for(Int_t module=0; module<5; module++) {
52     for(Int_t column=0; column<56; column++) {
53       for(Int_t row=0; row<64; row++) {
54         char hname[128];
55         sprintf(hname,"mod%dcol%drow%d",module,column,row);
56         TH1F* hist = (TH1F*)fHistoFile->Get(hname);
57         if(hist) 
58           fAmpHisto[module][column][row]=hist;
59         else
60           fAmpHisto[module][column][row]=0;
61       }
62     }
63   }
64 }
65
66 //-----------------------------------------------------------------------------
67 AliPHOSCalibHistoProducer::AliPHOSCalibHistoProducer() :
68   fRawReader(0),fHistoFile(0),fUpdatingRate(0),fIsOldRCUFormat(kFALSE)
69 {
70   // Default constructor
71 }
72
73 //-----------------------------------------------------------------------------
74 AliPHOSCalibHistoProducer::~AliPHOSCalibHistoProducer()
75 {
76   // Destructor
77   if(fHistoFile) delete fHistoFile;
78 }
79
80 //-----------------------------------------------------------------------------
81 void AliPHOSCalibHistoProducer::Run()
82 {
83   // Reads raw data stream and fills amplitude histograms
84   // The histograms are written to file every fUpdatingRate events
85
86   Int_t iEvent = 0;
87   Int_t runNum = -111;
88   Int_t relId[4];
89   AliPHOSDigit* digit;
90   Double_t energy;
91
92   TClonesArray *digits = new TClonesArray("AliPHOSDigit",100);
93   AliPHOSGeometry* geo = AliPHOSGeometry::GetInstance("IHEP");
94
95   AliPHOSRawDecoder dc(fRawReader);
96   if(fIsOldRCUFormat)
97     dc.SetOldRCUFormat(kTRUE);
98   
99   AliPHOSRawDigiProducer producer;
100
101   // Read raw data event by event
102   
103   while (fRawReader->NextEvent()) {
104     runNum = fRawReader->GetRunNumber();
105     producer.MakeDigits(digits,&dc);
106     
107     for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
108       digit = (AliPHOSDigit*)digits->At(iDigit);
109       energy = digit->GetEnergy(); // no pedestal subtraction!
110       geo->AbsToRelNumbering(digit->GetId(),relId);         
111       Int_t mod = relId[0]-1;
112       Int_t col = relId[3]-1;
113       Int_t row = relId[2]-1;
114
115       if(fAmpHisto[mod][col][row]) {
116         fAmpHisto[mod][col][row]->Fill(energy);
117       }
118       else {
119         char hname[128];
120         sprintf(hname,"mod%dcol%drow%d",mod,col,row);
121         fAmpHisto[mod][col][row] = new TH1F(hname,hname,100,0.,1000.);
122       }
123
124       // update histograms in local file every 100th event
125       if(iEvent%fUpdatingRate == 0) {
126         AliInfo(Form("Updating histo file, event %d, run %d\n",iEvent,runNum));
127         UpdateHistoFile();
128       } 
129     }
130
131     iEvent++;
132   }
133
134   UpdateHistoFile();
135   digits->Delete();
136   delete digits;
137   delete geo;
138   AliInfo(Form("%d events of run %d processed.",iEvent,runNum));
139 }
140
141 //-----------------------------------------------------------------------------
142 void AliPHOSCalibHistoProducer::UpdateHistoFile()
143 {
144   // Write histograms to file
145
146   if(!fHistoFile) return;
147   if(!fHistoFile->IsOpen()) return;
148
149   TH1F* hist=0;
150
151   for(Int_t module=0; module<5; module++) {
152     for(Int_t column=0; column<56; column++) {
153       for(Int_t row=0; row<64; row++) {
154         hist = fAmpHisto[module][column][row]; 
155         if(hist) hist->Write(hist->GetName(),TObject::kWriteDelete);
156       }
157     }
158   }
159
160 }