]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSCpvPreprocessor.cxx
PHOS trigger information in the ESD.
[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   AddRunType("PHYSICS");
53   AddRunType("STANDALONE");
54
55 }
56
57 //_______________________________________________________________________________________
58 UInt_t AliPHOSCpvPreprocessor::Process(TMap* /*valueSet*/)
59 {
60   // process data retrieved by the Shuttle
61
62   // The fileName with the histograms which have been produced by
63   // AliPHOSCpvDAs.
64   // It is a responsibility of the SHUTTLE framework to form the fileName.
65   
66   TString runType = GetRunType();
67   Log(Form("Run type: %s",runType.Data()));
68
69   if(runType=="PHYSICS") {
70
71     gRandom->SetSeed(0); //the seed is set to the current  machine clock!
72     AliPHOSCpvCalibData calibData;
73
74     TList* list = GetFileSources(kDAQ, "AMPLITUDES");
75     if(!list) {
76       Log("Sources list not found, exit.");
77       return 1;
78     }
79
80     TIter iter(list);
81     TObjString *source;
82   
83     while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
84       AliInfo(Form("found source %s", source->String().Data()));
85
86       TString fileName = GetFile(kDAQ, "AMPLITUDES", source->GetName());
87       AliInfo(Form("Got filename: %s",fileName.Data()));
88
89       TFile f(fileName);
90
91       if(!f.IsOpen()) {
92         Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
93         return 1;
94       }
95
96       const Int_t nMod=5;   // 1:5 modules
97       const Int_t nCol=56;  // 1:56 columns in each CPV module (along the global Z axis)
98       const Int_t nRow=128; // 1:128 rows in each CPV module
99
100       Double_t coeff;
101       char hnam[80];
102       TH1F* histo=0;
103     
104       //Get the reference histogram
105       //(author: Gustavo Conesa Balbastre)
106
107       TList * keylist = f.GetListOfKeys();
108       Int_t nkeys   = f.GetNkeys();
109       Bool_t ok = kFALSE;
110       TKey  *key;
111       TString refHistoName= "";
112       Int_t ikey = 0;
113       Int_t counter = 0;
114       TH1F* hRef = 0;
115     
116       //Check if the file contains any histogram
117     
118       if(nkeys< 2){
119         Log(Form("Not enough histograms (%d) for calibration.",nkeys));
120         return 1;
121       }
122
123       while(!ok){
124         ikey = gRandom->Integer(nkeys);
125         key = (TKey*)keylist->At(ikey);
126         refHistoName = key->GetName();
127         hRef = (TH1F*)f.Get(refHistoName);
128         counter++;
129         // Check if the reference histogram has too little statistics
130         if(hRef->GetEntries()>2) ok=kTRUE;
131         if(!ok && counter >= nkeys){
132           Log("No histogram with enough statistics for reference.");
133           return 1;
134         }
135       }
136     
137       Log(Form("reference histogram %s, %.1f entries, mean=%.3f, rms=%.3f.",
138                hRef->GetName(),hRef->GetEntries(),
139                hRef->GetMean(),hRef->GetRMS()));
140
141       Double_t refMean=hRef->GetMean();
142     
143       // Calculates relative calibration coefficients for all non-zero channels
144       
145       for(Int_t mod=0; mod<nMod; mod++) {
146         for(Int_t col=0; col<nCol; col++) {
147           for(Int_t row=0; row<nRow; row++) {
148             snprintf(hnam,80,"%d_%d_%d",mod,row,col); // mod_X_Z
149             histo = (TH1F*)f.Get(hnam);
150             //TODO: dead channels exclusion!
151             if(histo) {
152               coeff = histo->GetMean()/refMean;
153               if(coeff>0)
154                 calibData.SetADCchannelCpv(mod+1,col+1,row+1,1./coeff);
155               AliInfo(Form("mod %d col %d row %d  coeff %f\n",mod,col,row,coeff));
156             }
157           }
158         }
159       }
160     
161       f.Close();
162     }
163   
164     //Store CPV calibration data
165   
166     AliCDBMetaData cpvMetaData;
167     Bool_t cpvOK = Store("Calib", "CpvGainPedestals", &calibData, &cpvMetaData);
168
169     if(cpvOK) return 0;
170     else
171       return 1;
172   }
173
174   Log(Form("Unknown or unused run type %s. Do nothing and return OK.",runType.Data()));
175   return 0;
176
177 }
178