1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
20 #include "AliOmegaDalitz.h"
26 #include <TParticle.h>
27 #include <TDatabasePDG.h>
28 #include <TLorentzVector.h>
29 #include <TClonesArray.h>
31 ClassImp(AliOmegaDalitz)
34 //-----------------------------------------------------------------------------
36 // Generate lepton-pair mass distributions for Dalitz decays according
37 // to the Kroll-Wada parametrization: N. Kroll, W. Wada: Phys. Rev 98(1955)1355
39 // For the electromagnetic form factor the parameterization from
40 // Lepton-G is used: L.G. Landsberg et al.: Phys. Rep. 128(1985)301
42 //-----------------------------------------------------------------------------
44 // Adaption for AliRoot
49 AliOmegaDalitz::AliOmegaDalitz():
58 void AliOmegaDalitz::Init()
61 Int_t idecay, ibin, nbins = 1000;
62 Double_t min, max, binwidth;
63 Double_t pmass, lmass, omass, emass, mmass;
64 Double_t epsilon, delta, mLL, q, kwHelp, krollWada, formFactor, weight;
66 // Get the particle masses
68 emass = (TDatabasePDG::Instance()->GetParticle(kElectron))->Mass();
70 mmass = (TDatabasePDG::Instance()->GetParticle(kMuonPlus))->Mass();
72 pmass = (TDatabasePDG::Instance()->GetParticle(223)) ->Mass();
74 omass = (TDatabasePDG::Instance()->GetParticle(kPi0)) ->Mass();
76 for ( idecay = 1; idecay<=2; idecay++ )
85 binwidth = (max - min) / (Double_t)nbins;
87 fEPMass = new TH1F("fEPMass", "Dalitz electron pair", nbins, min, max);
89 fMPMass = new TH1F("fMPMass", "Dalitz muon pair", nbins, min, max);
91 epsilon = (lmass / pmass) * (lmass / pmass);
92 delta = (omass / pmass) * (omass / pmass);
94 for ( ibin = 1; ibin <= nbins; ibin++ )
96 mLL = min + (Double_t)(ibin - 1) * binwidth + binwidth / 2.0;
97 q = (mLL / pmass) * (mLL / pmass);
98 if ( q <= 4.0 * epsilon )
100 AliFatal("Error in calculating Dalitz mass histogram binning!");
103 kwHelp = (1.0 + q / (1.0 - delta)) * (1.0 + q / (1.0 - delta))
104 - 4.0 * q / ((1.0 - delta) * (1.0 - delta));
107 AliFatal("Error in calculating Dalitz mass histogram binning!");
109 krollWada = (2.0 / mLL) * TMath::Exp(1.5 * TMath::Log(kwHelp))
110 * TMath::Sqrt(1.0 - 4.0 * epsilon / q)
111 * (1.0 + 2.0 * epsilon / q);
113 (TMath::Power(TMath::Power(0.6519,2),2))
114 / (TMath::Power(TMath::Power(0.6519,2)-TMath::Power(mLL, 2), 2)
115 + TMath::Power(0.04198, 2)*TMath::Power(0.6519, 2));
116 weight = krollWada * formFactor;
118 fEPMass->AddBinContent(ibin, weight);
120 fMPMass->AddBinContent(ibin, weight);
126 void AliOmegaDalitz::Decay(Int_t idlepton, TLorentzVector* pparent)
128 //-----------------------------------------------------------------------------
130 // Generate Omega Dalitz decay
132 //-----------------------------------------------------------------------------
139 Double_t pmass, lmass, omass, lpmass;
140 Double_t e1, p1, e3, p3;
141 Double_t betaSquare, lambda;
142 Double_t costheta, sintheta, cosphi, sinphi, phi;
144 // Get the particle masses
146 lmass = (TDatabasePDG::Instance()->GetParticle(idlepton))->Mass();
148 pmass = pparent->M();
150 omass = (TDatabasePDG::Instance()->GetParticle(kPi0)) ->Mass();
152 // Sample the lepton pair mass from a histogram
155 if ( idlepton == kElectron )
156 lpmass = fEPMass->GetRandom();
158 lpmass = fMPMass->GetRandom();
159 if ( pmass - omass > lpmass && lpmass / 2. > lmass ) break;
162 // lepton pair kinematics in virtual photon rest frame
164 p1 = TMath::Sqrt((e1 + lmass) * (e1 - lmass));
165 betaSquare = 1.0 - 4.0 * (lmass * lmass) / (lpmass * lpmass);
166 lambda = betaSquare / (2.0 - betaSquare);
167 costheta = (2.0 * gRandom->Rndm()) - 1.;
168 sintheta = TMath::Sqrt((1. + costheta) * (1. - costheta));
169 phi = 2.0 * TMath::ACos(-1.) * gRandom->Rndm();
170 sinphi = TMath::Sin(phi);
171 cosphi = TMath::Cos(phi);
172 // momentum vectors of leptons in virtual photon rest frame
173 Double_t pProd1[3] = {p1 * sintheta * cosphi,
174 p1 * sintheta * sinphi,
176 Double_t pProd2[3] = {-1.0 * p1 * sintheta * cosphi,
177 -1.0 * p1 * sintheta * sinphi,
178 -1.0 * p1 * costheta};
180 // pizero kinematics in omega rest frame
181 e3 = (pmass * pmass + omass * omass - lpmass * lpmass)/(2. * pmass);
182 p3 = TMath::Sqrt((e3 + omass) * (e3 - omass));
183 costheta = (2.0 * gRandom->Rndm()) - 1.;
184 sintheta = TMath::Sqrt((1. + costheta) * (1. - costheta));
185 phi = 2.0 * TMath::ACos(-1.) * gRandom->Rndm();
186 sinphi = TMath::Sin(phi);
187 cosphi = TMath::Cos(phi);
188 // pizero 4-vector in omega rest frame
189 fProducts[2].SetPx(p3 * sintheta * cosphi);
190 fProducts[2].SetPy(p3 * sintheta * sinphi);
191 fProducts[2].SetPz(p3 * costheta);
192 fProducts[2].SetE(e3);
194 // lepton 4-vectors in properly rotated virtual photon rest frame
195 Double_t pRot1[3] = {0.};
196 Rot(pProd1, pRot1, costheta, -sintheta, -cosphi, -sinphi);
197 Double_t pRot2[3] = {0.};
198 Rot(pProd2, pRot2, costheta, -sintheta, -cosphi, -sinphi);
199 fProducts[0].SetPx(pRot1[0]);
200 fProducts[0].SetPy(pRot1[1]);
201 fProducts[0].SetPz(pRot1[2]);
202 fProducts[0].SetE(e1);
203 fProducts[1].SetPx(pRot2[0]);
204 fProducts[1].SetPy(pRot2[1]);
205 fProducts[1].SetPz(pRot2[2]);
206 fProducts[1].SetE(e1);
208 // boost the dilepton into the omega's rest frame
209 Double_t eLPparent = TMath::Sqrt(p3 * p3 + lpmass * lpmass);
210 TVector3 boostPair( -1.0 * fProducts[2].Px() / eLPparent,
211 -1.0 * fProducts[2].Py() / eLPparent,
212 -1.0 * fProducts[2].Pz() / eLPparent);
213 fProducts[0].Boost(boostPair);
214 fProducts[1].Boost(boostPair);
216 // boost all decay products into the lab frame
217 TVector3 boostLab( pparent->Px() / pparent->E(),
218 pparent->Py() / pparent->E(),
219 pparent->Pz() / pparent->E());
220 fProducts[0].Boost(boostLab);
221 fProducts[1].Boost(boostLab);
222 fProducts[2].Boost(boostLab);
228 Int_t AliOmegaDalitz::ImportParticles(TClonesArray *particles)
230 // Import TParticles in array particles
233 TClonesArray &clonesParticles = *particles;
235 Int_t pdg [3] = {kElectron, -kElectron, kPi0};
236 if ( fProducts[1].M() > 0.001 )
241 Int_t parent[3] = {0, 0, -1};
242 Int_t d1 [3] = {-1, -1, 1};
243 Int_t d2 [3] = {-1, -1, 2};
244 for (Int_t i = 2; i > -1; i--) {
245 Double_t px = fProducts[i].Px();
246 Double_t py = fProducts[i].Py();
247 Double_t pz = fProducts[i].Pz();
248 Double_t e = fProducts[i].E();
250 new(clonesParticles[2 - i]) TParticle(pdg[i], 1, parent[i], -1, d1[i], d2[i], px, py, pz, e, 0., 0., 0., 0.);
256 void AliOmegaDalitz::
257 Rot(Double_t pin[3], Double_t pout[3], Double_t costheta, Double_t sintheta,
258 Double_t cosphi, Double_t sinphi) const
261 pout[0] = pin[0]*costheta*cosphi-pin[1]*sinphi+pin[2]*sintheta*cosphi;
262 pout[1] = pin[0]*costheta*sinphi+pin[1]*cosphi+pin[2]*sintheta*sinphi;
263 pout[2] = -1.0 * pin[0] * sintheta + pin[2] * costheta;
267 void AliOmegaDalitz::Decay(TClonesArray* array)
270 // Replace all omega dalitz decays with the correct matrix element decays
272 printf("-->Correcting Dalitz \n");
273 Int_t nt = array->GetEntriesFast();
275 for (Int_t i = 0; i < nt; i++) {
276 TParticle* part = (TParticle*) (array->At(i));
277 if (part->GetPdgCode() != 223) continue;
279 Int_t fd = part->GetFirstDaughter() - 1;
280 Int_t ld = part->GetLastDaughter() - 1;
282 if (fd < 0) continue;
283 if ((ld - fd) != 2) continue;
285 for (Int_t j = 0; j < 3; j++) dp[j] = (TParticle*) (array->At(fd+j));
286 if ((dp[0]->GetPdgCode() != kPi0) ||
287 ((TMath::Abs(dp[1]->GetPdgCode()) != kElectron) &&
288 (TMath::Abs(dp[1]->GetPdgCode()) != kMuonPlus))) continue;
289 TLorentzVector omega(part->Px(), part->Py(), part->Pz(), part->Energy());
290 if ( TMath::Abs(dp[1]->GetPdgCode()) == kElectron )
291 Decay(kElectron, &omega);
293 Decay(kMuonPlus, &omega);
294 for (Int_t j = 0; j < 3; j++) dp[j]->SetMomentum(fProducts[2-j]);
295 printf("original %13.3f %13.3f %13.3f %13.3f \n",
296 part->Px(), part->Py(), part->Pz(), part->Energy());
297 printf("new %13.3f %13.3f %13.3f %13.3f \n",
298 dp[0]->Px()+dp[1]->Px()+dp[2]->Px(),
299 dp[0]->Py()+dp[1]->Py()+dp[2]->Py(),
300 dp[0]->Pz()+dp[1]->Pz()+dp[2]->Pz(),
301 dp[0]->Energy()+dp[1]->Energy()+dp[2]->Energy());
306 AliOmegaDalitz& AliOmegaDalitz::operator=(const AliOmegaDalitz& rhs)
308 // Assignment operator
313 void AliOmegaDalitz::Copy(TObject&) const
318 Fatal("Copy","Not implemented!\n");
321 AliOmegaDalitz::AliOmegaDalitz(const AliOmegaDalitz &dalitz)