]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONPedestal.cxx
initialisation bug fixed
[u/mrichter/AliRoot.git] / MUON / AliMUONPedestal.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 #include "AliMUONPedestal.h"
19 #include "AliMUONErrorCounter.h"
20 #include "AliMUONVStore.h"
21 #include "AliMUON2DMap.h"
22 #include "AliMUONCalibParamND.h"
23
24 #include <TString.h>
25 #include <THashTable.h>
26 #include <THashList.h>
27 #include <TTimeStamp.h>
28 #include <TMath.h>
29 #include <TTree.h>
30 #include <TFile.h>
31 #include <TH1F.h>
32 #include <Riostream.h>
33
34 #include <sstream>
35
36 //-----------------------------------------------------------------------------
37 /// \class AliMUONPedestal
38 ///
39 /// Implementation of the pedestal computing
40 ///
41 /// add
42 /// 
43 ///
44 /// \author Alberto Baldisseri, JL Charvet (05/05/2009)
45 //-----------------------------------------------------------------------------
46
47 /// \cond CLASSIMP
48 ClassImp(AliMUONPedestal)
49 /// \endcond
50
51 //______________________________________________________________________________
52 AliMUONPedestal::AliMUONPedestal()
53 : TObject(),
54 fN(0),
55 fNEvents(0),
56 fRunNumber(0),
57 fNChannel(0),
58 fNManu(0),
59 fErrorBuspatchTable(new THashTable(100,2)),
60 fManuBuspatchTable(new THashTable(100,2)),
61 fDate(new TTimeStamp()),
62 fFilcout(0),
63 fPedestalStore(new AliMUON2DMap(kFALSE)),
64 fIndex(-1)
65 {
66 /// Default constructor
67   sprintf(fHistoFileName," ");
68   sprintf(fprefixDA," "); 
69 }
70 //  AliMUONPedestal& operator=(const AliMUONPedestal& other); Copy ctor
71
72 //______________________________________________________________________________
73 AliMUONPedestal::~AliMUONPedestal()
74 {
75 /// Destructor
76   delete fPedestalStore;
77 }
78
79 //______________________________________________________________________________
80 void AliMUONPedestal::MakePed(Int_t busPatchId, Int_t manuId, Int_t channelId, Int_t charge)
81 {
82   /// Compute pedestals values
83   AliMUONVCalibParam* ped = 
84     static_cast<AliMUONVCalibParam*>(fPedestalStore ->FindObject(busPatchId, manuId));
85
86   if (!ped) {
87     fNManu++;
88     ped = new AliMUONCalibParamND(2, kNChannels,busPatchId, manuId, -1.); // put default wise -1, not connected channel
89     fPedestalStore ->Add(ped);  
90   }
91
92   // Initialization for the first value
93   if (ped->ValueAsDouble(channelId, 0) == -1) ped->SetValueAsDouble(channelId, 0, 0.);
94   if (ped->ValueAsDouble(channelId, 1) == -1) ped->SetValueAsDouble(channelId, 1, 0.);
95
96   Double_t pedMean  = ped->ValueAsDouble(channelId, 0) + (Double_t) charge;
97   Double_t pedSigma = ped->ValueAsDouble(channelId, 1) + (Double_t) charge*charge;
98
99   ped->SetValueAsDouble(channelId, 0, pedMean);
100   ped->SetValueAsDouble(channelId, 1, pedSigma);
101
102   char bpmanuname[256];
103   AliMUONErrorCounter* manuCounter;
104   sprintf(bpmanuname,"bp%dmanu%d",busPatchId,manuId);           
105   if (!(manuCounter = (AliMUONErrorCounter*)GetManuBuspatchTable()->FindObject(bpmanuname)))
106     {
107       // New (buspatch,manu)
108       manuCounter = new AliMUONErrorCounter(busPatchId,manuId);
109       manuCounter->SetName(bpmanuname);
110       GetManuBuspatchTable()->Add(manuCounter);
111     }
112   else
113     {
114       // Existing buspatch
115       manuCounter->Increment();
116     }   
117 }
118
119 //______________________________________________________________________________
120 void AliMUONPedestal::Finalize()
121 {
122   Double_t pedMean;
123   Double_t pedSigma;
124   Int_t busPatchId;
125   Int_t manuId;
126   Int_t channelId;
127
128   // print in logfile
129   if (fErrorBuspatchTable->GetSize())
130   {
131     cout<<"\n* Buspatches with less statistics (due to parity errors)"<<endl;
132     (*fFilcout)<<"\n* Buspatches with less statistics (due to parity errors)"<<endl;
133     TIterator* iter = fErrorBuspatchTable->MakeIterator();
134     AliMUONErrorCounter* parityerror;
135     while((parityerror = (AliMUONErrorCounter*) iter->Next()))
136     {
137       cout<<"  bp "<<parityerror->BusPatch()<<": events used = "<<fNEvents-parityerror->Events()<<endl;
138       (*fFilcout)<<"  bp "<<parityerror->BusPatch()<<": events used = "<<fNEvents-parityerror->Events()<<endl;
139     }
140
141   }
142
143 // iterator over pedestal
144   TIter next(fPedestalStore ->CreateIterator());
145   AliMUONVCalibParam* ped;
146
147   while ( ( ped = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
148   {
149     busPatchId              = ped->ID0();
150     manuId                  = ped->ID1();
151     if(manuId==0)
152       {
153         cout << " !!! BIG WARNING: ManuId = " << manuId << " !!! in  BP = " << busPatchId << endl;
154         (*fFilcout) << " !!! BIG WARNING: ManuId = " << manuId << " !!! in  BP = " << busPatchId << endl;
155       }
156     Int_t eventCounter;
157     // Correct the number of events for buspatch with errors
158     char bpname[256];
159     AliMUONErrorCounter* errorCounter;
160     sprintf(bpname,"bp%d",busPatchId);
161     if ((errorCounter = (AliMUONErrorCounter*)fErrorBuspatchTable->FindObject(bpname)))
162     {
163       eventCounter = fNEvents - errorCounter->Events();
164     }
165     else
166     {
167       eventCounter = fNEvents;
168     }
169
170     Int_t occupancy;
171     // value of (buspatch, manu) occupancy
172     char bpmanuname[256];
173     AliMUONErrorCounter* manuCounter;
174     sprintf(bpmanuname,"bp%dmanu%d",busPatchId,manuId);
175     manuCounter = (AliMUONErrorCounter*)fManuBuspatchTable->FindObject(bpmanuname);
176     occupancy = manuCounter->Events()/64/eventCounter;
177     if(occupancy>1)
178     {
179         cout << " !!! BIG WARNING: ManuId = " << manuId << " !!! in  BP = " << busPatchId << " occupancy (>1) = " << occupancy << endl;
180         (*fFilcout) << " !!! BIG WARNING: ManuId = " << manuId << " !!! in  BP = " << busPatchId << " occupancy (>1) = " << occupancy <<endl;
181     }
182
183     for (channelId = 0; channelId < ped->Size() ; ++channelId) 
184     {
185       pedMean  = ped->ValueAsDouble(channelId, 0);
186
187       if (pedMean > 0) // connected channels
188       {
189           ped->SetValueAsDouble(channelId, 0, pedMean/(Double_t)eventCounter);
190           pedMean  = ped->ValueAsDouble(channelId, 0);
191           pedSigma = ped->ValueAsDouble(channelId, 1);
192           ped->SetValueAsDouble(channelId, 1, TMath::Sqrt(TMath::Abs(pedSigma/(Double_t)eventCounter - pedMean*pedMean)));
193           if(manuId == 0 || occupancy>1)
194           {
195               ped->SetValueAsDouble(channelId, 0, kADCMax);
196               ped->SetValueAsDouble(channelId, 1, kADCMax);
197           }
198       }
199       else
200       {
201           ped->SetValueAsDouble(channelId, 0, kADCMax);
202           ped->SetValueAsDouble(channelId, 1, kADCMax);
203       }
204     }
205   }
206 }
207
208 //______________________________________________________________________________
209 void AliMUONPedestal::MakeASCIIoutput(ostream& out) const
210 {
211 /// put pedestal store in the output stream
212
213   out<<"//===========================================================================" << endl;
214   out<<"//                       Pedestal file calculated by MUONTRKda"<<endl;
215   out<<"//===========================================================================" << endl;
216   out<<"//       * Run           : " << fRunNumber << endl; 
217   out<<"//       * Date          : " << fDate->AsString("l") <<endl;
218   out<<"//       * Statictics    : " << fNEvents << endl;
219   out<<"//       * # of MANUS    : " << fNManu << endl;
220   out<<"//       * # of channels : " << fNChannel << endl;
221   if (fErrorBuspatchTable->GetSize())
222   {
223     out<<"//"<<endl;
224     out<<"//       * Buspatches with less statistics (due to parity errors)"<<endl;
225     TIterator* iter = fErrorBuspatchTable->MakeIterator();
226     AliMUONErrorCounter* parityerror;
227     while((parityerror = (AliMUONErrorCounter*) iter->Next()))
228     {
229       out<<"//         bp "<<parityerror->BusPatch()<<" events used "<<fNEvents-parityerror->Events()<<endl;
230     }
231   }  
232   out<<"//"<<endl;
233   out<<"//---------------------------------------------------------------------------" << endl;
234   out<<"//---------------------------------------------------------------------------" << endl;
235   out<<"//      BP     MANU     CH.      MEAN    SIGMA"<<endl;
236   out<<"//---------------------------------------------------------------------------" << endl;
237
238   // iterator over pedestal
239   TIter next(fPedestalStore ->CreateIterator());
240   AliMUONVCalibParam* ped;
241   
242   while ( ( ped = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
243   {
244     Int_t busPatchId = ped->ID0();
245     Int_t manuId = ped->ID1();
246     
247     for ( Int_t channelId = 0; channelId < ped->Size(); ++channelId ) 
248     {
249       Double_t pedMean  = ped->ValueAsDouble(channelId, 0);
250       Double_t pedSigma = ped->ValueAsDouble(channelId, 1);
251         
252       out << "\t" << busPatchId << "\t" << manuId <<"\t"<< channelId << "\t" << pedMean <<"\t"<< pedSigma << endl;        
253     }
254   }
255 }
256
257 //______________________________________________________________________________
258 void AliMUONPedestal::MakeControlHistos()
259 {
260
261   if (fIndex>=0) return; // Pedestal run (fIndex=-1)
262
263   Double_t pedMean;
264   Double_t pedSigma;
265   Int_t busPatchId;
266   Int_t manuId;
267   Int_t channelId;
268
269 // histo
270   TFile*  histoFile = 0;
271   TTree* tree = 0;
272   TH1F* pedMeanHisto = 0;
273   TH1F* pedSigmaHisto = 0;
274     
275   sprintf(fHistoFileName,"%s.root",fprefixDA);
276   histoFile = new TFile(fHistoFileName,"RECREATE","MUON Tracking pedestals");
277
278   Char_t name[255];
279   Char_t title[255];
280   sprintf(name,"pedmean_allch");
281   sprintf(title,"Pedestal mean all channels");
282   Int_t nx = kADCMax+1;
283   Int_t xmin = 0;
284   Int_t xmax = kADCMax; 
285   pedMeanHisto = new TH1F(name,title,nx,xmin,xmax);
286   pedMeanHisto->SetDirectory(histoFile);
287
288   sprintf(name,"pedsigma_allch");
289   sprintf(title,"Pedestal sigma all channels");
290   nx = 201;
291   xmin = 0;
292   xmax = 200; 
293   pedSigmaHisto = new TH1F(name,title,nx,xmin,xmax);
294   pedSigmaHisto->SetDirectory(histoFile);
295
296   tree = new TTree("t","Pedestal tree");
297   tree->Branch("bp",&busPatchId,"bp/I");
298   tree->Branch("manu",&manuId,",manu/I");
299   tree->Branch("channel",&channelId,",channel/I");
300   tree->Branch("pedMean",&pedMean,",pedMean/D");
301   tree->Branch("pedSigma",&pedSigma,",pedSigma/D");
302
303   // iterator over pedestal
304   TIter next(fPedestalStore ->CreateIterator());
305   AliMUONVCalibParam* ped;
306   
307   while ( ( ped = dynamic_cast<AliMUONVCalibParam*>(next() ) ) )
308   {
309     busPatchId = ped->ID0();
310     manuId = ped->ID1();
311     
312     for ( channelId = 0; channelId < ped->Size(); ++channelId ) 
313     {
314       pedMean  = ped->ValueAsDouble(channelId, 0);
315       pedSigma = ped->ValueAsDouble(channelId, 1);
316           
317       pedMeanHisto->Fill(pedMean);
318       pedSigmaHisto->Fill(pedSigma);
319       tree->Fill();  
320     }
321   }
322     
323   histoFile->Write();  
324   histoFile->Close(); 
325
326 }