Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MuonicAtomDecay.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//---------------------------------------------------------------------
28//
29// GEANT4 Class
30//
31// File name: G4MuonicAtomDecay
32//
33// 20170522 K L Genser first implementation based on code by
34// V.Ivantchenko & G4HadronicProcess & G4Decay
35//
36// Class Description:
37//
38// MuonicAtom Process where Muon either decays in orbit or is captured by the nucleus
39//
40// Modifications:
41//
42//
43//------------------------------------------------------------------------
44
45#include "G4MuonicAtomDecay.hh"
48#include "G4Nucleus.hh"
49#include "G4ProcessManager.hh"
50#include "G4HadFinalState.hh"
51#include "G4HadProjectile.hh"
52#include "G4HadSecondary.hh"
53#include "G4ForceCondition.hh"
54#include "G4MuonicAtom.hh"
55#include "G4MuonicAtomHelper.hh"
56#include "G4VDecayChannel.hh"
57#include "G4DecayTable.hh"
58#include "G4DecayProducts.hh"
59#include "G4CascadeInterface.hh"
61#include "G4RandomDirection.hh"
63#include "G4SystemOfUnits.hh"
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
68 const G4String& name)
70 fMuMass(G4MuonMinus::MuonMinus()->GetPDGMass()),
71 cmptr(hiptr),
72 verboseLevel(0)
73{
74 // This is not a hadronic process; assume it is a kind of decay
75 enableAtRestDoIt = true;
76 enablePostStepDoIt = true; // it is a streach; fixme
77 theProcessSubType = 221; // (see G4DecayProcessType.hh) fixme
78 if (!cmptr) {
79 // cmptr = new G4CascadeInterface(); // Bertini - Pointer owned by InteractionRegistry
80 cmptr = new G4MuMinusCapturePrecompound(); // Precompound
81 }
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85
87//{delete theTotalResult;}
88{}
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91
93{
94 return ( a.GetParticleType() == "MuonicAtom" );
95}
96
97//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
98
99// void
100// G4MuonicAtomDecay::PreparePhysicsTable(const G4ParticleDefinition& p)
101// {
102// G4HadronicProcessStore::Instance()->RegisterParticleForExtraProcess(this,&p);
103// }
104
105// //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
106
107// void G4MuonicAtomDecay::BuildPhysicsTable(const G4ParticleDefinition& p)
108// {
109// G4HadronicProcessStore::Instance()->PrintInfo(&p);
110// }
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
113
115 const G4Track& aTrack, G4ForceCondition* condition)
116{
118 // check if this is the beginning of tracking
121 }
123}
124
125//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
126
129{
131 return DBL_MAX; // this will need to be changed in future; fixme
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
135
138{
139 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
140 G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
141 G4MuonicAtom* muatom = static_cast<G4MuonicAtom*>(aParticleDef);
142 G4double meanlife = muatom->GetPDGLifeTime();
143#ifdef G4VERBOSE
144 if (GetVerboseLevel()>1) {
145 G4cout << "mean life time: "<< meanlife/ns << "[ns]" << G4endl;
146 }
147#endif
148 return meanlife;
149}
150
151
152G4VParticleChange* G4MuonicAtomDecay::DecayIt(const G4Track& aTrack,
153 const G4Step&)
154{
155
156 // mainly based on G4HadronStoppingProcess & G4Decay
157 // if primary is not Alive then do nothing
158 theTotalResult.Clear(); // G4ParticleChange*
159 theTotalResult.Initialize(aTrack);
160 theTotalResult.ProposeWeight(aTrack.GetWeight());
161 if(aTrack.GetTrackStatus() != fAlive &&
162 aTrack.GetTrackStatus() != fStopButAlive) {
163 return &theTotalResult;
164 }
165
166 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
167 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
168 G4MuonicAtom const* muatom = static_cast<G4MuonicAtom const*>(aParticleDef);
169 G4Ions const* baseion = muatom->GetBaseIon();
170 G4int Z = baseion->GetAtomicNumber();
171 G4double Zd = Z;
172 G4double KEnergy = G4MuonicAtomHelper::GetKShellEnergy(Zd); // fixme check
173 G4HadProjectile thePro(aTrack); // G4HadProjectile, here the muonic atom
174 thePro.SetBoundEnergy(KEnergy);
175
176 G4ForceCondition* condition = nullptr; // it is unused in the following call anyway
177 G4double meanlife = GetMeanLifeTime(aTrack, condition);
178
179 G4HadFinalState* result = nullptr; // result before converting to G4VParticleChange*
180 // G4int nSecondaries = 0;
181 // save track time and start from zero time
182 // G4double time0 = aTrack.GetGlobalTime(); FillResult does it
183 // see G4Decay DecayIt
184 // see time0 down below
185 thePro.SetGlobalTime(0.0);
186
187 // do we need G4double fRemainderLifeTime; ???
188
189 G4double maDTime = theNumberOfInteractionLengthLeft*meanlife; //fixme check
190#ifdef G4VERBOSE
191 if (GetVerboseLevel()>1) {
192 G4cout << "G4MuonicAtomDecay::DecayIt time set to: "<< maDTime/ns << "[ns]" << G4endl;
193 }
194#endif
195
196 // decide on DIO or Capture
197
198 G4double lambdac = 1./muatom->GetDIOLifeTime();
199 G4double lambdad = 1./muatom->GetNCLifeTime();
200 G4double lambda = lambdac + lambdad;
201
202 if ( G4UniformRand()*lambda < lambdac) {
203 // if ( false ) { // force NC for testing
204
205 // DIO
206 // result = dmptr->ApplyYourself(thePro, *nucleus); // not quite the reaction;
207 // using G4PhaseSpaceDecayChannel
208
209#ifdef G4VERBOSE
210 if (GetVerboseLevel()>0) {
211 G4cout << "G4MuonicAtomDecay::DecayIt: selected DIO mode" << G4endl;
212 }
213#endif
214
215 // decay table; we use it only for the DIO which is more of a decay
216 // code mostly copied from G4Decay
217
218 G4DecayProducts* products = nullptr;
219 G4DecayTable *decaytable = aParticleDef->GetDecayTable();
220 G4VDecayChannel* decaychannel = nullptr;
221 G4double massParent = aParticle->GetMass();
222 decaychannel = decaytable->SelectADecayChannel(massParent);
223 if ( decaychannel ==0) {
224 // decay channel not found
226 ed << "Can not determine decay channel for "
227 << aParticleDef->GetParticleName() << G4endl
228 << " mass of dynamic particle: " << massParent/GeV << " (GEV)" << G4endl
229 << " dacay table has " << decaytable->entries() << " entries" << G4endl;
230 G4double checkedmass=massParent;
231 if (massParent < 0.) {
232 checkedmass=aParticleDef->GetPDGMass();
233 ed << "Using PDG mass ("<<checkedmass/GeV << "(GeV)) in IsOKWithParentMass" << G4endl;
234 }
235 for (G4int ic =0;ic <decaytable->entries();++ic) {
236 G4VDecayChannel * dc= decaytable->GetDecayChannel(ic);
237 ed << ic << ": BR " << dc->GetBR() << ", IsOK? "
238 << dc->IsOKWithParentMass(checkedmass)
239 << ", --> ";
240 G4int ndaughters=dc->GetNumberOfDaughters();
241 for (G4int id=0;id<ndaughters;++id) {
242 if (id>0) ed << " + "; // seperator, except for first
243 ed << dc->GetDaughterName(id);
244 }
245 ed << G4endl;
246 }
247 G4Exception("G4MuonicAtomDecay::DecayIt", "DECAY003", FatalException,ed);
248 } else {
249 // execute DecayIt()
250#ifdef G4VERBOSE
251 G4int temp = decaychannel->GetVerboseLevel();
252 if (GetVerboseLevel()>1) {
253 G4cout << "G4MuonicAtomDecay::DecayIt : selected decay channel addr:"
254 << decaychannel <<G4endl;
255 decaychannel->SetVerboseLevel(GetVerboseLevel());
256 }
257#endif
258 products = decaychannel->DecayIt(aParticle->GetMass());
259 if(!products) {
261 ed << "No products are generated for "
262 << aParticleDef->GetParticleName();
263 G4Exception("G4MuonicAtomDecay::DecayIt","DECAY003",FatalException,ed);
264 }
265#ifdef G4VERBOSE
266 if (GetVerboseLevel()>1) {
267 decaychannel->SetVerboseLevel(temp);
268 }
269 if (GetVerboseLevel()>2 && products) {
270 if (! products->IsChecked() ) products->DumpInfo();
271 }
272#endif
273 }
274
275 // get parent particle information ...................................
276 G4double ParentEnergy = aParticle->GetTotalEnergy();
277 G4double ParentMass = aParticle->GetMass();
278 if (ParentEnergy < ParentMass) {
279 if (GetVerboseLevel()>0) {
280 G4cout << "G4MuonicAtomDecay::DecayIt : Total Energy is less than its mass" << G4endl;
281 G4cout << " Particle: " << aParticle->GetDefinition()->GetParticleName();
282 G4cout << " Energy:" << ParentEnergy/MeV << "[MeV]";
283 G4cout << " Mass:" << ParentMass/MeV << "[MeV]";
284 G4cout << G4endl;
285 }
286 G4Exception( "G4MuonicAtomDecay::DecayIt ",
287 "DECAY102",JustWarning,
288 "Total Energy is less than its mass");
289 ParentEnergy = ParentMass;
290 }
291
292 G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
293
294 //boost all decay products to laboratory frame
295 G4double energyDeposit = 0.0;
296 G4double finalGlobalTime = aTrack.GetGlobalTime();
297 G4double finalLocalTime = aTrack.GetLocalTime();
298 if (aTrack.GetTrackStatus() == fStopButAlive ){
299 // AtRest case
300 finalGlobalTime += maDTime;
301 finalLocalTime += maDTime;
302 energyDeposit += aParticle->GetKineticEnergy();
303 } else {
304 // PostStep case
305 products->Boost( ParentEnergy, ParentDirection);
306 }
307 // G4ParticleChangeForDecay fParticleChangeForDecay; // is it equivalent to G4ParticleChange* theTotalResult;
308
309 //add products in theTotalResult
310 G4int numberOfSecondaries = products->entries();
311 theTotalResult.SetNumberOfSecondaries(numberOfSecondaries);
312#ifdef G4VERBOSE
313 if (GetVerboseLevel()>1) {
314 G4cout << "G4MuonicAtomDecay::DecayIt : Decay vertex :";
315 G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
316 G4cout << " X:" << (aTrack.GetPosition()).x() /cm << "[cm]";
317 G4cout << " Y:" << (aTrack.GetPosition()).y() /cm << "[cm]";
318 G4cout << " Z:" << (aTrack.GetPosition()).z() /cm << "[cm]";
319 G4cout << G4endl;
320 G4cout << "G4MuonicAtomDecay::DecayIt : decay products in Lab. Frame" << G4endl;
321 products->DumpInfo();
322 }
323#endif
324 G4int index;
325 G4ThreeVector currentPosition;
326 const G4TouchableHandle thand = aTrack.GetTouchableHandle();
327 for (index=0; index < numberOfSecondaries; index++)
328 {
329 // get current position of the track
330 currentPosition = aTrack.GetPosition();
331 // create a new track object
332 G4Track* secondary = new G4Track( products->PopProducts(),
333 finalGlobalTime ,
334 currentPosition );
335 // switch on good for tracking flag
336 secondary->SetGoodForTrackingFlag();
337 secondary->SetTouchableHandle(thand);
338 // add the secondary track in the List
339 theTotalResult.AddSecondary(secondary);
340 }
341 delete products;
342
343 // Kill the parent particle
344 theTotalResult.ProposeTrackStatus( fStopAndKill ) ;
345 theTotalResult.ProposeLocalEnergyDeposit(energyDeposit);
346 theTotalResult.ProposeLocalTime( finalLocalTime );
347
348 // Clear NumberOfInteractionLengthLeft
350
351 return &theTotalResult ;
352
353 } else { //either or
354
355 // nuclearCapture
356
357 // model
358 // need to be able to choose between preco or bertini; no good way to do it?
359 // hardcoded in the constructor for now
360
361#ifdef G4VERBOSE
362 if (GetVerboseLevel()>0) {
363 G4cout << "G4MuonicAtomDecay::DecayIt: selected NC mode" << G4endl;
364 }
365#endif
366
367 G4int A = baseion->GetAtomicMass();
368 // G4Nucleus* nucleus = GetTargetNucleusPointer(); // from G4HadronicProcess
369 G4Nucleus nucleus;
370 nucleus.SetParameters(A, Z);
371
372 // we define a local projectile here which will be the orbiting muon
373 // we shall assume it is at rest; fixme
374
375 // G4HadProjectile, here the muon
377 G4ThreeVector(0.,0.,0.)));
378 theMuPro.SetBoundEnergy(KEnergy);
379 theMuPro.SetGlobalTime(0.0);
380
381 G4int reentryCount = 0; // this may be in the model already; check fixme <---
382 do {
383 // sample final state
384 // nuclear interaction should keep G4HadFinalState object
385 // model should define time of each secondary particle
386 try {
387 result = cmptr->ApplyYourself(theMuPro, nucleus); // muon and muonic atom nucleus
388 ++reentryCount;
389 }
390 catch(G4HadronicException & aR) {
392 ed << "Call for " << cmptr->GetModelName() << G4endl;
393 ed << " Z= "
394 << nucleus.GetZ_asInt()
395 << " A= " << nucleus.GetA_asInt() << G4endl;
396 DumpState(aTrack,"ApplyYourself",ed);
397 ed << " ApplyYourself failed" << G4endl;
398 G4Exception("G4MuonicAtomDecay::DecayIt", "HAD_MAD_101",
399 FatalException, ed);
400 }
401
402 // Check the result for catastrophic energy non-conservation
403 // result = CheckResult(theMuPro, nucleus, result);
404
405 if(reentryCount>100) {
407 ed << "Call for " << cmptr->GetModelName() << G4endl;
408 ed << " Z= "
409 << nucleus.GetZ_asInt()
410 << " A= " << nucleus.GetA_asInt() << G4endl;
411 DumpState(aTrack,"ApplyYourself",ed);
412 ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
413 G4Exception("G4MuonicAtomDecay::DecayIt", "HAD_MAD_102",
414 FatalException, ed);
415 }
416 // Loop checking, 06-Aug-2015, Vladimir Ivanchenko
417 } while(!result);
418
419 // add delay time of capture (inter + intra)
420 G4int nsec = result->GetNumberOfSecondaries();
421 for(G4int i=0; i<nsec; ++i) {
422 G4HadSecondary* sec = result->GetSecondary(i);
423 G4double ctime = sec->GetTime();
424 sec->SetTime(maDTime + ctime); // we add time0 in the next stage
425#ifdef G4VERBOSE
426 if (GetVerboseLevel()>1) {
427 G4cout << "G4MuonicAtomDecay::DecayIt time set to: "
428 << (maDTime + ctime)/ns << "[ns]" << G4endl;
429 }
430#endif
431 }
432
433 FillResult(result,aTrack);
434
435 // delete result;// causes bad free check fixme; move to the class members?
436
438 return &theTotalResult;
439
440 }
441
442
443}
444
445//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
446
447void G4MuonicAtomDecay::ProcessDescription(std::ostream& outFile) const
448{
449 outFile << "MuonicAtom process where Muon decays in orbit or is captured by the nucleus." <<G4endl;
450}
451
452//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
453
454void G4MuonicAtomDecay::FillResult(G4HadFinalState * aR, const G4Track & aT)
455{
456 // based on G4HadronicProcess::FillResult
457
459
460 G4double rotation = CLHEP::twopi*G4UniformRand();
461 G4ThreeVector it(0., 0., 1.);
462
463 G4double efinal = aR->GetEnergyChange();
464 if(efinal < 0.0) { efinal = 0.0; }
465
466 // check status of primary
467 if(aR->GetStatusChange() == stopAndKill) {
468 theTotalResult.ProposeTrackStatus(fStopAndKill);
469 theTotalResult.ProposeEnergy( 0.0 );
470
471 // check its final energy
472 } else if(0.0 == efinal) {
473 theTotalResult.ProposeEnergy( 0.0 );
475 ->GetAtRestProcessVector()->size() > 0)
476 { theTotalResult.ProposeTrackStatus(fStopButAlive); }
477 else { theTotalResult.ProposeTrackStatus(fStopAndKill); } // check fixme
478
479 // primary is not killed apply rotation and Lorentz transformation
480 } else {
481 theTotalResult.ProposeTrackStatus(fAlive);
483 G4double newE = efinal + mass;
484 G4double newP = std::sqrt(efinal*(efinal + 2*mass));
485 G4ThreeVector newPV = newP*aR->GetMomentumChange();
486 G4LorentzVector newP4(newE, newPV);
487 newP4.rotate(rotation, it);
488 newP4 *= aR->GetTrafoToLab();
489 theTotalResult.ProposeMomentumDirection(newP4.vect().unit());
490 newE = newP4.e() - mass;
491#ifdef G4VERBOSE
492 if (GetVerboseLevel()>1 && newE <= 0.0) {
494 DumpState(aT,"Primary has zero energy after interaction",ed);
495 G4Exception("G4MuonicAtomDecay::FillResults", "HAD_MAD_103", JustWarning, ed);
496 }
497#endif
498 if(newE < 0.0) { newE = 0.0; }
499 theTotalResult.ProposeEnergy( newE );
500 }
501 //G4cout << "FillResult: Efinal= " << efinal << " status= "
502 // << theTotalResult.GetTrackStatus()
503 // << " fKill= " << fStopAndKill << G4endl;
504
505 // check secondaries: apply rotation and Lorentz transformation
506 G4int nSec = aR->GetNumberOfSecondaries();
507 theTotalResult.SetNumberOfSecondaries(nSec);
508 G4double weight = aT.GetWeight();
509
510 if (nSec > 0) {
511 G4double time0 = aT.GetGlobalTime();
512 for (G4int i = 0; i < nSec; ++i) {
514 theM.rotate(rotation, it);
515 theM *= aR->GetTrafoToLab();
516 aR->GetSecondary(i)->GetParticle()->Set4Momentum(theM);
517
518 // time of interaction starts from zero
519 G4double time = aR->GetSecondary(i)->GetTime();
520 if (time < 0.0) { time = 0.0; }
521
522 // take into account global time
523 time += time0;
524
525 G4Track* track = new G4Track(aR->GetSecondary(i)->GetParticle(),
526 time, aT.GetPosition());
528 G4double newWeight = weight*aR->GetSecondary(i)->GetWeight();
529 // G4cout << "#### ParticleDebug "
530 // <<GetProcessName()<<" "
531 //<<aR->GetSecondary(i)->GetParticle()->GetDefinition()->GetParticleName()<<" "
532 // <<aScaleFactor<<" "
533 // <<XBiasSurvivalProbability()<<" "
534 // <<XBiasSecondaryWeight()<<" "
535 // <<aT.GetWeight()<<" "
536 // <<aR->GetSecondary(i)->GetWeight()<<" "
537 // <<aR->GetSecondary(i)->GetParticle()->Get4Momentum()<<" "
538 // <<G4endl;
539 track->SetWeight(newWeight);
541 theTotalResult.AddSecondary(track);
542#ifdef G4VERBOSE
543 if (GetVerboseLevel()>1) {
544 G4double e = track->GetKineticEnergy();
545 if (e <= 0.0) {
547 DumpState(aT,"Secondary has zero energy",ed);
548 ed << "Secondary " << track->GetDefinition()->GetParticleName()
549 << G4endl;
550 G4Exception("G4MuonicAtomDecay::FillResults", "HAD_MAD_103",
551 JustWarning,ed);
552 }
553 }
554#endif
555 }
556 }
557 aR->Clear();
558}
559
560//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
561
562void G4MuonicAtomDecay::DumpState(const G4Track& aTrack,
563 const G4String& method,
565{
566 ed << "Unrecoverable error in the method " << method << " of "
567 << GetProcessName() << G4endl;
568 ed << "TrackID= "<< aTrack.GetTrackID() << " ParentID= "
569 << aTrack.GetParentID()
570 << " " << aTrack.GetParticleDefinition()->GetParticleName()
571 << G4endl;
572 ed << "Ekin(GeV)= " << aTrack.GetKineticEnergy()/CLHEP::GeV
573 << "; direction= " << aTrack.GetMomentumDirection() << G4endl;
574 ed << "Position(mm)= " << aTrack.GetPosition()/CLHEP::mm << ";";
575
576 if (aTrack.GetMaterial()) {
577 ed << " material " << aTrack.GetMaterial()->GetName();
578 }
579 ed << G4endl;
580
581 if (aTrack.GetVolume()) {
582 ed << "PhysicalVolume <" << aTrack.GetVolume()->GetName()
583 << ">" << G4endl;
584 }
585}
586
587//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
588
590{
591 // based on G4Decay::GetMeanFreePath check; fixme
592
593 // get particle
594 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
595 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
596 G4double aMass = aParticle->GetMass();
597 G4double aLife = aParticleDef->GetPDGLifeTime();
598
599 // returns the mean free path in GEANT4 internal units
600 G4double pathlength;
601 G4double aCtau = c_light * aLife;
602
603 // check if the particle is stable?
604 if (aParticleDef->GetPDGStable()) {
605 pathlength = DBL_MAX;
606
607 //check if the particle has very short life time ?
608 } else if (aCtau < DBL_MIN) {
609 pathlength = DBL_MIN;
610
611 } else {
612 //calculate the mean free path
613 // by using normalized kinetic energy (= Ekin/mass)
614 G4double rKineticEnergy = aParticle->GetKineticEnergy()/aMass;
615 const G4double HighestValue = 20.0; //
616 if ( rKineticEnergy > HighestValue) {
617 // gamma >> 1
618 pathlength = ( rKineticEnergy + 1.0)* aCtau;
619 } else if ( rKineticEnergy < DBL_MIN ) {
620 // too slow particle
621#ifdef G4VERBOSE
622 if (GetVerboseLevel()>1) {
623 G4cout << "G4MuonicAtomDecay::GetMeanFreePath() !!particle stops!!";
624 G4cout << aParticleDef->GetParticleName() << G4endl;
625 G4cout << "KineticEnergy:" << aParticle->GetKineticEnergy()/GeV <<"[GeV]";
626 }
627#endif
628 pathlength = DBL_MIN;
629 } else {
630 // beta <1
631 pathlength = (aParticle->GetTotalMomentum())/aMass*aCtau ;
632 }
633 }
634 return pathlength;
635}
double A(double temperature)
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ NotForced
@ stopAndKill
@ fDecay
CLHEP::Hep3Vector G4ThreeVector
@ fAlive
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
HepLorentzVector & rotate(double, const Hep3Vector &)
void DumpInfo() const
G4int entries() const
G4DynamicParticle * PopProducts()
G4bool IsChecked() const
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
G4VDecayChannel * GetDecayChannel(G4int index) const
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
Definition: G4DecayTable.cc:82
G4int entries() const
G4double GetMass() const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
void Set4Momentum(const G4LorentzVector &momentum)
G4double GetTotalMomentum() const
G4double GetEnergyChange() const
const G4LorentzRotation & GetTrafoToLab() const
G4HadFinalStateStatus GetStatusChange() const
G4double GetLocalEnergyDeposit() const
const G4ThreeVector & GetMomentumChange() const
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
G4DynamicParticle * GetParticle()
G4double GetWeight() const
G4int GetCreatorModelType() const
void SetTime(G4double aT)
G4double GetTime() const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
const G4String & GetModelName() const
Definition: G4Ions.hh:52
const G4String & GetName() const
Definition: G4Material.hh:175
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
virtual G4bool IsApplicable(const G4ParticleDefinition &)
G4int GetVerboseLevel() const
G4MuonicAtomDecay(G4HadronicInteraction *hiptr=nullptr, const G4String &processName="MuonicAtomDecay")
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
virtual void ProcessDescription(std::ostream &outFile) const
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
virtual ~G4MuonicAtomDecay()
static G4double GetKShellEnergy(G4double A)
G4double GetNCLifeTime() const
G4Ions const * GetBaseIon() const
Definition: G4MuonicAtom.hh:96
G4double GetDIOLifeTime() const
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:212
void AddSecondary(G4Track *aSecondary)
void ProposeLocalTime(G4double t)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual void Initialize(const G4Track &)
G4ProcessManager * GetProcessManager() const
G4bool GetPDGStable() const
G4int GetAtomicNumber() const
const G4String & GetParticleType() const
G4int GetAtomicMass() const
G4DecayTable * GetDecayTable() const
G4double GetPDGLifeTime() const
const G4String & GetParticleName() const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
Definition: G4Step.hh:62
G4TrackStatus GetTrackStatus() const
G4int GetTrackID() const
const G4ParticleDefinition * GetParticleDefinition() const
G4VPhysicalVolume * GetVolume() const
G4double GetWeight() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4double GetLocalTime() const
G4Material * GetMaterial() const
G4ParticleDefinition * GetDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4int GetParentID() const
void SetCreatorModelIndex(G4int idx)
void SetGoodForTrackingFlag(G4bool value=true)
G4double GetBR() const
void SetVerboseLevel(G4int value)
G4int GetVerboseLevel() const
G4int GetNumberOfDaughters() const
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
virtual G4bool IsOKWithParentMass(G4double parentMass)
const G4String & GetDaughterName(G4int anIndex) const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
const G4String & GetName() const
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:80
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:424
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:359
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:361
G4int theProcessSubType
Definition: G4VProcess.hh:349
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
#define DBL_MIN
Definition: templates.hh:54
#define DBL_MAX
Definition: templates.hh:62
#define ns
Definition: xmlparse.cc:614