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