ATLAS Offline Software
Loading...
Searching...
No Matches
ParticleDecayer.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// SUMMARY: This code implements a "particle decayer" to allow us to augment the standard
6// "particle gun" with correlated particles, such as those from a decay. We assume that
7// some previous code has generated a scalar or a boson according to some preferred angular and momentum
8// distributions using a "geantino" (PDG ID 999). ParticleDecayer will then find this input
9// geantino assigns a given mass and (using the user defined parameters provided) select a decay vertex location where
10// the "geantino" is forced to decay into a particle/antiparticle pair.
11// The processes that the user can generate are
12//
13// s -> gd + gd
14// gd -> l + l
15//
16// or
17//
18// gd -> l + l
19//
20//
21// See the jobOption in share folder for more information.
22//
23
25#include "HepPDT/ParticleDataTable.hh"
27#include "CLHEP/Random/RandomEngine.h"
28#include <cmath> //M_PI
29
31//function to generate a lifetime (actually decay-length) according to the proper lifetime of the particle particle
32double ParticleDecayer::rnd_ExpLifetime(CLHEP::HepRandomEngine* engine, double ct) {
33 double r = engine->flat(); //< Return random num in [0,1]
34 return ((-ct)*std::log(1.-r));
35}
36
38//function to generate a uniform distribution in (a,b)
39double ParticleDecayer::rnd_DoubleRange(CLHEP::HepRandomEngine* engine, double a, double b) {
40 double r = engine->flat(); //< Return random num in [0,1]
41 return a + r*(b-a);
42}
43
44
46//function to generate a cos theta valua according to angular distribution
47double ParticleDecayer::cosgen(CLHEP::HepRandomEngine* engine, int itype){
48 double x,fx,hit;
49 x=0;
50
51 if(itype==0){ // flat
52 x=rnd_DoubleRange(engine,-1.,1.);
53 return x;
54 }
55
56 if(itype==1){ // f(cost)=1.5*cost^2
57 hit=1.5;
58 fx=0;
59 while(hit>fx){
60 x=rnd_DoubleRange(engine,-1.,1.);
61 hit=rnd_DoubleRange(engine,0.,1.5);
62 fx=1.5*x*x;
63 }
64 return x;
65 }
66
67 if(itype==2){ // f(cost)=0.375(1+cost^2)
68 hit=0.75;
69 fx=0;
70 while(hit>fx){
71 x=rnd_DoubleRange(engine,-1.,1.);
72 hit=rnd_DoubleRange(engine,0.,0.75);
73 fx=0.375*(1+x*x);
74 }
75 return x;
76 }
77
78 if(itype==3){ // f(cost)=0.75(1-cost^2)
79 hit=1.5;
80 fx=0;
81 while(hit>fx){
82 x=rnd_DoubleRange(engine,-1.,1.);
83 hit=rnd_DoubleRange(engine,0.,1.5);
84 fx=0.75*(1-x*x);
85 }
86 return x;
87 }
88
89 return x;
90
91}
92
94StatusCode ParticleDecayer::changeMass( HepMC::GenParticlePtr genpart, double newMass )
95{
96 double e = genpart->momentum().e();
97 double theta = genpart->momentum().theta();
98 double phi = genpart->momentum().phi();
99 //Sanity check
100 double p2 = e*e - newMass*newMass;
101 if ( p2 < 0 ) {
102 ATH_MSG_FATAL("ParticleDecayer::fillEvt: -- you have generated a tachyon!");
103 return StatusCode::FAILURE;
104 }
105 //At this point, we have e, theta, and phi. Put them together to get the four-momentum.
106 double p = std::sqrt(p2);
107 double px = p*std::sin(theta)*std::cos(phi);
108 double py = p*std::sin(theta)*std::sin(phi);
109 double pz = p*std::cos(theta);
110 //Fill the four-momentum
111 const CLHEP::HepLorentzVector updatedLV(px,py,pz,e);
112 genpart->set_momentum(HepMC::FourVector(updatedLV.x(),updatedLV.y(),updatedLV.z(),updatedLV.e()));
113 genpart->set_generated_mass(newMass);
114 return StatusCode::SUCCESS;
115}
116
118StatusCode ParticleDecayer::setDecayPosition( CLHEP::HepRandomEngine* engine, HepMC::GenParticlePtr genpart, HepMC::GenEvent* event, bool doScalarDecay )
119{
120 HepMC::GenVertexPtr vtxp = genpart->production_vertex();
121 if(!vtxp) {
122 ATH_MSG_FATAL("ParticleDecayer::fillEvt: -- no production vertex position found!");
123 return StatusCode::FAILURE;
124 }
125 if ( doScalarDecay ) // scalar decay is prompt. decay vtx == scalar production vertex
126 {
128 end_vtx->set_position(vtxp->position());
129 end_vtx->add_particle_in(std::move(genpart));
130 event->add_vertex(std::move(end_vtx));
131 return StatusCode::SUCCESS;
132 }
133
134 //lifetime handling
135 double gamma = genpart->momentum().e()/genpart->momentum().m();
136 double ctau = 0.;
137 double theta = genpart->momentum().theta();
138
139 // Make sure theta is between 0 and pi
140 while( theta > M_PI )
141 {
142 theta -= M_PI;
143 }
144 while( theta < 0 )
145 {
146 theta += M_PI;
147 }
148
150 {
151 if(m_expDecayDoVariableLifetime) // Variable proper lifetime distribution, such that fixed fraction of events decays within detector, independent of boost
152 {
153 double distanceToEdge = -999.;
154 if ( theta < m_thetaEndCapBarrel || theta > ( M_PI - m_thetaEndCapBarrel) ) // Particle escapes through endcap
155 {
156 distanceToEdge = std::abs(m_endCapDistance/std::cos(theta));
157 }
158 else // Particle escapes through barrel
159 {
160 distanceToEdge = m_barrelRadius/std::sin(theta);
161 }
162 if ( gamma < 1. )
163 {
164 ATH_MSG_FATAL("Gamma factor cannot be smaller than 1" );
165 return StatusCode::FAILURE;
166 }
167 double Limit = distanceToEdge / gamma; // ctau is enhanced by factor gamma in lab frame
168 double lambda = -1.*Limit/std::log(1. - m_expDecayFractionToKeep);
169 ctau = rnd_ExpLifetime(engine,lambda);
171 {
172 while ( ctau > Limit ) // If decay is outside detector, let's sample again (Produces a truncated exponential)
173 {
174 ctau = rnd_ExpLifetime(engine, lambda);
175 }
176 }
177 }
178 else if( m_particleLifeTime > 0 ) // Fixed proper lifetime distribution. Note that if m_particleLifeTime <= 0. ctau remains 0
179 {
180 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- set lifetime of the dark photon according to ctau = " << m_particleLifeTime);
181 ctau = rnd_ExpLifetime(engine,m_particleLifeTime);
182 }
183 }else if (m_doUniformDecay)
184 {
185 double decayRadius = -999.;
186 if ( theta < m_thetaEndCapBarrel || theta > ( M_PI - m_thetaEndCapBarrel) ) // Particle escapes through endcap
187 {
188 double outerLength = std::abs(m_endCapDistance/std::cos(theta));
189 double outerRadius = outerLength*std::sin(theta);
190 decayRadius = rnd_DoubleRange(engine,0., std::min(outerRadius, std::abs(m_barrelRadius)) );
191 }else // Particle escapes through barrel
192 {
193 decayRadius = rnd_DoubleRange(engine,0., std::abs(m_barrelRadius));
194 }
195
196 double decayLength = decayRadius/std::sin(theta);
197 ctau = decayLength/gamma;
198 }else
199 {
200 ATH_MSG_FATAL("have to pick uniform or exponential decay distance");
201 return StatusCode::FAILURE;
202 }
203
204 double ctg = gamma * ctau;
205 double px = genpart->momentum().px();
206 double py = genpart->momentum().py();
207 double pz = genpart->momentum().pz();
208 double p = std::sqrt(px*px + py*py + pz*pz);
209
210 const CLHEP::HepLorentzVector posLV(((ctg*px/p)+(vtxp->position().x())), ((ctg*py/p)+(vtxp->position().y())), ((ctg*pz/p)+(vtxp->position().z())), ((ctg)+(vtxp->position().t())));
211
212 //set the decay vertex position of the particle
213 //Create a HepMC vertex at the decay position of the particle
214 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- set the decay vertex");
216 end_vtx->set_position(HepMC::FourVector(posLV.x(),posLV.y(),posLV.z(),posLV.t()));
217 end_vtx->add_particle_in(std::move(genpart));
218 event->add_vertex(std::move(end_vtx));
219 return StatusCode::SUCCESS;
220}
221
223ParticleDecayer::ParticleDecayer(const std::string& name, ISvcLocator* pSvcLocator) :
224 GenModule(name, pSvcLocator),
225 m_truthParticleContainerName("GEN_EVENT"),
227{
228 declareProperty("McEventCollection", m_truthParticleContainerName);
229 declareProperty("LJType", m_LJType = 1); //1 = one dark photon per LeptonJet, 2 = two dark photons per LeptonJet
230 declareProperty("ScalarPDGID", m_scalarPDGID); //new PDG ID of the scalar
231 declareProperty("ScalarMass", m_scalarMass); //new mass of the scalar in MeV
232 declareProperty("ParticleID", m_particleID = 999); //PDG ID of the geantino
233 declareProperty("ParticleMass", m_particleMass = 400); //new mass of the dark photon in MeV (default 400 MeV)
234 declareProperty("ParticleLifeTime", m_particleLifeTime = 0); //ctau of the dark photon in mm (default prompt)
235 declareProperty("ParticlePolarization", m_particlePolarization = 0); //polarization of the dark photon (default 0 : isotropic decay, -1 : transverse, 1 : longitudinal)
236 declareProperty("OppositePolarization", m_oppositePolarization = false);//In case of LJType == 2 and opposite polarization for the two dark photons in the event
237 declareProperty("ParticlePDGID", m_particlePDGID = 700022); //new PDG ID of the dark photon (default 700022)
238 declareProperty("DecayOnlyStable", m_decayOnlyStable = false);
239 declareProperty("DecayBRElectrons", m_BRElectron = 1.); //BR of dark photon decay to electrons
240 declareProperty("DecayBRMuons", m_BRMuon = 0.); //BR of dark photon decay to muons
241 declareProperty("DecayBRPions", m_BRPion = 0.); //BR of dark photon decay to pions
242
243 //Choose between the different types of displaced decays
244 declareProperty("DoUniformDecay", m_doUniformDecay = false);
245 declareProperty("DoExponentialDecay", m_doExponentialDecay = true );
246 declareProperty("ExpDecayDoVariableLifetime", m_expDecayDoVariableLifetime = false );
247 declareProperty("ExpDecayPercentageToKeep", m_expDecayFractionToKeep = 0.8 );
248 declareProperty("ExpDecayDoTruncateLongDecays", m_expDecayDoTruncateLongDecays = false );
249
250 //Dimensions within which to decay in case of non-prompt events
251 declareProperty("BarrelRadius", m_barrelRadius = 10.e3 ); // outer limit for decay radius
252 declareProperty("EndCapDistance", m_endCapDistance = 15.e3 ); // outer limit along z for decays in endcap
253 declareProperty("ThetaEndCapBarrel", m_thetaEndCapBarrel = 0.439067982); // theta if where to switch between barrel and endcap (default: eta = 1.5)
254}
255
257
259{
260 // Check here if the tool has been configured correctly
262 {
263 ATH_MSG_FATAL("Cannot configure exponential and uniform decay at the same time.");
264 return StatusCode::FAILURE;
265 }
266
267 double TOTBR = m_BRElectron + m_BRMuon + m_BRPion;
268 if (TOTBR>1.0000001) {
269 ATH_MSG_FATAL("ParticleDecayer::genInitialize: -- Branching Ratio sum is larger than 1!! Please check the values in your jobOption.");
270 ATH_MSG_FATAL("ParticleDecayer::genInitialize: -- Branching Ratio Electrons " << m_BRElectron);
271 ATH_MSG_FATAL("ParticleDecayer::genInitialize: -- Branching Ratio Muons " << m_BRMuon);
272 ATH_MSG_FATAL("ParticleDecayer::genInitialize: -- Branching Ratio Pions " << m_BRPion);
273 ATH_MSG_FATAL("ParticleDecayer::genInitialize: -- Total Branching Ratio " << TOTBR);
274 return StatusCode::FAILURE;
275 }
276
277 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- daughters randomly produced: " << m_BRElectron << "% e, " << m_BRMuon << "% mu, " << m_BRPion << "% pi");
278
279 if( m_LJType != 1 && m_LJType != 2 )
280 {
281 ATH_MSG_FATAL("LJType is set to " << m_LJType);
282 ATH_MSG_FATAL("LJType has to be set to either 1 or 2");
283 return StatusCode::FAILURE;
284 }
285
286
287 return StatusCode::SUCCESS;
288}
289
291StatusCode ParticleDecayer::fillEvt(HepMC::GenEvent* event) {
292
293 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- EventCounter: " << m_eventCounter);
294
295 const EventContext& ctx = Gaudi::Hive::currentContext();
296 CLHEP::HepRandomEngine* engine = this->getRandomEngine("ParticleDecayer", ctx);
297
298 StatusCode status = StatusCode::SUCCESS;
299
300 // Extract the event from the TES
301 McEventCollection* mcEvtColl = nullptr;
303 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- Found an McEventCollection for ParticleDecayer");
304 } else {
305 ATH_MSG_FATAL("TES is empty!!!!");
306 status = StatusCode::FAILURE;
307 return status;
308 }
309
310 event = mcEvtColl->back();
311
312 if (event == NULL) {
313 ATH_MSG_FATAL("ParticleDecayer::fillEvt: -- McEvent was not successfully created");
314 status = StatusCode::FAILURE;
315 return status;
316 }
317
318 for ( auto genpart : *event) {
319
320 const bool hasRightPDG = (genpart->pdg_id() == m_particleID);
321 const bool isStable = (genpart->status() == 1);
322 const bool hasEndVtx = (genpart->end_vertex() != nullptr);
323
324 if (!hasRightPDG) {
325 continue;
326 }
327
328 if (m_decayOnlyStable) {
329 if (!isStable || hasEndVtx) {
330 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- skip particle with PDG ID = "
331 << genpart->pdg_id()
332 << " status = " << genpart->status()
333 << " end_vertex = " << (hasEndVtx ? "yes" : "no"));
334 continue;
335 }
336 }
337
339 // only one dark photon per LeptonJet //
341 if (m_LJType == 1) {
342 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- only one dark photon per LeptonJet");
343 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- found MC particle with PDG ID = " << genpart->pdg_id());
344 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- assign the new mass of the dark photon, m = " << m_particleMass);
345
346 //Change the mass of the parent particle ( set by user input + command )
347 //Changes the magnitude of the spatial part of the 4-momentum such that the new 4-momentum has the desired inv mass
348 CHECK( changeMass( genpart, m_particleMass ) );
349
350 //Add decay position to the event
351 CHECK( setDecayPosition( engine, genpart, event ) );
352
353 //assign the new PDG_ID of the particle
354 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- set the new PDG ID = " << m_particlePDGID);
355 genpart->set_pdg_id(m_particlePDGID);
356
357 //set the new status (decayed) of the particle
358 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- set the new status = 2");
359 genpart->set_status(2);
360
361 //set the new momentum of the particle
362 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- set the new momentum");
363
365 CHECK( DFTwoBodyDecay( engine, std::move(genpart), m_particlePolarization ) );
366
367 }else if (m_LJType == 2) {
369 // two dark photons per LeptonJet //
371
372 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- two dark photons per LeptonJet");
373 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- found MC particle with PDG ID = " << genpart->pdg_id());
374 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- assign the new mass of the dark scalar, m = " << m_scalarMass);
375
376 //Get the mass of the parent particle ( set by user input + command )
377 //Change the mass of the parent particle ( set by user input + command )
378 //Changes the magnitude of the spatial part of the 4-momentum such that the new 4-momentum has the desired inv mass
379 CHECK( changeMass( genpart, m_scalarMass ) );
380
381 //set the new PDG_ID of the scalar
382 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- set the new PDG ID = " << m_scalarPDGID);
383 genpart->set_pdg_id(m_scalarPDGID);
384
385 //set the new status (decayed) of the scalar
386 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- set the new status = 2");
387 genpart->set_status(2);
388
389 //set the decay vertex position of the scalar
390 //Create a HepMC vertex at the decay position of the scalar
391 CHECK( setDecayPosition( engine, genpart, event, true ) );
392
394 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- allow the two-body decay of the dark scalar to dark photons...");
395 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- dark photon has PDG ID = " << m_particlePDGID);
396
397 std::vector<CLHEP::HepLorentzVector> darkPhotonLVs;
398 CHECK( getDecayProducts( engine, CLHEP::HepLorentzVector( genpart->momentum().px(), genpart->momentum().py(), genpart->momentum().pz(), genpart->momentum().e() ),
400 darkPhotonLVs) );
401
402 //add dark photon 1
403 auto v0=darkPhotonLVs.at(0).vect();
404 addParticle( genpart->end_vertex(), m_particlePDGID, HepMC::FourVector(v0.x(),v0.y(),v0.z(),0.0) , 2);
405 //add dark photon 2
406 auto v1=darkPhotonLVs.at(1).vect();
407 addParticle( genpart->end_vertex(), m_particlePDGID, HepMC::FourVector(v1.x(),v1.y(),v1.z(),0.0), 2);
408
409 //lifetime handling of the dark photons
410 int polarizationSwitch = 1;
411#ifdef HEPMC3
412 const std::vector<HepMC::GenParticlePtr>& particlesOut = genpart->end_vertex()->particles_out();
413 std::vector<HepMC::GenParticlePtr>::const_iterator pItBegin = particlesOut.begin();
414 std::vector<HepMC::GenParticlePtr>::const_iterator pItEnd = particlesOut.end();
415#else
416 HepMC::GenVertex::particles_out_const_iterator pItBegin = genpart->end_vertex()->particles_out_const_begin();
417 HepMC::GenVertex::particles_out_const_iterator pItEnd = genpart->end_vertex()->particles_out_const_end();
418#endif
419 for ( auto pIt=pItBegin ; pIt != pItEnd; ++pIt )
420 {
421
422 //Add decay position to the event
423 CHECK( setDecayPosition( engine, *pIt, event ) );
424 //And perform two-body decay
425 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- Now allow the two-body decay of the dark photons");
426 CHECK( DFTwoBodyDecay( engine, *pIt, polarizationSwitch*m_particlePolarization ) );
428 {
429 polarizationSwitch = -polarizationSwitch;
430 }
431 }
432
433 }else
434 {
435 ATH_MSG_FATAL("LJType set to " << m_LJType );
436 ATH_MSG_FATAL("LJType must be 1 or 2" );
437 return StatusCode::FAILURE;
438 }
439 }
440
441 //set the event number
442 event->set_event_number(m_eventCounter);
443
444 //Add a unit entry to the event weight vector if it's currently empty
445 if (event->weights().empty()) {
446 event->weights().push_back(1.);
447 }
448
449 //increment the event counter
451
452 return status;
453}
454
455
457void ParticleDecayer::addParticle(HepMC::GenVertexPtr prod_vtx, int pdg, HepMC::FourVector momentum, int statusCode) {
458
459 double mass = 0.;
460 if( pdg == m_particlePDGID)
461 {
462 mass = m_particleMass;
463 }else
464 {
465 mass = getParticleMass(pdg);
466 }
467 double energy=std::sqrt(std::pow(momentum.x(),2)+std::pow(momentum.y(),2)+std::pow(momentum.z(),2)+mass*mass);
468HepMC::GenParticlePtr aParticle = HepMC::newGenParticlePtr (HepMC::FourVector(momentum.x(), momentum.y(), momentum.z(), energy),
469 pdg, statusCode);
470
471 prod_vtx->add_particle_out(std::move(aParticle));
472}
473
474
477 SmartIF<IPartPropSvc> partPropSvc(Gaudi::svcLocator()->service("PartPropSvc"));
478 if (!partPropSvc) throw GaudiException("PartPropSvc error", "I_ParticleDecayer", StatusCode::FAILURE);
480
481 const HepPDT::ParticleData* particle = m_particleTable->particle(HepPDT::ParticleID(std::abs(pid)));
482 return particle->mass().value();
483}
484
485
487StatusCode ParticleDecayer::DFTwoBodyDecay( CLHEP::HepRandomEngine* engine, HepMC::GenParticlePtr genpart, int Polarization ) {
488
489
490 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- allow the two-body decay of the dark photon...");
491
492 //create heplorentzvector from genpart's fourMomentum, as it is more useful
493 CLHEP::HepLorentzVector boostDF( genpart->momentum().px(), genpart->momentum().py(), genpart->momentum().pz(), genpart->momentum().e() );
494
495 //Given branching fractions, pick decay mode
496 int ModeOfDecay;
497 double unif = rnd_DoubleRange(engine, 0.,1.);
498 if (unif<=m_BRElectron) {
499 ModeOfDecay = 11;
500 } else if (unif<=(m_BRElectron+m_BRMuon) && unif>m_BRElectron) {
501 ModeOfDecay = 13;
502 } else {
503 ModeOfDecay = 211;
504 }
505 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- decayMode = " << ModeOfDecay);
506
507 //Now that we have a decay mode. get the associated particle mass
508 double decayPartMass = getParticleMass(ModeOfDecay);
509
510 //Choose how to decay
511 //angular distribution handling, see pag.6 of http://arxiv.org/pdf/hep-ph/0605296v2.pdf
512 int type = 0; //isotropic by default
513 if(Polarization==0) {
514 type = 0; //isotropic
515 } else if(Polarization==-1 && (ModeOfDecay==11 || ModeOfDecay==13)) {
516 type = 2; //transverse polarization of the dark photon + fermions
517 } else if(Polarization==1 && (ModeOfDecay==11 || ModeOfDecay==13)) {
518 type = 3; //longitudinal polarization of the dark photon + fermions
519 } else if(Polarization==-1 && ModeOfDecay==211) {
520 type = 3; //transverse polarization of the dark photon + scalars
521 } else if(Polarization==1 && ModeOfDecay==211) {
522 type = 1; //longitudinal polarization of the dark photon + scalars
523 } else {
524 ATH_MSG_FATAL("ParticleDecayer::fillEvt: -- wrong polarization value... please check!");
525 return StatusCode::FAILURE;
526 }
527
528 //Now we are ready to decay the dark photon
529 std::vector<CLHEP::HepLorentzVector> daughterLVs;
530 CHECK( getDecayProducts( engine, boostDF, decayPartMass, daughterLVs, type) );
531
532 //Add the daughters to the pool file
533 ATH_MSG_DEBUG("ParticleDecayer::fillEvt: -- Add the daughters to the pool file");
534 auto end_vtx = genpart->end_vertex();
535 auto v0=daughterLVs.at(0).vect();
536 addParticle(end_vtx, ModeOfDecay, HepMC::FourVector(v0.x(),v0.y(),v0.z(),0.0), 1);
537 auto v1=daughterLVs.at(1).vect();
538 addParticle(std::move(end_vtx), -ModeOfDecay, HepMC::FourVector(v1.x(),v1.y(),v1.z(),0.0), 1);
539
540 return StatusCode::SUCCESS;
541}
542
543
544StatusCode ParticleDecayer::getDecayProducts( CLHEP::HepRandomEngine* engine, CLHEP::HepLorentzVector parentLV,
545 double decayPartMass,
546 std::vector<CLHEP::HepLorentzVector>& daughterLVs,
547 int decayType)
548{
549 double parentMass = parentLV.m();
550 CLHEP::Hep3Vector boostVec = parentLV.boostVector();
551
552 if( decayPartMass > parentMass/2.)
553 {
554 ATH_MSG_FATAL("Decay particle has more than half the mass of parent.");
555 return StatusCode::FAILURE;
556 }
557 //Get the angles in the rest frame
558 double phi_rf = rnd_DoubleRange(engine, -M_PI, M_PI);
559 double ct_rf = cosgen(engine, decayType);
560 double theta_rf = std::acos(ct_rf);
561
562 //construct p1 particle momentum in rest-frame (_rf)
563 double p1_rf = std::sqrt(parentMass*parentMass/4. - decayPartMass*decayPartMass);
564 double px_rf = p1_rf*std::cos(phi_rf)*std::sin(theta_rf);
565 double py_rf = p1_rf*std::sin(phi_rf)*std::sin(theta_rf);
566 double pz_rf = p1_rf*ct_rf;
567 CLHEP::HepLorentzVector hlv1( px_rf, py_rf, pz_rf, parentMass/2.);
568 CLHEP::HepLorentzVector hlv2( -px_rf, -py_rf, -pz_rf, parentMass/2.);
569
570 //Transform from the rotated frame to the (x,y,z) frame.
571 hlv1.rotateUz((parentLV.vect()).unit());
572 hlv2.rotateUz((parentLV.vect()).unit());
573
574 //Now boost back to the lab-frame
575 hlv1.boost(boostVec);
576 hlv2.boost(boostVec);
577
578 daughterLVs.push_back(hlv1);
579 daughterLVs.push_back(hlv2);
580
581 return StatusCode::SUCCESS;
582}
#define M_PI
Scalar phi() const
phi method
Scalar theta() const
theta method
#define ATH_MSG_FATAL(x)
#define ATH_MSG_DEBUG(x)
#define CHECK(...)
Evaluate an expression and check for errors.
static Double_t a
#define x
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
const T * back() const
Access the last element in the collection as an rvalue.
const ServiceHandle< IPartPropSvc > partPropSvc() const
Access the particle property service.
Definition GenBase.h:113
GenModule(const std::string &name, ISvcLocator *pSvcLocator)
Constructor.
Definition GenModule.cxx:14
CLHEP::HepRandomEngine * getRandomEngine(const std::string &streamName, const EventContext &ctx) const
Definition GenModule.cxx:34
This defines the McEventCollection, which is really just an ObjectVector of McEvent objectsFile: Gene...
double getParticleMass(int pdgID)
StatusCode changeMass(HepMC::GenParticlePtr, double)
bool m_expDecayDoVariableLifetime
HepPDT::ParticleDataTable * m_particleTable
StatusCode DFTwoBodyDecay(CLHEP::HepRandomEngine *engine, HepMC::GenParticlePtr, int)
StatusCode getDecayProducts(CLHEP::HepRandomEngine *engine, CLHEP::HepLorentzVector, double, std::vector< CLHEP::HepLorentzVector > &, int decayType=0)
StatusCode setDecayPosition(CLHEP::HepRandomEngine *engine, HepMC::GenParticlePtr, HepMC::GenEvent *, bool doScalarDecay=false)
double rnd_DoubleRange(CLHEP::HepRandomEngine *engine, double a, double b)
std::string m_truthParticleContainerName
void addParticle(HepMC::GenVertexPtr, int pdg, HepMC::FourVector, int statusCode)
ParticleDecayer(const std::string &name, ISvcLocator *pSvcLocator)
double cosgen(CLHEP::HepRandomEngine *engine, int itype)
double m_expDecayFractionToKeep
StatusCode fillEvt(HepMC::GenEvent *)
For filling the HepMC event object.
bool m_expDecayDoTruncateLongDecays
StatusCode genInitialize()
For initializing the generator, if required.
double rnd_ExpLifetime(CLHEP::HepRandomEngine *engine, double ct)
int r
Definition globals.cxx:22
bool contains(const std::string &s, const std::string &regx)
does a string contain the substring
Definition hcg.cxx:114
HepMC::GenVertex * GenVertexPtr
Definition GenVertex.h:59
GenVertexPtr newGenVertexPtr(const HepMC::FourVector &pos=HepMC::FourVector(0.0, 0.0, 0.0, 0.0), const int i=0)
Definition GenVertex.h:64
GenParticlePtr newGenParticlePtr(const HepMC::FourVector &mom=HepMC::FourVector(0.0, 0.0, 0.0, 0.0), int pid=0, int status=0)
Definition GenParticle.h:39
GenParticle * GenParticlePtr
Definition GenParticle.h:37