]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSCpvPreprocessor.cxx
Restored compilation on Windows/Cygwin
[u/mrichter/AliRoot.git] / PHOS / AliPHOSCpvPreprocessor.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2008, 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
17 ///////////////////////////////////////////////////////////////////////////////
18 // CPV Preprocessor class. It runs by Shuttle at the end of the run,
19 // calculates calibration coefficients and dead/bad channels
20 // to be posted in OCDB.
21 //
22 // Author: Boris Polichtchouk, 25 January 2008
23 ///////////////////////////////////////////////////////////////////////////////
24
25 #include "AliPHOSCpvPreprocessor.h"
26 #include "AliLog.h"
27 #include "AliCDBMetaData.h"
28 #include "AliPHOSCpvCalibData.h"
29 #include "TFile.h"
30 #include "TH1.h"
31 #include "TMap.h"
32 #include "TRandom.h"
33 #include "TKey.h"
34 #include "TList.h"
35 #include "TObjString.h"
36
37 ClassImp(AliPHOSCpvPreprocessor)
38
39 //_______________________________________________________________________________________
40 AliPHOSCpvPreprocessor::AliPHOSCpvPreprocessor() :
41 AliPreprocessor("CPV",0)
42 {
43   //default constructor
44 }
45
46 //_______________________________________________________________________________________
47 AliPHOSCpvPreprocessor::AliPHOSCpvPreprocessor(AliShuttleInterface* shuttle):
48 AliPreprocessor("CPV",shuttle)
49 {
50   // Constructor
51 }
52
53 //_______________________________________________________________________________________
54 UInt_t AliPHOSCpvPreprocessor::Process(TMap* /*valueSet*/)
55 {
56   // process data retrieved by the Shuttle
57
58   // The fileName with the histograms which have been produced by
59   // AliPHOSCpvDAs.
60   // It is a responsibility of the SHUTTLE framework to form the fileName.
61   
62   TString runType = GetRunType();
63   Log(Form("Run type: %s",runType.Data()));
64
65   gRandom->SetSeed(0); //the seed is set to the current  machine clock!
66   AliPHOSCpvCalibData calibData;
67
68   TList* list = GetFileSources(kDAQ, "AMPLITUDES");
69   if(!list) {
70     Log("Sources list not found, exit.");
71     return 1;
72   }
73
74   TIter iter(list);
75   TObjString *source;
76   
77   while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
78     AliInfo(Form("found source %s", source->String().Data()));
79
80     TString fileName = GetFile(kDAQ, "AMPLITUDES", source->GetName());
81     AliInfo(Form("Got filename: %s",fileName.Data()));
82
83     TFile f(fileName);
84
85     if(!f.IsOpen()) {
86       Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
87       return 1;
88     }
89
90     const Int_t nMod=5;   // 1:5 modules
91     const Int_t nCol=56;  // 1:56 columns in each CPV module (along the global Z axis)
92     const Int_t nRow=128; // 1:128 rows in each CPV module
93
94     Double_t coeff;
95     char hnam[80];
96     TH1F* histo=0;
97     
98     //Get the reference histogram
99     //(author: Gustavo Conesa Balbastre)
100
101     TList * keylist = f.GetListOfKeys();
102     Int_t nkeys   = f.GetNkeys();
103     Bool_t ok = kFALSE;
104     TKey  *key;
105     TString refHistoName= "";
106     Int_t ikey = 0;
107     Int_t counter = 0;
108     TH1F* hRef = 0;
109     
110     //Check if the file contains any histogram
111     
112     if(nkeys< 2){
113       Log(Form("Not enough histograms (%d) for calibration.",nkeys));
114       return 1;
115     }
116
117     while(!ok){
118       ikey = gRandom->Integer(nkeys);
119       key = (TKey*)keylist->At(ikey);
120       refHistoName = key->GetName();
121       hRef = (TH1F*)f.Get(refHistoName);
122       counter++;
123       // Check if the reference histogram has too little statistics
124       if(hRef->GetEntries()>2) ok=kTRUE;
125       if(!ok && counter >= nkeys){
126         Log("No histogram with enough statistics for reference.");
127         return 1;
128       }
129     }
130     
131     Log(Form("reference histogram %s, %.1f entries, mean=%.3f, rms=%.3f.",
132              hRef->GetName(),hRef->GetEntries(),
133              hRef->GetMean(),hRef->GetRMS()));
134
135     Double_t refMean=hRef->GetMean();
136     
137     // Calculates relative calibration coefficients for all non-zero channels
138     
139     for(Int_t mod=0; mod<nMod; mod++) {
140       for(Int_t col=0; col<nCol; col++) {
141         for(Int_t row=0; row<nRow; row++) {
142           sprintf(hnam,"%d_%d_%d",mod,row,col); // mod_X_Z
143           histo = (TH1F*)f.Get(hnam);
144           //TODO: dead channels exclusion!
145           if(histo) {
146             coeff = histo->GetMean()/refMean;
147             if(coeff>0)
148               calibData.SetADCchannelCpv(mod+1,col+1,row+1,1./coeff);
149             AliInfo(Form("mod %d col %d row %d  coeff %f\n",mod,col,row,coeff));
150           }
151         }
152       }
153     }
154     
155     f.Close();
156   }
157   
158   //Store CPV calibration data
159   
160   AliCDBMetaData cpvMetaData;
161   Bool_t cpvOK = Store("Calib", "CpvGainPedestals", &calibData, &cpvMetaData);
162
163   if(cpvOK) return 0;
164   else
165     return 1;
166
167 }
168