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