]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/TenderSupplies/AliVZEROTenderSupply.cxx
Update VZERO reco in order to deal with high lumi data. Weighted mean on V0 A and...
[u/mrichter/AliRoot.git] / ANALYSIS / TenderSupplies / AliVZEROTenderSupply.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
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 //  Recalculate VZERO timing and decision using the tender                   //
20 //  (in case the LHC phase drift is updated in OCDB)                         //
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23
24 #include <TList.h>
25 #include <TObjString.h>
26 #include <TChain.h>
27 #include <TF1.h>
28
29 #include <AliLog.h>
30 #include <AliESDEvent.h>
31 #include <AliESDVZERO.h>
32 #include <AliCDBId.h>
33 #include <AliCDBManager.h>
34 #include <AliCDBEntry.h>
35 #include <AliTender.h>
36 #include <AliLHCClockPhase.h>
37 #include <AliVZEROCalibData.h>
38 #include <AliVZEROTriggerMask.h>
39 #include <AliVZEROReconstructor.h>
40
41 #include "AliVZEROTenderSupply.h"
42
43 ClassImp(AliVZEROTenderSupply)
44
45 AliVZEROTenderSupply::AliVZEROTenderSupply() :
46   AliTenderSupply(),
47   fCalibData(NULL),
48   fTimeSlewing(NULL),
49   fRecoParam(NULL),
50   fLHCClockPhase(0),
51   fDebug(kFALSE)
52 {
53   //
54   // default ctor
55   //
56 }
57
58 //_____________________________________________________
59 AliVZEROTenderSupply::AliVZEROTenderSupply(const char *name, const AliTender *tender) :
60   AliTenderSupply(name,tender),
61   fCalibData(NULL),
62   fTimeSlewing(NULL),
63   fRecoParam(NULL),
64   fLHCClockPhase(0),
65   fDebug(kFALSE)
66 {
67   //
68   // named ctor
69   //
70 }
71
72 //_____________________________________________________
73 void AliVZEROTenderSupply::Init()
74 {
75   //
76   // Initialise VZERO tender
77   //
78 }
79
80 //_____________________________________________________
81 void AliVZEROTenderSupply::ProcessEvent()
82 {
83   //
84   // Reapply the LHC-clock phase drift
85   //
86
87   AliESDEvent *event=fTender->GetEvent();
88   if (!event) return;
89   
90   //load gain correction if run has changed
91   if (fTender->RunChanged()){
92     if (fDebug) printf("AliVZEROTenderSupply::ProcessEvent - Run Changed (%d)\n",fTender->GetRun());
93     GetPhaseCorrection();
94
95     AliCDBEntry *entryGeom = fTender->GetCDBManager()->Get("GRP/Geometry/Data",fTender->GetRun());
96     if (!entryGeom) {
97       AliError("No geometry entry is found");
98       return;
99     } else {
100       if (fDebug) printf("AliVZEROTenderSupply::Used geometry entry: %s\n",entryGeom->GetId().ToString().Data());
101     }
102
103     AliCDBEntry *entryCal = fTender->GetCDBManager()->Get("VZERO/Calib/Data",fTender->GetRun());
104     if (!entryCal) {
105       AliError("No VZERO calibration entry is found");
106       fCalibData = NULL;
107       return;
108     } else {
109       fCalibData = (AliVZEROCalibData*)entryCal->GetObject();
110       if (fDebug) printf("AliVZEROTenderSupply::Used VZERO calibration entry: %s\n",entryCal->GetId().ToString().Data());
111     }
112
113     AliCDBEntry *entrySlew = fTender->GetCDBManager()->Get("VZERO/Calib/TimeSlewing",fTender->GetRun());
114     if (!entrySlew) {
115       AliError("VZERO time slewing function is not found in OCDB !");
116       fTimeSlewing = NULL;
117       return;
118     } else {
119       fTimeSlewing = (TF1*)entrySlew->GetObject();
120       if (fDebug) printf("AliVZEROTenderSupply::Used VZERO time slewing entry: %s\n",entrySlew->GetId().ToString().Data());
121     }
122
123     AliCDBEntry *entryRecoParam = fTender->GetCDBManager()->Get("VZERO/Calib/RecoParam",fTender->GetRun());
124     if (!entryRecoParam) {
125       AliError("VZERO reco-param object is not found in OCDB !");
126       fRecoParam = NULL;
127       return;
128     } else {
129       fRecoParam = (AliVZERORecoParam*)entryRecoParam->GetObject();
130       if (fDebug) printf("AliVZEROTenderSupply::Used VZERO reco-param entry: %s\n",entryRecoParam->GetId().ToString().Data());
131     }
132   }
133
134   if (!fCalibData || !fTimeSlewing || !fRecoParam) {
135     AliWarning("VZERO calibration objects not found!");
136     return;
137   }
138   //
139   // correct VZERO time signals and decision
140   //
141   AliESDVZERO *esdVZERO = event->GetVZEROData();
142   if (!esdVZERO) {
143     AliError("No VZERO object is found inside ESD!");
144     return;
145   }
146   if (!esdVZERO->TestBit(AliESDVZERO::kDecisionFilled)) {
147     AliWarning("VZERO offline trigger decisions were not filled in ESD, the tender supply is disabled");
148     return;
149   }
150
151   if (fDebug) printf("LHC-clock phase correction: %f\n",fLHCClockPhase);
152
153   if (fDebug) printf("Original VZERO decision %d (%f ns) and %d (%f ns)\n",
154                      esdVZERO->GetV0ADecision(),esdVZERO->GetV0ATime(),
155                      esdVZERO->GetV0CDecision(),esdVZERO->GetV0CTime());
156   Float_t time[64];
157   for(Int_t i = 0; i < 64; ++i) {
158     time[i] = esdVZERO->GetTime(i);
159     if (time[i] > (AliVZEROReconstructor::kInvalidTime + 1e-6))
160       time[i] += fLHCClockPhase;
161   }
162   esdVZERO->SetTime(time);
163
164   {
165     AliVZEROTriggerMask triggerMask;
166     triggerMask.SetRecoParam(fRecoParam);
167     triggerMask.FillMasks(esdVZERO, fCalibData, fTimeSlewing);
168   }
169   if (fDebug) printf("Modified VZERO decision %d (%f ns) and %d (%f ns)\n",
170                      esdVZERO->GetV0ADecision(),esdVZERO->GetV0ATime(),
171                      esdVZERO->GetV0CDecision(),esdVZERO->GetV0CTime());
172
173 }
174
175 //_____________________________________________________
176 void AliVZEROTenderSupply::GetPhaseCorrection()
177 {
178   //
179   // Get Gain splines from OCDB
180   //
181
182   AliInfo("Get LHC-clock phase correction");
183
184   //
185   //find previous entry from the UserInfo
186   //
187   TTree *tree=((TChain*)fTender->GetInputData(0))->GetTree();
188   if (!tree) {
189     AliError("Tree not found in ESDhandler");
190     return;
191   }
192   
193   TList *userInfo=(TList*)tree->GetUserInfo();
194   if (!userInfo) {
195     AliError("No UserInfo found in tree");
196     return;
197   }
198
199   TList *cdbList=(TList*)userInfo->FindObject("cdbList");
200   if (!cdbList) {
201     AliError("No cdbList found in UserInfo");
202     if (AliLog::GetGlobalLogLevel()>=AliLog::kError) userInfo->Print();
203     return;
204   }
205
206   Float_t oldPhase = 0;
207
208   TIter nextCDB(cdbList);
209   TObjString *os=0x0;
210   while ( (os=(TObjString*)nextCDB()) ){
211     if (!(os->GetString().Contains("GRP/Calib/LHCClockPhase"))) continue;
212     AliCDBId *id=AliCDBId::MakeFromString(os->GetString());
213     
214     AliCDBEntry *entry=fTender->GetCDBManager()->Get(*id);
215     if (!entry) {
216       AliError("The previous LHC-clock phase entry is not found");
217       delete id;
218       return;
219     }
220     
221     if (fDebug) printf("AliVZEROTenderSupply::Used old LHC-clock phase entry: %s\n",entry->GetId().ToString().Data());
222     
223     AliLHCClockPhase *phase = (AliLHCClockPhase*)entry->GetObject();
224     if (!phase) {
225       AliError("Phase object is not found in the calibration entry");
226       delete id;
227       return;
228     }
229
230     oldPhase = phase->GetMeanPhase();
231
232     delete id;
233     break;
234   }
235
236   //
237   //new LHC-clock phase entry
238   //
239   Float_t newPhase = 0;
240   AliCDBEntry *entryNew=fTender->GetCDBManager()->Get("GRP/Calib/LHCClockPhase",fTender->GetRun());
241   if (!entryNew) {
242     AliError("No new LHC-clock phase calibration entry is found");
243     return;
244   }
245   if (fDebug) printf("AliVZEROTenderSupply::Used new LHC-clock phase entry: %s\n",entryNew->GetId().ToString().Data());
246   
247   AliLHCClockPhase *phase2 = (AliLHCClockPhase*)entryNew->GetObject();
248   if (!phase2) {
249     AliError("Phase object is not found in the calibration entry");
250     return;
251   }
252
253   newPhase = phase2->GetMeanPhase();
254
255   fLHCClockPhase = newPhase - oldPhase;
256 }