Line data Source code
1 : // AliAODDimuon: a class for AODs for the MUON Arm of the ALICE Experiment
2 : // Author: P. Cortese, Universita' del Piemonte Orientale in Alessandria and
3 : // INFN of Torino - Italy
4 : //
5 : // The class defines a dimuon pair object from two AliAODTrack objects.
6 : // AliAODDimuon objects are supposed to be added to the AliAODEvent structure
7 : // during analysis. They would then allow to calculate the dimuon-related
8 : // kinematic variables with a minimal disk occupancy.
9 : // The payload of the class has been reduced to two pointers to the two
10 : // tracks. An instance of this class has also to be added to the AliAODEvent
11 : // structure to provide additional information that is specific to MUON and
12 : // therefore has not been included into the AOD header.
13 : // Two transient data members are not stored on file as they can be recomputed
14 : // at runtime.
15 : //
16 :
17 : #include "AliAODDimuon.h"
18 : #include "TLorentzVector.h"
19 : #define AliAODDimuon_CXX
20 :
21 170 : ClassImp(AliAODDimuon)
22 :
23 : //______________________________________________________________________________
24 14 : AliAODDimuon::AliAODDimuon():AliVParticle(),fP(0),fMProton(0.93827231)
25 10 : {
26 : // default constructor
27 2 : fMu[0]=0;
28 2 : fMu[1]=0;
29 4 : }
30 :
31 : //______________________________________________________________________________
32 0 : AliAODDimuon::AliAODDimuon(const AliAODDimuon& dimu):
33 0 : AliVParticle(dimu),
34 0 : fP(0),fMProton(0.93827231)
35 0 : {
36 : // copy constructor
37 0 : fMu[0]=dimu.Mu(0);
38 0 : fMu[1]=dimu.Mu(1);
39 0 : }
40 :
41 : //______________________________________________________________________________
42 : AliAODDimuon &AliAODDimuon::operator=(const AliAODDimuon& dimu)
43 : {
44 : // assignment operator
45 0 : if(&dimu != this){
46 0 : delete fP;
47 0 : fP=0;
48 0 : fMProton=0.93827231;
49 0 : fMu[0]=dimu.Mu(0);
50 0 : fMu[1]=dimu.Mu(1);
51 0 : }
52 0 : return *this;
53 : }
54 :
55 : //______________________________________________________________________________
56 0 : AliAODDimuon::AliAODDimuon(TObject *mu0, TObject *mu1):
57 0 : fP(0),fMProton(0.93827231)
58 0 : {
59 : // Creates a dimuon pair from two tracks
60 :
61 0 : fMu[0]=mu0;
62 0 : fMu[1]=mu1;
63 0 : }
64 :
65 : //______________________________________________________________________________
66 : AliAODDimuon::~AliAODDimuon()
67 12 : {
68 : // destructor
69 2 : delete fP;
70 12 : }
71 :
72 : //______________________________________________________________________________
73 : void AliAODDimuon::Clear(Option_t*)
74 : {
75 : /// delete our internal memory
76 0 : delete fP;
77 0 : fP = 0x0;
78 0 : }
79 :
80 : //______________________________________________________________________________
81 : Double_t AliAODDimuon::Px() const {
82 : // Px of the dimuon
83 0 : if(CheckPointers())return -999999999;
84 0 : return ((AliAODTrack*)fMu[0].GetObject())->Px()+
85 0 : ((AliAODTrack*)fMu[1].GetObject())->Px();
86 0 : }
87 :
88 : //______________________________________________________________________________
89 : Double_t AliAODDimuon::Py() const {
90 : // Py of the dimuon
91 0 : if(CheckPointers())return -999999999;
92 0 : return ((AliAODTrack*)fMu[0].GetObject())->Py()+
93 0 : ((AliAODTrack*)fMu[1].GetObject())->Py();
94 0 : }
95 :
96 : //______________________________________________________________________________
97 : Double_t AliAODDimuon::Pz() const {
98 : // Pz of the dimuon
99 0 : if(CheckPointers())return -999999999;
100 0 : return ((AliAODTrack*)fMu[0].GetObject())->Pz()+
101 0 : ((AliAODTrack*)fMu[1].GetObject())->Pz();
102 0 : }
103 :
104 : //______________________________________________________________________________
105 : Double_t AliAODDimuon::Pt() const {
106 : // Pt of the dimuon
107 0 : if(CheckPointers())return -999999999;
108 0 : Double_t px=Px();
109 0 : Double_t py=Py();
110 0 : return TMath::Sqrt(px*px+py*py);
111 0 : }
112 :
113 : //______________________________________________________________________________
114 : Double_t AliAODDimuon::E() const
115 : {
116 : // Dimuon energy
117 :
118 0 : if(CheckPointers())return -999999999;
119 :
120 0 : return ((AliAODTrack*)fMu[0].GetObject())->E()+ ((AliAODTrack*)fMu[1].GetObject())->E();
121 0 : }
122 :
123 : //______________________________________________________________________________
124 : Double_t AliAODDimuon::P() const {
125 : // Dimuon momentum
126 0 : if(CheckPointers())return -999999999;
127 0 : return TLV()->P();
128 0 : }
129 :
130 : //______________________________________________________________________________
131 : Double_t AliAODDimuon::M() const {
132 : // Dimuon invariant mass
133 0 : if(CheckPointers())return -999999999;
134 0 : return TLV()->M();
135 0 : }
136 :
137 : //______________________________________________________________________________
138 : Double_t AliAODDimuon::Eta() const {
139 : // Dimuon pseudorapidity
140 0 : if(CheckPointers())return -999999999;
141 0 : return TLV()->Eta();
142 0 : }
143 :
144 : //______________________________________________________________________________
145 : Double_t AliAODDimuon::Phi() const {
146 : // Dimuon asimuthal angle
147 0 : if(CheckPointers())return -999999999;
148 0 : return TLV()->Phi();
149 0 : }
150 :
151 :
152 : //______________________________________________________________________________
153 : Double_t AliAODDimuon::Theta() const {
154 : // Dimuon polar angle
155 0 : if(CheckPointers())return -999999999;
156 0 : return TLV()->Theta();
157 0 : }
158 :
159 : //______________________________________________________________________________
160 : Double_t AliAODDimuon::Y() const {
161 : // Dimuon rapidity
162 0 : if(CheckPointers())return -999999999;
163 0 : return TLV()->Rapidity();
164 0 : }
165 :
166 : //______________________________________________________________________________
167 : Short_t AliAODDimuon::Charge() const {
168 : // Dimuon charge
169 0 : if(CheckPointers())return -999;
170 0 : return ((AliAODTrack*)fMu[0].GetObject())->Charge()+
171 0 : ((AliAODTrack*)fMu[1].GetObject())->Charge();
172 0 : }
173 :
174 : //______________________________________________________________________________
175 : Int_t AliAODDimuon::CheckPointers() const{
176 : // Checks if the track pointers have been initialized
177 0 : if(fMu[0]==0||fMu[1]==0){
178 0 : printf("Dimuon not initialized\n");
179 0 : return -999;
180 : }
181 0 : if((fMu[0].GetObject())==0||(fMu[1].GetObject())==0){
182 0 : printf("Can not get objects\n");
183 0 : return -999;
184 : }
185 0 : return 0;
186 0 : }
187 :
188 : //______________________________________________________________________________
189 : void AliAODDimuon::SetMu(Int_t imu, AliAODTrack *mu){
190 : // Assign a track pointer
191 0 : if (imu==0||imu==1){
192 0 : delete fP;
193 0 : fP=0;
194 0 : fMu[imu]=mu;
195 0 : }
196 0 : }
197 :
198 : //______________________________________________________________________________
199 : void AliAODDimuon::SetMuons(AliAODTrack *mu0, AliAODTrack *mu1){
200 : // Assign the track pointers
201 0 : fMu[0]=mu0;
202 0 : fMu[1]=mu1;
203 0 : delete fP;
204 0 : fP=0;
205 0 : }
206 :
207 : //______________________________________________________________________________
208 : Double_t AliAODDimuon::XF() {
209 : // Dimuon Feynman x
210 :
211 0 : if(CheckPointers())return -999999999;
212 :
213 : //Double_t ebeam=((AliAODEventInfo*)fEi.GetObject())->EBeam();
214 : Double_t ebeam = 3500.; // temporary
215 0 : if(ebeam<=0){
216 0 : printf("AliAODDimuon::xf: can not compute xf with EBeam=%f\n",ebeam);
217 0 : return -999999999;
218 : }
219 0 : Double_t mDimu=M();
220 0 : Double_t pMax=TMath::Sqrt(ebeam*ebeam-mDimu*mDimu);
221 0 : return Pz()/pMax;
222 0 : }
223 :
224 : //______________________________________________________________________________
225 : Double_t AliAODDimuon::CostCS(){
226 : // Calculation the Collins-Soper angle (adapted from code by R. Arnaldi)
227 0 : if(CheckPointers())return -999999999;
228 : Double_t ebeam=3500.; //temporary
229 0 : if(ebeam<=0){
230 0 : printf("Can not compute costCS with EBeam=%f\n",ebeam);
231 0 : return -999999999;
232 : }
233 0 : Double_t mp=fMProton;
234 0 : Double_t pbeam=TMath::Sqrt(ebeam*ebeam-mp*mp);
235 0 : Double_t pla10=((AliAODTrack*)fMu[0].GetObject())->Px();
236 0 : Double_t pla11=((AliAODTrack*)fMu[0].GetObject())->Py();
237 0 : Double_t pla12=((AliAODTrack*)fMu[0].GetObject())->Pz();
238 0 : Double_t e1=((AliAODTrack*)fMu[0].GetObject())->E();
239 0 : Double_t mu1Charge=((AliAODTrack*)fMu[0].GetObject())->Charge();
240 0 : Double_t pla20=((AliAODTrack*)fMu[1].GetObject())->Px();
241 0 : Double_t pla21=((AliAODTrack*)fMu[1].GetObject())->Py();
242 0 : Double_t pla22=((AliAODTrack*)fMu[1].GetObject())->Pz();
243 0 : Double_t e2=((AliAODTrack*)fMu[1].GetObject())->E();
244 0 : Double_t mu2Charge=((AliAODTrack*)fMu[1].GetObject())->Charge();
245 :
246 : // Fill the Lorentz vector for projectile and target
247 : // For the moment we do not consider the crossing angle
248 : // Projectile runs towards the MUON arm
249 0 : TLorentzVector pProjCM(0.,0.,-pbeam,ebeam); // projectile
250 0 : TLorentzVector pTargCM(0.,0., pbeam,ebeam); // target
251 : //
252 : // --- Get the muons parameters in the CM frame
253 : //
254 0 : TLorentzVector pMu1CM(pla10,pla11,pla12,e1);
255 0 : TLorentzVector pMu2CM(pla20,pla21,pla22,e2);
256 : //
257 : // --- Obtain the dimuon parameters in the CM frame
258 : //
259 0 : TLorentzVector pDimuCM=pMu1CM+pMu2CM;
260 : //
261 : // --- Translate the dimuon parameters in the dimuon rest frame
262 : //
263 0 : TVector3 beta=(-1./pDimuCM.E())*pDimuCM.Vect();
264 0 : TLorentzVector pMu1Dimu=pMu1CM;
265 0 : TLorentzVector pMu2Dimu=pMu2CM;
266 0 : TLorentzVector pProjDimu=pProjCM;
267 0 : TLorentzVector pTargDimu=pTargCM;
268 0 : pMu1Dimu.Boost(beta);
269 0 : pMu2Dimu.Boost(beta);
270 0 : pProjDimu.Boost(beta);
271 0 : pTargDimu.Boost(beta);
272 : //
273 : // --- Determine the z axis for the CS angle
274 : //
275 0 : TVector3 zaxisCS=(((pProjDimu.Vect()).Unit())-((pTargDimu.Vect()).Unit())).Unit();
276 : //
277 : // --- Determine the CS angle (angle between mu+ and the z axis defined above)
278 : //
279 : Double_t cost;
280 0 : if(mu1Charge > 0) {
281 0 : cost = zaxisCS.Dot((pMu1Dimu.Vect()).Unit());
282 : // Theta CS is not properly defined for Like-Sign muons
283 0 : if(mu2Charge > 0 && cost<0) cost=-cost;
284 : } else {
285 : // Theta CS is not properly defined for Like-Sign muons
286 0 : cost = zaxisCS.Dot((pMu2Dimu.Vect()).Unit());
287 0 : if(mu2Charge < 0 && cost<0) cost=-cost;
288 : }
289 : return cost;
290 0 : }
291 :
292 : //______________________________________________________________________________
293 : Double_t AliAODDimuon::CostHe(){
294 : // Calculation the Helicity polarization angle (adapted from code by R. Arnaldi)
295 0 : if(CheckPointers())return -999999999;
296 : Double_t ebeam=3500; //temporary
297 0 : if(ebeam<=0){
298 0 : printf("Can not compute costCS with EBeam=%f\n",ebeam);
299 0 : return -999999999;
300 : }
301 0 : Double_t pla10=((AliAODTrack*)fMu[0].GetObject())->Px();
302 0 : Double_t pla11=((AliAODTrack*)fMu[0].GetObject())->Py();
303 0 : Double_t pla12=((AliAODTrack*)fMu[0].GetObject())->Pz();
304 0 : Double_t e1=((AliAODTrack*)fMu[0].GetObject())->E();
305 0 : Double_t mu1Charge=((AliAODTrack*)fMu[0].GetObject())->Charge();
306 0 : Double_t pla20=((AliAODTrack*)fMu[1].GetObject())->Px();
307 0 : Double_t pla21=((AliAODTrack*)fMu[1].GetObject())->Py();
308 0 : Double_t pla22=((AliAODTrack*)fMu[1].GetObject())->Pz();
309 0 : Double_t e2=((AliAODTrack*)fMu[1].GetObject())->E();
310 0 : Double_t mu2Charge=((AliAODTrack*)fMu[1].GetObject())->Charge();
311 : //
312 : // --- Get the muons parameters in the CM frame
313 : //
314 0 : TLorentzVector pMu1CM(pla10,pla11,pla12,e1);
315 0 : TLorentzVector pMu2CM(pla20,pla21,pla22,e2);
316 : //
317 : // --- Obtain the dimuon parameters in the CM frame
318 : //
319 0 : TLorentzVector pDimuCM=pMu1CM+pMu2CM;
320 : //
321 : // --- Translate the dimuon parameters in the dimuon rest frame
322 : //
323 0 : TVector3 beta=(-1./pDimuCM.E())*pDimuCM.Vect();
324 0 : TLorentzVector pMu1Dimu=pMu1CM;
325 0 : TLorentzVector pMu2Dimu=pMu2CM;
326 0 : pMu1Dimu.Boost(beta);
327 0 : pMu2Dimu.Boost(beta);
328 : //
329 : // --- Determine the z axis for the calculation of the polarization angle
330 : // (i.e. the direction of the dimuon in the CM system)
331 : //
332 0 : TVector3 zaxis;
333 0 : zaxis=(pDimuCM.Vect()).Unit();
334 : //
335 : // --- Calculation of the polarization angle (Helicity)
336 : // (angle between mu+ and the z axis defined above)
337 : //
338 : Double_t cost;
339 0 : if(mu1Charge > 0) {
340 0 : cost = zaxis.Dot((pMu1Dimu.Vect()).Unit());
341 : // Theta Helicity is not properly defined for Like-Sign muons
342 0 : if(mu2Charge > 0 && cost<0) cost=-cost;
343 : } else {
344 0 : cost = zaxis.Dot((pMu2Dimu.Vect()).Unit());
345 : // Theta Helicity is not properly defined for Like-Sign muons
346 0 : if(mu2Charge < 0 && cost<0) cost=-cost;
347 : }
348 : return cost;
349 0 : }
350 :
351 : //________________________________________________________________________
352 : Double_t AliAODDimuon::PhiCS(){
353 : // Cosinus of the Collins-Soper polar decay angle
354 0 : if(CheckPointers())return -999999999;
355 : Double_t ebeam=3500.; //temporary
356 0 : if(ebeam<=0){
357 0 : printf("Can not compute phiCS with EBeam=%f\n",ebeam);
358 0 : return -999999999;
359 : }
360 0 : Double_t mp=fMProton;
361 0 : Double_t pbeam=TMath::Sqrt(ebeam*ebeam-mp*mp);
362 0 : Double_t pla10=((AliAODTrack*)fMu[0].GetObject())->Px();
363 0 : Double_t pla11=((AliAODTrack*)fMu[0].GetObject())->Py();
364 0 : Double_t pla12=((AliAODTrack*)fMu[0].GetObject())->Pz();
365 0 : Double_t e1=((AliAODTrack*)fMu[0].GetObject())->E();
366 0 : Double_t mu1Charge=((AliAODTrack*)fMu[0].GetObject())->Charge();
367 0 : Double_t pla20=((AliAODTrack*)fMu[1].GetObject())->Px();
368 0 : Double_t pla21=((AliAODTrack*)fMu[1].GetObject())->Py();
369 0 : Double_t pla22=((AliAODTrack*)fMu[1].GetObject())->Pz();
370 0 : Double_t e2=((AliAODTrack*)fMu[1].GetObject())->E();
371 :
372 : // Fill the Lorentz vector for projectile and target
373 : // For the moment we do not consider the crossing angle
374 : // Projectile runs towards the MUON arm
375 0 : TLorentzVector pProjCM(0.,0.,-pbeam,ebeam); // projectile
376 0 : TLorentzVector pTargCM(0.,0., pbeam,ebeam); // target
377 : //
378 : // --- Get the muons parameters in the CM frame
379 : //
380 0 : TLorentzVector pMu1CM(pla10,pla11,pla12,e1);
381 0 : TLorentzVector pMu2CM(pla20,pla21,pla22,e2);
382 : //
383 : // --- Obtain the dimuon parameters in the CM frame
384 : //
385 0 : TLorentzVector pDimuCM=pMu1CM+pMu2CM;
386 : //
387 : // --- Translate the dimuon parameters in the dimuon rest frame
388 : //
389 0 : TVector3 beta=(-1./pDimuCM.E())*pDimuCM.Vect();
390 0 : TLorentzVector pMu1Dimu=pMu1CM;
391 0 : TLorentzVector pMu2Dimu=pMu2CM;
392 0 : TLorentzVector pProjDimu=pProjCM;
393 0 : TLorentzVector pTargDimu=pTargCM;
394 0 : pMu1Dimu.Boost(beta);
395 0 : pMu2Dimu.Boost(beta);
396 0 : pProjDimu.Boost(beta);
397 0 : pTargDimu.Boost(beta);
398 : //
399 : // --- Determine the z axis for the CS angle
400 : //
401 0 : TVector3 zaxisCS=(((pProjDimu.Vect()).Unit())-((pTargDimu.Vect()).Unit())).Unit();
402 : //
403 : // --- Determine the CS angle (angle between mu+ and the z axis defined above)
404 : //
405 0 : TVector3 yaxisCS=(((pProjDimu.Vect()).Unit()).Cross((pTargDimu.Vect()).Unit())).Unit();
406 0 : TVector3 xaxisCS=(yaxisCS.Cross(zaxisCS)).Unit();
407 :
408 : Double_t phi;
409 0 : if(mu1Charge>0) phi = TMath::ATan2((pMu1Dimu.Vect()).Dot(yaxisCS),((pMu1Dimu.Vect()).Dot(xaxisCS)));
410 0 : else phi = TMath::ATan2((pMu2Dimu.Vect()).Dot(yaxisCS),((pMu2Dimu.Vect()).Dot(xaxisCS)));
411 :
412 : return phi;
413 0 : }
414 :
415 : //______________________________________________________________________________
416 : Double_t AliAODDimuon::PhiHe(){
417 : // Calculation the Helicity aimuthal angle (adapted from code by R. Arnaldi)
418 0 : if(CheckPointers())return -999999999;
419 : Double_t ebeam=3500; //temporary
420 0 : if(ebeam<=0){
421 0 : printf("Can not compute phiHE with EBeam=%f\n",ebeam);
422 0 : return -999999999;
423 : }
424 0 : Double_t pbeam=TMath::Sqrt(ebeam*ebeam-fMProton*fMProton);
425 0 : Double_t pla10=((AliAODTrack*)fMu[0].GetObject())->Px();
426 0 : Double_t pla11=((AliAODTrack*)fMu[0].GetObject())->Py();
427 0 : Double_t pla12=((AliAODTrack*)fMu[0].GetObject())->Pz();
428 0 : Double_t e1=((AliAODTrack*)fMu[0].GetObject())->E();
429 0 : Double_t mu1Charge=((AliAODTrack*)fMu[0].GetObject())->Charge();
430 0 : Double_t pla20=((AliAODTrack*)fMu[1].GetObject())->Px();
431 0 : Double_t pla21=((AliAODTrack*)fMu[1].GetObject())->Py();
432 0 : Double_t pla22=((AliAODTrack*)fMu[1].GetObject())->Pz();
433 0 : Double_t e2=((AliAODTrack*)fMu[1].GetObject())->E();
434 :
435 : // Fill the Lorentz vector for projectile and target
436 : // For the moment we consider no crossing angle
437 : // Projectile runs towards the MUON arm
438 0 : TLorentzVector pProjCM(0.,0.,-pbeam,ebeam); // projectile
439 0 : TLorentzVector pTargCM(0.,0., pbeam,ebeam); // target
440 : //
441 : // --- Get the muons parameters in the CM frame
442 : //
443 0 : TLorentzVector pMu1CM(pla10,pla11,pla12,e1);
444 0 : TLorentzVector pMu2CM(pla20,pla21,pla22,e2);
445 : //
446 : // --- Obtain the dimuon parameters in the CM frame
447 : //
448 0 : TLorentzVector pDimuCM=pMu1CM+pMu2CM;
449 : //
450 : // --- Translate the muon parameters in the dimuon rest frame
451 : //
452 0 : TVector3 zaxis=(pDimuCM.Vect()).Unit();
453 : //
454 : // --- Translate the dimuon parameters in the dimuon rest frame
455 : //
456 0 : TVector3 beta=(-1./pDimuCM.E())*pDimuCM.Vect();
457 0 : TLorentzVector pMu1Dimu=pMu1CM;
458 0 : TLorentzVector pMu2Dimu=pMu2CM;
459 0 : pMu1Dimu.Boost(beta);
460 0 : pMu2Dimu.Boost(beta);
461 :
462 0 : TLorentzVector pProjDimu=pProjCM;
463 0 : TLorentzVector pTargDimu=pTargCM;
464 0 : pProjDimu.Boost(beta);
465 0 : pTargDimu.Boost(beta);
466 :
467 0 : TVector3 yaxis=((pProjDimu.Vect()).Cross(pTargDimu.Vect())).Unit();
468 0 : TVector3 xaxis=(yaxis.Cross(zaxis)).Unit();
469 : //
470 : // --- Calculation of the azimuthal angle (Helicity)
471 : //
472 : Double_t phi;
473 0 : if(mu1Charge>0) phi = TMath::ATan2((pMu1Dimu.Vect()).Dot(yaxis),(pMu1Dimu.Vect()).Dot(xaxis));
474 0 : else phi = TMath::ATan2((pMu2Dimu.Vect()).Dot(yaxis),(pMu2Dimu.Vect()).Dot(xaxis));
475 :
476 : return phi;
477 0 : }
478 :
479 : //______________________________________________________________________________
480 : Int_t AliAODDimuon::AnyPt(){
481 : // Test if the two muons match two trigger tracks
482 0 : if(CheckPointers())return 0;
483 0 : return (((AliAODTrack*)fMu[0].GetObject())->MatchTrigger())&&
484 0 : (((AliAODTrack*)fMu[1].GetObject())->MatchTrigger());
485 0 : }
486 :
487 : //______________________________________________________________________________
488 : Int_t AliAODDimuon::LowPt(){
489 : // Test if the two muons match two trigger tracks with a "Low Pt" cut
490 0 : if(CheckPointers())return 0;
491 0 : return (((AliAODTrack*)fMu[0].GetObject())->MatchTriggerLowPt())&&
492 0 : (((AliAODTrack*)fMu[1].GetObject())->MatchTriggerLowPt());
493 0 : }
494 :
495 : //______________________________________________________________________________
496 : Int_t AliAODDimuon::HighPt(){
497 : // Test if the two muons match two trigger tracks with a "High Pt" cut
498 0 : if(CheckPointers())return 0;
499 0 : return (((AliAODTrack*)fMu[0].GetObject())->MatchTriggerHighPt())&&
500 0 : (((AliAODTrack*)fMu[1].GetObject())->MatchTriggerHighPt());
501 0 : }
502 :
503 : //______________________________________________________________________________
504 : Double_t AliAODDimuon::MaxChi2Match(){
505 : // Maximum matching Chi2 between track and trigger track
506 0 : if(CheckPointers())return -999999999;
507 0 : return TMath::Max((((AliAODTrack*)fMu[0].GetObject())->GetChi2MatchTrigger()),
508 0 : (((AliAODTrack*)fMu[1].GetObject())->GetChi2MatchTrigger()));
509 0 : }
510 :
511 : //______________________________________________________________________________
512 : TLorentzVector* AliAODDimuon::TLV() const
513 : {
514 0 : delete fP;
515 0 : fP = new TLorentzVector(Px(),Py(),Pz(),E());
516 0 : return fP;
517 0 : }
|