1 /*****************************************************************
3 General utility class for doing pp vs multiplicity studies
5 --- Functionality being added on demand.
7 --- Please report any bugs, comments, suggestions to:
8 david.dobrigkeit.chinellato@cern.ch
10 *****************************************************************/
12 #include "AliVEvent.h"
13 #include "AliESDEvent.h"
14 #include "AliAODEvent.h"
15 #include "AliVVertex.h"
17 #include "AliAODVertex.h"
18 #include "AliVTrack.h"
19 #include "AliVEvent.h"
20 #include <TMatrixDSym.h>
22 #include "AliESDUtils.h"
23 #include "AliPPVsMultUtils.h"
27 ClassImp(AliPPVsMultUtils)
29 //______________________________________________________________________
30 AliPPVsMultUtils::AliPPVsMultUtils():TObject(),
32 fCalibrationLoaded(0),
33 fBoundaryHisto_V0M(0),
34 fBoundaryHisto_V0A(0),
35 fBoundaryHisto_V0C(0),
36 fBoundaryHisto_V0MEq(0),
37 fBoundaryHisto_V0AEq(0),
38 fBoundaryHisto_V0CEq(0)
44 //_____________________________________________________________________________
45 AliPPVsMultUtils::AliPPVsMultUtils(const AliPPVsMultUtils &c) : TObject(c),
47 fCalibrationLoaded(0),
48 fBoundaryHisto_V0M(0),
49 fBoundaryHisto_V0A(0),
50 fBoundaryHisto_V0C(0),
51 fBoundaryHisto_V0MEq(0),
52 fBoundaryHisto_V0AEq(0),
53 fBoundaryHisto_V0CEq(0)
56 // copy constructor - untested
58 ((AliPPVsMultUtils &) c).Copy(*this);
61 //_____________________________________________________________________________
62 AliPPVsMultUtils &AliPPVsMultUtils::operator=(const AliPPVsMultUtils &c)
65 // Assignment operator - untested
68 if (this != &c) ((AliPPVsMultUtils &) c).Copy(*this);
72 //______________________________________________________________________
73 Float_t AliPPVsMultUtils::GetMultiplicityPercentile(AliVEvent *event, TString lMethod)
74 // Function to get multiplicity quantiles
76 Int_t lRequestedRunNumber = event->GetRunNumber();
77 if( lRequestedRunNumber != fRunNumber ) fCalibrationLoaded = LoadCalibration( lRequestedRunNumber );
79 if ( !fCalibrationLoaded ){
80 return -1000; //Return absurd value (hopefully noone will use this...)
82 if ( !fBoundaryHisto_V0M || !fBoundaryHisto_V0A || !fBoundaryHisto_V0C ||
83 !fBoundaryHisto_V0MEq || !fBoundaryHisto_V0AEq || !fBoundaryHisto_V0CEq ){
84 return -1000; //Return absurd value (hopefully noone will use this...)
87 Float_t lreturnval = -1;
89 //Get VZERO Information for multiplicity later
90 AliVVZERO* esdV0 = event->GetVZEROData();
92 AliError("AliVVZERO not available");
97 Float_t multV0A = 0; // multiplicity from V0 reco side A
98 Float_t multV0C = 0; // multiplicity from V0 reco side C
99 Float_t multV0AEq = 0; // multiplicity from V0 reco side A
100 Float_t multV0CEq = 0; // multiplicity from V0 reco side C
101 Float_t multV0ACorr = 0; // multiplicity from V0 reco side A
102 Float_t multV0CCorr = 0; // multiplicity from V0 reco side C
104 //Non-Equalized Signal: copy of multV0ACorr and multV0CCorr from AliCentralitySelectionTask
105 //Getters for uncorrected multiplicity
106 multV0A=esdV0->GetMTotV0A();
107 multV0C=esdV0->GetMTotV0C();
109 //Getting around to the SPD vertex -> typecast to ESD/AOD
110 const AliVVertex *lPrimarySPDVtx = NULL;
111 /* get ESD vertex SPD */
112 if (event->InheritsFrom("AliESDEvent")) {
113 AliESDEvent *esdevent = dynamic_cast<AliESDEvent *>(event);
114 if (!esdevent) return kFALSE;
115 lPrimarySPDVtx = esdevent->GetPrimaryVertexSPD();
117 /* get AOD vertex SPD */
118 else if (event->InheritsFrom("AliAODEvent")) {
119 AliAODEvent *aodevent = dynamic_cast<AliAODEvent *>(event);
120 if (!aodevent) return kFALSE;
121 lPrimarySPDVtx = aodevent->GetPrimaryVertexSPD();
124 //Get Z vertex position of SPD vertex (why not Tracking if available?)
125 Float_t zvtx = lPrimarySPDVtx->GetZ();
127 //Acquire Corrected multV0A
128 multV0ACorr = AliESDUtils::GetCorrV0A(multV0A,zvtx);
129 multV0CCorr = AliESDUtils::GetCorrV0C(multV0C,zvtx);
131 // Equalized signals // From AliCentralitySelectionTask // Updated
132 for(Int_t iCh = 32; iCh < 64; ++iCh) {
133 Double_t mult = event->GetVZEROEqMultiplicity(iCh);
136 for(Int_t iCh = 0; iCh < 32; ++iCh) {
137 Double_t mult = event->GetVZEROEqMultiplicity(iCh);
141 if ( lMethod == "V0M" ) lreturnval =
142 fBoundaryHisto_V0M -> GetBinContent( fBoundaryHisto_V0M->FindBin(multV0ACorr+multV0CCorr) );
143 if ( lMethod == "V0A" ) lreturnval =
144 fBoundaryHisto_V0A -> GetBinContent( fBoundaryHisto_V0A->FindBin(multV0ACorr) );
145 if ( lMethod == "V0C" ) lreturnval =
146 fBoundaryHisto_V0C -> GetBinContent( fBoundaryHisto_V0C->FindBin(multV0CCorr) );
148 if ( lMethod == "V0MEq" ) lreturnval =
149 fBoundaryHisto_V0MEq -> GetBinContent( fBoundaryHisto_V0MEq->FindBin(multV0AEq+multV0CEq) );
150 if ( lMethod == "V0AEq" ) lreturnval =
151 fBoundaryHisto_V0AEq -> GetBinContent( fBoundaryHisto_V0AEq->FindBin(multV0AEq) );
152 if ( lMethod == "V0CEq" ) lreturnval =
153 fBoundaryHisto_V0CEq -> GetBinContent( fBoundaryHisto_V0CEq->FindBin(multV0CEq) );
158 //______________________________________________________________________
159 Bool_t AliPPVsMultUtils::LoadCalibration(Int_t lLoadThisCalibration)
160 //To be called if starting analysis on a new run
162 //If Histograms exist, de-allocate to prevent memory leakage
163 if( fBoundaryHisto_V0M ) {fBoundaryHisto_V0M->Delete(); fBoundaryHisto_V0M = 0x0; }
164 if( fBoundaryHisto_V0A ) {fBoundaryHisto_V0A->Delete(); fBoundaryHisto_V0A = 0x0; }
165 if( fBoundaryHisto_V0C ) {fBoundaryHisto_V0C->Delete(); fBoundaryHisto_V0C = 0x0; }
166 if( fBoundaryHisto_V0MEq ) {fBoundaryHisto_V0MEq->Delete(); fBoundaryHisto_V0MEq = 0x0; }
167 if( fBoundaryHisto_V0AEq ) {fBoundaryHisto_V0AEq->Delete(); fBoundaryHisto_V0AEq = 0x0; }
168 if( fBoundaryHisto_V0CEq ) {fBoundaryHisto_V0CEq->Delete(); fBoundaryHisto_V0CEq = 0x0; }
170 AliInfo(Form( "Loading calibration file for run %i",lLoadThisCalibration) );
171 TFile *lCalibFile_V0M = 0x0;
172 TFile *lCalibFile_V0A = 0x0;
173 TFile *lCalibFile_V0C = 0x0;
174 TFile *lCalibFile_V0MEq = 0x0;
175 TFile *lCalibFile_V0AEq = 0x0;
176 TFile *lCalibFile_V0CEq = 0x0;
178 //AliInfo("Calling TFile::Open");
179 lCalibFile_V0M = TFile::Open("$ALICE_ROOT/PWGLF/STRANGENESS/Cascades/corrections/calibration_V0M.root");
180 lCalibFile_V0A = TFile::Open("$ALICE_ROOT/PWGLF/STRANGENESS/Cascades/corrections/calibration_V0A.root");
181 lCalibFile_V0C = TFile::Open("$ALICE_ROOT/PWGLF/STRANGENESS/Cascades/corrections/calibration_V0C.root");
182 lCalibFile_V0MEq = TFile::Open("$ALICE_ROOT/PWGLF/STRANGENESS/Cascades/corrections/calibration_V0MEq.root");
183 lCalibFile_V0AEq = TFile::Open("$ALICE_ROOT/PWGLF/STRANGENESS/Cascades/corrections/calibration_V0AEq.root");
184 lCalibFile_V0CEq = TFile::Open("$ALICE_ROOT/PWGLF/STRANGENESS/Cascades/corrections/calibration_V0CEq.root");
186 //AliInfo("Casting");
187 fBoundaryHisto_V0M = dynamic_cast<TH1F *>(lCalibFile_V0M -> Get(Form("histocalib%i",lLoadThisCalibration)) );
188 fBoundaryHisto_V0A = dynamic_cast<TH1F *>(lCalibFile_V0A -> Get(Form("histocalib%i",lLoadThisCalibration)) );
189 fBoundaryHisto_V0C = dynamic_cast<TH1F *>(lCalibFile_V0C -> Get(Form("histocalib%i",lLoadThisCalibration)) );
190 fBoundaryHisto_V0MEq = dynamic_cast<TH1F *>(lCalibFile_V0MEq -> Get(Form("histocalib%i",lLoadThisCalibration)) );
191 fBoundaryHisto_V0AEq = dynamic_cast<TH1F *>(lCalibFile_V0AEq -> Get(Form("histocalib%i",lLoadThisCalibration)) );
192 fBoundaryHisto_V0CEq = dynamic_cast<TH1F *>(lCalibFile_V0CEq -> Get(Form("histocalib%i",lLoadThisCalibration)) );
194 if ( !fBoundaryHisto_V0M || !fBoundaryHisto_V0A || !fBoundaryHisto_V0C ||
195 !fBoundaryHisto_V0MEq || !fBoundaryHisto_V0AEq || !fBoundaryHisto_V0CEq ){
196 AliInfo(Form("No calibration for run %i exists at the moment!",lLoadThisCalibration));
197 if( lCalibFile_V0M ) {lCalibFile_V0M->Close(); lCalibFile_V0M->Delete(); delete lCalibFile_V0M; }
198 if( lCalibFile_V0A ) {lCalibFile_V0A->Close(); lCalibFile_V0A->Delete(); delete lCalibFile_V0A; }
199 if( lCalibFile_V0C ) {lCalibFile_V0C->Close(); lCalibFile_V0C->Delete(); delete lCalibFile_V0C; }
200 if( lCalibFile_V0MEq ) {lCalibFile_V0MEq->Close(); lCalibFile_V0MEq->Delete(); delete lCalibFile_V0MEq; }
201 if( lCalibFile_V0AEq ) {lCalibFile_V0AEq->Close(); lCalibFile_V0AEq->Delete(); delete lCalibFile_V0AEq; }
202 if( lCalibFile_V0CEq ) {lCalibFile_V0CEq->Close(); lCalibFile_V0CEq->Delete(); delete lCalibFile_V0CEq; }
203 fRunNumber = lLoadThisCalibration;
204 return kFALSE; //return denial
207 //Careful with manual cleanup if needed: to be implemented
208 fBoundaryHisto_V0M->SetDirectory(0);
209 fBoundaryHisto_V0A->SetDirectory(0);
210 fBoundaryHisto_V0C->SetDirectory(0);
211 fBoundaryHisto_V0MEq->SetDirectory(0);
212 fBoundaryHisto_V0AEq->SetDirectory(0);
213 fBoundaryHisto_V0CEq->SetDirectory(0);
215 //AliInfo("Closing");
216 if( lCalibFile_V0M ) {lCalibFile_V0M->Close(); lCalibFile_V0M->Delete(); delete lCalibFile_V0M; }
217 if( lCalibFile_V0A ) {lCalibFile_V0A->Close(); lCalibFile_V0A->Delete(); delete lCalibFile_V0A; }
218 if( lCalibFile_V0C ) {lCalibFile_V0C->Close(); lCalibFile_V0C->Delete(); delete lCalibFile_V0C; }
219 if( lCalibFile_V0MEq ) {lCalibFile_V0MEq->Close(); lCalibFile_V0MEq->Delete(); delete lCalibFile_V0MEq; }
220 if( lCalibFile_V0AEq ) {lCalibFile_V0AEq->Close(); lCalibFile_V0AEq->Delete(); delete lCalibFile_V0AEq; }
221 if( lCalibFile_V0CEq ) {lCalibFile_V0CEq->Close(); lCalibFile_V0CEq->Delete(); delete lCalibFile_V0CEq; }
223 fRunNumber = lLoadThisCalibration; //Loaded!
224 AliInfo(Form("Finished loading calibration for run %i",lLoadThisCalibration));