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