]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/RESONANCES/AliRsnMother.cxx
Added an experimental class for single-track monitors.
[u/mrichter/AliRoot.git] / PWG2 / RESONANCES / AliRsnMother.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 //  This class implements a candidate resonance. It has two pointers to its
19 //  two candidate daughters, whose 4-momenta are combined to obtain the mother
20 //  invariant mass and other kinematical quantities.
21 //  This class contains also some methods used to compute kinematical relations
22 //  between the candidate resonance and other particles.
23 //
24 //  authors: A. Pulvirenti (alberto.pulvirenti@ct.infn.it)
25 //           M. Vala (martin.vala@cern.ch)
26 //
27 ////////////////////////////////////////////////////////////////////////////////
28
29 #include <Riostream.h>
30 #include <TVector3.h>
31
32 #include "AliAODMCParticle.h"
33 #include "AliMCParticle.h"
34 #include "AliRsnEvent.h"
35
36 #include "AliRsnMother.h"
37
38 ClassImp(AliRsnMother)
39
40 //__________________________________________________________________________________________________
41 AliRsnMother::AliRsnMother(const AliRsnMother &obj) :
42    TObject(obj),
43    fRefEvent(obj.fRefEvent),
44    fSum(obj.fSum),
45    fSumMC(obj.fSumMC),
46    fRef(obj.fRef),
47    fRefMC(obj.fRefMC)
48 {
49 //
50 // Copy constructor.
51 // Does not duplicate pointers.
52 //
53
54    Int_t i;
55    for (i = 0; i < 2; i++) fDaughter[i] = obj.fDaughter[i];
56 }
57
58 //__________________________________________________________________________________________________
59 AliRsnMother& AliRsnMother::operator=(const AliRsnMother &obj)
60 {
61 //
62 // Assignment operator.
63 // Does not duplicate pointers.
64 //
65
66    fSum = obj.fSum;
67    fRef = obj.fRef;
68    fSumMC = obj.fSumMC;
69    fRefMC = obj.fRefMC;
70    fRefEvent = obj.fRefEvent;
71    fDaughter[0] = obj.fDaughter[0];
72    fDaughter[1] = obj.fDaughter[1];
73    
74    return (*this);
75 }
76
77 //__________________________________________________________________________________________________
78 AliRsnMother::~AliRsnMother()
79 {
80 //
81 // Desctructor.
82 // Does nothing, since pointers are not created in this class.
83 //
84 }
85
86 //_______________________________________________________________________________________________________________________
87 void AliRsnMother::Reset()
88 {
89 //
90 // Resets the mother, zeroing all data members.
91 //
92
93    Int_t i;
94    for (i = 0; i < 2; i++) fDaughter[i] = 0x0;
95    fRefEvent = 0x0;
96    fSum.SetXYZM(0.0, 0.0, 0.0, 0.0);
97    fRef.SetXYZM(0.0, 0.0, 0.0, 0.0);
98 }
99
100 //__________________________________________________________________________________________________
101 Int_t AliRsnMother::CommonMother() const
102 {
103 //
104 // If MC info is available, checks if the two tracks in the pair have the same mother.
105 // If the mother label is the same, the function returns the PDG code of mother,
106 // otherwise it returns 0.
107 // The two arguments passed by reference contain the GEANT labels of the mother
108 // of the two particles to which the two daughters point. This is for being able 
109 // to check if they are really coming from a resonance (indexes >= 0) or not.
110 //
111
112    Int_t m1  = fDaughter[0]->GetMother();
113    Int_t m2  = fDaughter[1]->GetMother();
114    Int_t out = 0;
115    
116    // a true mother makes sense only if both mothers
117    // are not-negative and equal
118    if (m1 >= 0 && m2 >= 0 && m1 == m2) {
119       out = TMath::Abs(fDaughter[0]->GetMotherPDG());
120       AliDebugClass(1, Form("Mothers are: %d %d --> EQUAL (PDG = %d)", m1, m2, out));
121    } 
122    
123    return out;
124 }
125
126 //__________________________________________________________________________________________________
127 Double_t AliRsnMother::AngleToLeading(Bool_t &success)
128 {
129 //
130 // Compute the angle betwee this and the leading particls
131 // of the reference event (if this was set properly).
132 // In case one of the two daughters is the leading, return
133 // a meaningless value, in order to skip this pair.
134 // if second argument is kTRUE, use MC values.
135 //
136
137    if (!fRefEvent) {
138       success = kFALSE;
139       return -99.0;
140    }
141    
142    Int_t id1 = fDaughter[0]->GetID();
143    Int_t id2 = fDaughter[1]->GetID();
144    Int_t idL = fRefEvent->GetLeadingIndex();
145    
146    if (id1 == idL || id2 == idL) {
147       success = kFALSE;
148       return -99.0;
149    }
150    
151    AliRsnDaughter leading = fRefEvent->GetDaughter(idL, kFALSE);
152    AliVParticle  *ref     = leading.GetRef();
153    Double_t       angle   = ref->Phi() - fSum.Phi();
154    
155    //return angle w.r.t. leading particle in the range -pi/2, 3/2pi
156    while (angle >= 1.5 * TMath::Pi()) angle -= 2 * TMath::Pi();
157    while (angle < -0.5 * TMath::Pi()) angle += 2 * TMath::Pi();
158    success = kTRUE;
159    
160    return angle;
161 }
162
163 //__________________________________________________________________________________________________
164 void AliRsnMother::ComputeSum(Double_t mass0, Double_t mass1, Double_t motherMass)
165 {
166 //
167 // Sets the masses for the 4-momenta of the daughters and then
168 // sums them, taking into account that the space part is set to
169 // each of them when the reference object is set (see AliRsnDaughter::SetRef)
170 //
171
172    fDaughter[0]->FillP(mass0);
173    fDaughter[1]->FillP(mass1);
174
175    // sum
176    fSum   = fDaughter[0]->Prec() + fDaughter[1]->Prec();
177    fSumMC = fDaughter[0]->Psim() + fDaughter[1]->Psim();
178    
179    // reference
180    fRef.SetXYZM(fSum.X(), fSum.Y(), fSum.Z(), motherMass);
181    fRefMC.SetXYZM(fSumMC.X(), fSumMC.Y(), fSumMC.Z(), motherMass);
182 }
183
184 //__________________________________________________________________________________________________
185 Double_t AliRsnMother::CosThetaStar(Bool_t first, Bool_t useMC)
186 {
187 //
188 // Computes the cosine of theta*, which is the angle of one of the daughters
189 // with respect to the total momentum of the resonance, in its rest frame.
190 // The arguments are needed to choose which of the daughters one want to use
191 // and if reconstructed or MC momentum must be used.
192 // [Contribution from Z. Feckova]
193 //
194
195    TLorentzVector &mother    = Sum(useMC);
196    TLorentzVector &daughter0 = (first ? fDaughter[0]->P(useMC) : fDaughter[1]->P(useMC));
197    TLorentzVector &daughter1 = (first ? fDaughter[1]->P(useMC) : fDaughter[0]->P(useMC));
198    TVector3 momentumM(mother.Vect());
199    TVector3 normal(mother.Y() / momentumM.Mag(), -mother.X() / momentumM.Mag(), 0.0);
200
201    // Computes first the invariant mass of the mother
202    Double_t mass0      = fDaughter[0]->P(useMC).M();
203    Double_t mass1      = fDaughter[1]->P(useMC).M();
204    Double_t p0         = daughter0.Vect().Mag();
205    Double_t p1         = daughter1.Vect().Mag();
206    Double_t E0         = TMath::Sqrt(mass0 * mass0 + p0 * p0);
207    Double_t E1         = TMath::Sqrt(mass1 * mass1 + p1 * p1);
208    Double_t MotherMass = TMath::Sqrt((E0 + E1) * (E0 + E1) - (p0 * p0 + 2.0 * daughter0.Vect().Dot(daughter1.Vect()) + p1 * p1));
209    MotherMass = fSum.M();
210
211    // Computes components of beta
212    Double_t betaX = -mother.X() / mother.E();
213    Double_t betaY = -mother.Y() / mother.E();
214    Double_t betaZ = -mother.Z() / mother.E();
215
216    // Computes Lorentz transformation of the momentum of the first daughter
217    // into the rest frame of the mother and theta*
218    daughter0.Boost(betaX, betaY, betaZ);
219    TVector3 momentumD = daughter0.Vect();
220
221    Double_t cosThetaStar = normal.Dot(momentumD) / momentumD.Mag();
222
223    return cosThetaStar;
224 }
225
226 //__________________________________________________________________________________________________
227 void AliRsnMother::PrintInfo(const Option_t * /*option*/) const
228 {
229 //
230 // Print some info of the pair.
231 // The options are passed to the AliRsnDaughter::Print() method
232 //
233
234    AliInfo("======== BEGIN PAIR INFO ===========");
235    AliInfo("Track #1");
236    fDaughter[0]->Print();
237    AliInfo("Track #2");
238    fDaughter[1]->Print();
239    AliInfo("========= END PAIR INFO ===========");
240 }
241
242 //__________________________________________________________________________________________________
243 Bool_t AliRsnMother::CheckPair(Bool_t checkMC) const
244 {
245 //
246 // Checks that the pair is well initialized:
247 // - both daughters are good pointers
248 // - if MC is required, both daughters have a MC reference
249 //
250
251    if (!fDaughter[0] || !fDaughter[1]) {
252       AliError("One of the two tracks is NULL in this pair!");
253       return kFALSE;
254    }
255
256    if (checkMC) {
257       if (fDaughter[0]->GetRefMC() == 0x0 || fDaughter[1]->GetRefMC() == 0x0) {
258          AliError("Required MC info but not all MC refs are available");
259          return kFALSE;
260       }
261    }
262
263    return kTRUE;
264 }