Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HadronicProcess.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 source file
30//
31// G4HadronicProcess
32//
33// original by H.P.Wellisch
34// J.L. Chuma, TRIUMF, 10-Mar-1997
35//
36// Modifications:
37// 05-Jul-2010 V.Ivanchenko cleanup commented lines
38// 20-Jul-2011 M.Kelsey -- null-pointer checks in DumpState()
39// 24-Sep-2011 M.Kelsey -- Use envvar G4HADRONIC_RANDOM_FILE to save random
40// engine state before each model call
41// 18-Oct-2011 M.Kelsey -- Handle final-state cases in conservation checks.
42// 14-Mar-2012 G.Folger -- enhance checks for conservation of energy, etc.
43// 28-Jul-2012 M.Maire -- add function GetTargetDefinition()
44// 14-Sep-2012 Inherit from RestDiscrete, use subtype code (now in ctor) to
45// configure base-class
46// 28-Sep-2012 Restore inheritance from G4VDiscreteProcess, remove enable-flag
47// changing, remove warning message from original ctor.
48// 21-Aug-2019 V.Ivanchenko leave try/catch only for ApplyYourself(..), cleanup
49
50#include "G4HadronicProcess.hh"
51
52#include "G4Types.hh"
53#include "G4SystemOfUnits.hh"
54#include "G4HadProjectile.hh"
55#include "G4ElementVector.hh"
56#include "G4Track.hh"
57#include "G4Step.hh"
58#include "G4Element.hh"
59#include "G4ParticleChange.hh"
60#include "G4ProcessVector.hh"
61#include "G4ProcessManager.hh"
62#include "G4NucleiProperties.hh"
63
67
68#include "G4NistManager.hh"
71#include "G4Exp.hh"
72
73#include <typeinfo>
74#include <sstream>
75#include <iostream>
76
77// File-scope variable to capture environment variable at startup
78
79static const char* G4Hadronic_Random_File = std::getenv("G4HADRONIC_RANDOM_FILE");
80
81//////////////////////////////////////////////////////////////////
82
84 G4ProcessType procType)
85 : G4VDiscreteProcess(processName, procType)
86{
87 SetProcessSubType(fHadronInelastic); // Default unless subclass changes
88 InitialiseLocal();
89}
90
91//////////////////////////////////////////////////////////////////
92
94 G4HadronicProcessType aHadSubType)
95 : G4VDiscreteProcess(processName, fHadronic)
96{
97 SetProcessSubType(aHadSubType);
98 InitialiseLocal();
99}
100
102{
103 theProcessStore->DeRegister(this);
104 delete theTotalResult;
105 delete theCrossSectionDataStore;
106}
107
108void G4HadronicProcess::InitialiseLocal() {
111 theInteraction = nullptr;
112 theCrossSectionDataStore = new G4CrossSectionDataStore();
113 theProcessStore = G4HadronicProcessStore::Instance();
114 theProcessStore->Register(this);
115 theInitialNumberOfInteractionLength = 0.0;
116 aScaleFactor = 1.0;
117 fWeight = 1.0;
118 nMatWarn = nKaonWarn = 0;
119 useIntegralXS = true;
120 theLastCrossSection = 0.0;
121 nICelectrons = 0;
122 idxIC = -1;
123 G4HadronicProcess_debug_flag = false;
124 levelsSetByProcess = false;
125 epReportLevel = 0;
126 epCheckLevels.first = DBL_MAX;
127 epCheckLevels.second = DBL_MAX;
128 GetEnergyMomentumCheckEnvvars();
129}
130
131void G4HadronicProcess::GetEnergyMomentumCheckEnvvars() {
132 if ( std::getenv("G4Hadronic_epReportLevel") ) {
133 epReportLevel = std::strtol(std::getenv("G4Hadronic_epReportLevel"),0,10);
134 }
135 if ( std::getenv("G4Hadronic_epCheckRelativeLevel") ) {
136 epCheckLevels.first = std::strtod(std::getenv("G4Hadronic_epCheckRelativeLevel"),0);
137 }
138 if ( std::getenv("G4Hadronic_epCheckAbsoluteLevel") ) {
139 epCheckLevels.second = std::strtod(std::getenv("G4Hadronic_epCheckAbsoluteLevel"),0);
140 }
141}
142
144{
145 if(!a) { return; }
146 theEnergyRangeManager.RegisterMe( a );
148}
149
152 const G4Element * elm,
153 const G4Material* mat)
154{
155 if(!mat)
156 {
157 static const G4int nmax = 5;
158 if(nMatWarn < nmax) {
159 ++nMatWarn;
161 ed << "Cannot compute Element x-section for " << GetProcessName()
162 << " because no material defined \n"
163 << " Please, specify material pointer or define simple material"
164 << " for Z= " << elm->GetZasInt();
165 G4Exception("G4HadronicProcess::GetElementCrossSection", "had066",
166 JustWarning, ed);
167 }
168 }
169 return
170 std::max(theCrossSectionDataStore->GetCrossSection(part, elm, mat),0.0);
171}
172
174{
175 if(std::getenv("G4HadronicProcess_debug")) {
176 G4HadronicProcess_debug_flag = true;
177 }
178 theProcessStore->RegisterParticle(this, &p);
179}
180
182{
183 theCrossSectionDataStore->BuildPhysicsTable(p);
184 theEnergyRangeManager.BuildPhysicsTable(p);
186}
187
190{
191 //G4cout << "GetMeanFreePath " << aTrack.GetDefinition()->GetParticleName()
192 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
193 theLastCrossSection = aScaleFactor*theCrossSectionDataStore
195 G4double res = (theLastCrossSection>0.0) ? 1.0/theLastCrossSection : DBL_MAX;
196 //G4cout << " xsection= " << theLastCrossSection << G4endl;
197 return res;
198}
199
202{
203 //G4cout << "PostStepDoIt " << aTrack.GetDefinition()->GetParticleName()
204 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
205 // if primary is not Alive then do nothing
207 theTotalResult->Initialize(aTrack);
208 fWeight = aTrack.GetWeight();
210 if(aTrack.GetTrackStatus() != fAlive) { return theTotalResult; }
211
212 // Find cross section at end of step and check if <= 0
213 //
214 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
215 const G4Material* aMaterial = aTrack.GetMaterial();
216
217 // check only for charged particles
218 if(aParticle->GetDefinition()->GetPDGCharge() != 0.0) {
219 G4double xs = aScaleFactor*
220 theCrossSectionDataStore->ComputeCrossSection(aParticle,aMaterial);
221 if(xs <= 0.0 || xs < theLastCrossSection*G4UniformRand()) {
222 // No interaction
223 return theTotalResult;
224 }
225 }
226
227 const G4Element* anElement =
228 theCrossSectionDataStore->SampleZandA(aParticle,aMaterial,targetNucleus);
229
230 // Next check for illegal track status
231 //
232 if (aTrack.GetTrackStatus() != fAlive &&
233 aTrack.GetTrackStatus() != fSuspend) {
234 if (aTrack.GetTrackStatus() == fStopAndKill ||
238 ed << "G4HadronicProcess: track in unusable state - "
239 << aTrack.GetTrackStatus() << G4endl;
240 ed << "G4HadronicProcess: returning unchanged track " << G4endl;
241 DumpState(aTrack,"PostStepDoIt",ed);
242 G4Exception("G4HadronicProcess::PostStepDoIt", "had004", JustWarning, ed);
243 }
244 // No warning for fStopButAlive which is a legal status here
245 return theTotalResult;
246 }
247
248 // Initialize the hadronic projectile from the track
249 thePro.Initialise(aTrack);
250
251 theInteraction = ChooseHadronicInteraction(thePro, targetNucleus,
252 aMaterial, anElement);
253 if(!theInteraction) {
255 ed << "Target element "<<anElement->GetName()<<" Z= "
256 << targetNucleus.GetZ_asInt() << " A= "
257 << targetNucleus.GetA_asInt() << G4endl;
258 DumpState(aTrack,"ChooseHadronicInteraction",ed);
259 ed << " No HadronicInteraction found out" << G4endl;
260 G4Exception("G4HadronicProcess::PostStepDoIt", "had005", FatalException, ed);
261 return theTotalResult;
262 }
263
264 G4HadFinalState* result = nullptr;
265 G4int reentryCount = 0;
266 /*
267 G4cout << "### " << aParticle->GetDefinition()->GetParticleName()
268 << " Ekin(MeV)= " << aParticle->GetKineticEnergy()
269 << " Z= " << targetNucleus.GetZ_asInt()
270 << " A= " << targetNucleus.GetA_asInt()
271 << " by " << theInteraction->GetModelName()
272 << G4endl;
273 */
274 do
275 {
276 try
277 {
278 // Save random engine if requested for debugging
279 if (G4Hadronic_Random_File) {
280 CLHEP::HepRandom::saveEngineStatus(G4Hadronic_Random_File);
281 }
282 // Call the interaction
283 result = theInteraction->ApplyYourself( thePro, targetNucleus);
284 ++reentryCount;
285 }
286 catch(G4HadronicException & aR)
287 {
289 aR.Report(ed);
290 ed << "Call for " << theInteraction->GetModelName() << G4endl;
291 ed << "Target element "<<anElement->GetName()<<" Z= "
292 << targetNucleus.GetZ_asInt()
293 << " A= " << targetNucleus.GetA_asInt() << G4endl;
294 DumpState(aTrack,"ApplyYourself",ed);
295 ed << " ApplyYourself failed" << G4endl;
296 G4Exception("G4HadronicProcess::PostStepDoIt", "had006", FatalException,
297 ed);
298 }
299
300 // Check the result for catastrophic energy non-conservation
301 result = CheckResult(thePro, targetNucleus, result);
302
303 if(reentryCount>100) {
305 ed << "Call for " << theInteraction->GetModelName() << G4endl;
306 ed << "Target element "<<anElement->GetName()<<" Z= "
307 << targetNucleus.GetZ_asInt()
308 << " A= " << targetNucleus.GetA_asInt() << G4endl;
309 DumpState(aTrack,"ApplyYourself",ed);
310 ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
311 G4Exception("G4HadronicProcess::PostStepDoIt", "had006", FatalException,
312 ed);
313 }
314 }
315 while(!result); /* Loop checking, 30-Oct-2015, G.Folger */
316
317 // Check whether kaon0 or anti_kaon0 are present between the secondaries:
318 // if this is the case, transform them into either kaon0S or kaon0L,
319 // with equal, 50% probability, keeping their dynamical masses (and
320 // the other kinematical properties).
321 // When this happens - very rarely - a "JustWarning" exception is thrown.
322 G4int nSec = result->GetNumberOfSecondaries();
323 if ( nSec > 0 ) {
324 for ( G4int i = 0; i < nSec; ++i ) {
325 G4DynamicParticle* dynamicParticle = result->GetSecondary(i)->GetParticle();
326 const G4ParticleDefinition* particleDefinition =
327 dynamicParticle->GetParticleDefinition();
328 if ( particleDefinition == G4KaonZero::Definition() ||
329 particleDefinition == G4AntiKaonZero::Definition() ) {
330 G4ParticleDefinition* newPart;
331 if( G4UniformRand() > 0.5 ) { newPart = G4KaonZeroShort::Definition(); }
332 else { newPart = G4KaonZeroLong::Definition(); }
333 dynamicParticle->SetDefinition( newPart );
334 if(nKaonWarn < 5) {
335 ++nKaonWarn;
337 ed << " Hadronic model " << theInteraction->GetModelName() << G4endl;
338 ed << " created " << particleDefinition->GetParticleName() << G4endl;
339 ed << " -> forced to be " << newPart->GetParticleName() << G4endl;
340 G4Exception( "G4HadronicProcess::PostStepDoIt", "had007", JustWarning, ed );
341 }
342 }
343 }
344 }
345
347
349
350 FillResult(result, aTrack);
351
352 if (epReportLevel != 0) {
353 CheckEnergyMomentumConservation(aTrack, targetNucleus);
354 }
355 //G4cout << "PostStepDoIt done nICelectrons= " << nICelectrons << G4endl;
356 return theTotalResult;
357}
358
359void G4HadronicProcess::ProcessDescription(std::ostream& outFile) const
360{
361 outFile << "The description for this process has not been written yet.\n";
362}
363
364
365G4double G4HadronicProcess::XBiasSurvivalProbability()
366{
368 G4double biasedProbability = 1.-G4Exp(-nLTraversed);
369 G4double realProbability = 1-G4Exp(-nLTraversed/aScaleFactor);
370 G4double result = (biasedProbability-realProbability)/biasedProbability;
371 return result;
372}
373
374G4double G4HadronicProcess::XBiasSecondaryWeight()
375{
377 G4double result =
378 1./aScaleFactor*G4Exp(-nLTraversed/aScaleFactor*(1-1./aScaleFactor));
379 return result;
380}
381
382void
384{
386 const G4ThreeVector& dir = aT.GetMomentumDirection();
387
388 G4double efinal = std::max(aR->GetEnergyChange(), 0.0);
389
390 // check status of primary
391 if(aR->GetStatusChange() == stopAndKill) {
394
395 // check its final energy
396 } else if(0.0 == efinal) {
399 ->GetAtRestProcessVector()->size() > 0)
402
403 // primary is not killed apply rotation and Lorentz transformation
404 } else {
406 G4ThreeVector newDir = aR->GetMomentumChange();
407 newDir.rotateUz(dir);
410 }
411 //G4cout << "FillResult: Efinal= " << efinal << " status= "
412 // << theTotalResult->GetTrackStatus()
413 // << " fKill= " << fStopAndKill << G4endl;
414
415 // check secondaries
416 nICelectrons = 0;
417 if(idxIC == -1) {
418 G4int idx = G4PhysicsModelCatalog::GetIndex("e-InternalConvertion");
419 idxIC = -1 == idx ? -2 : idx;
420 }
421 G4int nSec = aR->GetNumberOfSecondaries();
423 G4double time0 = aT.GetGlobalTime();
424
425 for (G4int i = 0; i < nSec; ++i) {
426 G4DynamicParticle* dynParticle = aR->GetSecondary(i)->GetParticle();
427
428 // apply rotation
429 G4ThreeVector newDir = dynParticle->GetMomentumDirection();
430 newDir.rotateUz(dir);
431 dynParticle->SetMomentumDirection(newDir);
432
433 // check if secondary is on the mass shell
434 const G4ParticleDefinition* part = dynParticle->GetDefinition();
435 G4double mass = part->GetPDGMass();
436 G4double dmass= dynParticle->GetMass();
437 const G4double delta_mass_lim = 1.0*CLHEP::keV;
438 const G4double delta_ekin = 0.001*CLHEP::eV;
439 if(std::abs(dmass - mass) > delta_mass_lim) {
440 G4double e = std::max(dynParticle->GetKineticEnergy() + dmass - mass, delta_ekin);
441 if(G4HadronicProcess_debug_flag) {
443 ed << "TrackID= "<< aT.GetTrackID()
444 << " " << aT.GetParticleDefinition()->GetParticleName()
445 << " Target Z= " << targetNucleus.GetZ_asInt() << " A= "
446 << targetNucleus.GetA_asInt()
447 << " Ekin(GeV)= " << aT.GetKineticEnergy()/CLHEP::GeV
448 << "\n Secondary is out of mass shell: " << part->GetParticleName()
449 << " EkinNew(MeV)= " << e
450 << " DeltaMass(MeV)= " << dmass - mass << G4endl;
451 G4Exception("G4HadronicProcess::FillResults", "had012", JustWarning, ed);
452 }
453 dynParticle->SetKineticEnergy(e);
454 dynParticle->SetMass(mass);
455 }
456 G4int idxModel = aR->GetSecondary(i)->GetCreatorModelType();
457 //if(idxIC == idxModel) { ++nICelectrons; }
458 if(part->GetPDGEncoding() == 11) { ++nICelectrons; }
459
460 // time of interaction starts from zero + global time
461 G4double time = std::max(aR->GetSecondary(i)->GetTime(), 0.0) + time0;
462
463 G4Track* track = new G4Track(dynParticle, time, aT.GetPosition());
464 track->SetCreatorModelIndex(idxModel);
465 G4double newWeight = fWeight*aR->GetSecondary(i)->GetWeight();
466 track->SetWeight(newWeight);
469 if (G4HadronicProcess_debug_flag) {
470 G4double e = dynParticle->GetKineticEnergy();
471 if (e == 0.0) {
473 DumpState(aT,"Secondary has zero energy",ed);
474 ed << "Secondary " << part->GetParticleName()
475 << G4endl;
476 G4Exception("G4HadronicProcess::FillResults", "had011",
477 JustWarning,ed);
478 }
479 }
480 }
481 aR->Clear();
482 // G4cout << "FillResults done nICe= " << nICelectrons << G4endl;
483}
484
486{
488}
489
491{
492 if (aScale <= 0.0) {
494 ed << " Wrong biasing factor " << aScale << " for " << GetProcessName();
495 G4Exception("G4HadronicProcess::BiasCrossSectionByFactor", "had010",
496 JustWarning, ed, "Cross-section bias is ignored");
497 } else {
498 aScaleFactor = aScale;
499 }
500}
501
503 const G4Nucleus &aNucleus,
504 G4HadFinalState * result)
505{
506 // check for catastrophic energy non-conservation
507 // to re-sample the interaction
508
510 G4double nuclearMass(0);
511 if (theModel) {
512
513 // Compute final-state total energy
514 G4double finalE(0.);
515 G4int nSec = result->GetNumberOfSecondaries();
516
517 nuclearMass = G4NucleiProperties::GetNuclearMass(aNucleus.GetA_asInt(),
518 aNucleus.GetZ_asInt());
519 if (result->GetStatusChange() != stopAndKill) {
520 // Interaction didn't complete, returned "do nothing" state
521 // and reset nucleus or the primary survived the interaction
522 // (e.g. electro-nuclear ) => keep nucleus
523 finalE=result->GetLocalEnergyDeposit() +
524 aPro.GetDefinition()->GetPDGMass() + result->GetEnergyChange();
525 if( nSec == 0 ){
526 // Since there are no secondaries, there is no recoil nucleus.
527 // To check energy balance we must neglect the initial nucleus too.
528 nuclearMass=0.0;
529 }
530 }
531 for (G4int i = 0; i < nSec; i++) {
532 G4DynamicParticle *pdyn=result->GetSecondary(i)->GetParticle();
533 finalE += pdyn->GetTotalEnergy();
534 G4double mass_pdg=pdyn->GetDefinition()->GetPDGMass();
535 G4double mass_dyn=pdyn->GetMass();
536 if ( std::abs(mass_pdg - mass_dyn) > 0.1*mass_pdg + 1.*MeV ) {
537 // If it is shortlived, then a difference less than 3 times the width is acceptable
538 if ( pdyn->GetDefinition()->IsShortLived() &&
539 std::abs(mass_pdg - mass_dyn) < 3.0*pdyn->GetDefinition()->GetPDGWidth() ) {
540 continue;
541 }
542 result->Clear();
543 result = nullptr;
545 desc << "Warning: Secondary with off-shell dynamic mass detected: "
546 << G4endl
547 << " " << pdyn->GetDefinition()->GetParticleName()
548 << ", PDG mass: " << mass_pdg << ", dynamic mass: "
549 << mass_dyn << G4endl
550 << (epReportLevel<0 ? "abort the event"
551 : "re-sample the interaction") << G4endl
552 << " Process / Model: " << GetProcessName()<< " / "
553 << theModel->GetModelName() << G4endl
554 << " Primary: " << aPro.GetDefinition()->GetParticleName()
555 << " (" << aPro.GetDefinition()->GetPDGEncoding() << "), "
556 << " E= " << aPro.Get4Momentum().e()
557 << ", target nucleus (" << aNucleus.GetZ_asInt() << ", "
558 << aNucleus.GetA_asInt() << ")" << G4endl;
559 G4Exception("G4HadronicProcess:CheckResult()", "had012",
561 // must return here.....
562 return result;
563 }
564 }
565 G4double deltaE= nuclearMass + aPro.GetTotalEnergy() - finalE;
566
567 std::pair<G4double, G4double> checkLevels =
568 theModel->GetFatalEnergyCheckLevels(); // (relative, absolute)
569 if (std::abs(deltaE) > checkLevels.second &&
570 std::abs(deltaE) > checkLevels.first*aPro.GetKineticEnergy()){
571 // do not delete result, this is a pointer to a data member;
572 result->Clear();
573 result = nullptr;
575 desc << "Warning: Bad energy non-conservation detected, will "
576 << (epReportLevel<0 ? "abort the event"
577 : "re-sample the interaction") << G4endl
578 << " Process / Model: " << GetProcessName()<< " / "
579 << theModel->GetModelName() << G4endl
580 << " Primary: " << aPro.GetDefinition()->GetParticleName()
581 << " (" << aPro.GetDefinition()->GetPDGEncoding() << "), "
582 << " E= " << aPro.Get4Momentum().e()
583 << ", target nucleus (" << aNucleus.GetZ_asInt() << ", "
584 << aNucleus.GetA_asInt() << ")" << G4endl
585 << " E(initial - final) = " << deltaE << " MeV." << G4endl;
586 G4Exception("G4HadronicProcess:CheckResult()", "had012",
588 }
589 }
590 return result;
591}
592
593void
595 const G4Nucleus& aNucleus)
596{
597 G4int target_A=aNucleus.GetA_asInt();
598 G4int target_Z=aNucleus.GetZ_asInt();
599 G4double targetMass = G4NucleiProperties::GetNuclearMass(target_A,target_Z);
600 G4LorentzVector target4mom(0, 0, 0, targetMass
601 + nICelectrons*CLHEP::electron_mass_c2);
602
603 G4LorentzVector projectile4mom = aTrack.GetDynamicParticle()->Get4Momentum();
604 G4int track_A = aTrack.GetDefinition()->GetBaryonNumber();
605 G4int track_Z = G4lrint(aTrack.GetDefinition()->GetPDGCharge());
606
607 G4int initial_A = target_A + track_A;
608 G4int initial_Z = target_Z + track_Z - nICelectrons;
609
610 G4LorentzVector initial4mom = projectile4mom + target4mom;
611
612 // Compute final-state momentum for scattering and "do nothing" results
613 G4LorentzVector final4mom;
614 G4int final_A(0), final_Z(0);
615
617 if (theTotalResult->GetTrackStatus() != fStopAndKill) { // If it is Alive
618 // Either interaction didn't complete, returned "do nothing" state
619 // or the primary survived the interaction (e.g. electro-nucleus )
620
621 // Interaction didn't complete, returned "do nothing" state
622 // - or suppressed recoil (e.g. Neutron elastic )
623 final4mom = initial4mom;
624 final_A = initial_A;
625 final_Z = initial_Z;
626 if (nSec > 0) {
627 // The primary remains in final state (e.g. electro-nucleus )
628 // Use the final energy / momentum
631 G4double mass = aTrack.GetDefinition()->GetPDGMass();
632 G4double ptot = std::sqrt(ekin*(ekin + 2*mass));
633 final4mom.set(ptot*v.x(), ptot*v.y(), ptot*v.z(), mass + ekin);
634 final_A = track_A;
635 final_Z = track_Z;
636 // Expect that the target nucleus will have interacted,
637 // and its products, including recoil, will be included in secondaries.
638 }
639 }
640 if( nSec > 0 ) {
641 G4Track* sec;
642
643 for (G4int i = 0; i < nSec; i++) {
645 final4mom += sec->GetDynamicParticle()->Get4Momentum();
646 final_A += sec->GetDefinition()->GetBaryonNumber();
647 final_Z += G4lrint(sec->GetDefinition()->GetPDGCharge());
648 }
649 }
650
651 // Get level-checking information (used to cut-off relative checks)
652 G4String processName = GetProcessName();
654 G4String modelName("none");
655 if (theModel) modelName = theModel->GetModelName();
656 std::pair<G4double, G4double> checkLevels = epCheckLevels;
657 if (!levelsSetByProcess) {
658 if (theModel) checkLevels = theModel->GetEnergyMomentumCheckLevels();
659 checkLevels.first= std::min(checkLevels.first, epCheckLevels.first);
660 checkLevels.second=std::min(checkLevels.second, epCheckLevels.second);
661 }
662
663 // Compute absolute total-energy difference, and relative kinetic-energy
664 G4bool checkRelative = (aTrack.GetKineticEnergy() > checkLevels.second);
665
666 G4LorentzVector diff = initial4mom - final4mom;
667 G4double absolute = diff.e();
668 G4double relative = checkRelative ? absolute/aTrack.GetKineticEnergy() : 0.;
669
670 G4double absolute_mom = diff.vect().mag();
671 G4double relative_mom = checkRelative ? absolute_mom/aTrack.GetMomentum().mag() : 0.;
672
673 // Evaluate relative and absolute conservation
674 G4bool relPass = true;
675 G4String relResult = "pass";
676 if ( std::abs(relative) > checkLevels.first
677 || std::abs(relative_mom) > checkLevels.first) {
678 relPass = false;
679 relResult = checkRelative ? "fail" : "N/A";
680 }
681
682 G4bool absPass = true;
683 G4String absResult = "pass";
684 if ( std::abs(absolute) > checkLevels.second
685 || std::abs(absolute_mom) > checkLevels.second ) {
686 absPass = false ;
687 absResult = "fail";
688 }
689
690 G4bool chargePass = true;
691 G4String chargeResult = "pass";
692 if ( (initial_A-final_A)!=0
693 || (initial_Z-final_Z)!=0 ) {
694 chargePass = checkLevels.second < DBL_MAX ? false : true;
695 chargeResult = "fail";
696 }
697
698 G4bool conservationPass = (relPass || absPass) && chargePass;
699
700 std::stringstream Myout;
701 G4bool Myout_notempty(false);
702 // Options for level of reporting detail:
703 // 0. off
704 // 1. report only when E/p not conserved
705 // 2. report regardless of E/p conservation
706 // 3. report only when E/p not conserved, with model names, process names, and limits
707 // 4. report regardless of E/p conservation, with model names, process names, and limits
708 // negative -1.., as above, but send output to stderr
709
710 if( std::abs(epReportLevel) == 4
711 || ( std::abs(epReportLevel) == 3 && ! conservationPass ) ){
712 Myout << " Process: " << processName << " , Model: " << modelName << G4endl;
713 Myout << " Primary: " << aTrack.GetParticleDefinition()->GetParticleName()
714 << " (" << aTrack.GetParticleDefinition()->GetPDGEncoding() << "),"
715 << " E= " << aTrack.GetDynamicParticle()->Get4Momentum().e()
716 << ", target nucleus (" << aNucleus.GetZ_asInt() << ","
717 << aNucleus.GetA_asInt() << ")" << G4endl;
718 Myout_notempty=true;
719 }
720 if ( std::abs(epReportLevel) == 4
721 || std::abs(epReportLevel) == 2
722 || ! conservationPass ){
723
724 Myout << " "<< relResult <<" relative, limit " << checkLevels.first << ", values E/T(0) = "
725 << relative << " p/p(0)= " << relative_mom << G4endl;
726 Myout << " "<< absResult << " absolute, limit (MeV) " << checkLevels.second/MeV << ", values E / p (MeV) = "
727 << absolute/MeV << " / " << absolute_mom/MeV << " 3mom: " << (diff.vect())*1./MeV << G4endl;
728 Myout << " "<< chargeResult << " charge/baryon number balance " << (initial_Z-final_Z) << " / " << (initial_A-final_A) << " "<< G4endl;
729 Myout_notempty=true;
730
731 }
732 Myout.flush();
733 if ( Myout_notempty ) {
734 if (epReportLevel > 0) G4cout << Myout.str()<< G4endl;
735 else if (epReportLevel < 0) G4cerr << Myout.str()<< G4endl;
736 }
737}
738
740 const G4String& method,
742{
743 ed << "Unrecoverable error in the method " << method << " of "
744 << GetProcessName() << G4endl;
745 ed << "TrackID= "<< aTrack.GetTrackID() << " ParentID= "
746 << aTrack.GetParentID()
747 << " " << aTrack.GetParticleDefinition()->GetParticleName()
748 << G4endl;
749 ed << "Ekin(GeV)= " << aTrack.GetKineticEnergy()/CLHEP::GeV
750 << "; direction= " << aTrack.GetMomentumDirection() << G4endl;
751 ed << "Position(mm)= " << aTrack.GetPosition()/CLHEP::mm << ";";
752
753 if (aTrack.GetMaterial()) {
754 ed << " material " << aTrack.GetMaterial()->GetName();
755 }
756 ed << G4endl;
757
758 if (aTrack.GetVolume()) {
759 ed << "PhysicalVolume <" << aTrack.GetVolume()->GetName()
760 << ">" << G4endl;
761 }
762}
763
765{
766 theCrossSectionDataStore->DumpPhysicsTable(p);
767}
768
770{
771 theCrossSectionDataStore->AddDataSet(aDataSet);
772}
773
774std::vector<G4HadronicInteraction*>&
776{
777 return theEnergyRangeManager.GetHadronicInteractionList();
778}
779
782{
783 std::vector<G4HadronicInteraction*>& list
784 = theEnergyRangeManager.GetHadronicInteractionList();
785 for (size_t li=0; li<list.size(); li++) {
786 if (list[li]->GetModelName() == modelName) return list[li];
787 }
788 return nullptr;
789}
@ JustWarning
@ FatalException
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4ForceCondition
@ stopAndKill
G4HadronicProcessType
@ fHadronInelastic
G4ProcessType
@ fHadronic
@ fKillTrackAndSecondaries
@ fSuspend
@ fAlive
@ fStopAndKill
@ fStopButAlive
@ fPostponeToNextEvent
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
double x() const
double y() const
double mag() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
Hep3Vector vect() const
void set(double x, double y, double z, double t)
static void saveEngineStatus(const char filename[]="Config.conf")
Definition: Random.cc:278
static G4AntiKaonZero * Definition()
void BuildPhysicsTable(const G4ParticleDefinition &)
void AddDataSet(G4VCrossSectionDataSet *)
void DumpPhysicsTable(const G4ParticleDefinition &)
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
const G4Element * SampleZandA(const G4DynamicParticle *, const G4Material *, G4Nucleus &target)
G4double GetMass() const
void SetMomentumDirection(const G4ThreeVector &aDirection)
const G4ThreeVector & GetMomentumDirection() const
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMass(G4double mass)
const G4ParticleDefinition * GetParticleDefinition() const
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
void SetKineticEnergy(G4double aEnergy)
const G4String & GetName() const
Definition: G4Element.hh:126
G4int GetZasInt() const
Definition: G4Element.hh:131
void RegisterMe(G4HadronicInteraction *a)
void BuildPhysicsTable(const G4ParticleDefinition &)
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList()
G4double GetEnergyChange() const
G4HadFinalStateStatus GetStatusChange() const
void SetTrafoToLab(const G4LorentzRotation &aT)
G4double GetLocalEnergyDeposit() const
const G4ThreeVector & GetMomentumChange() const
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
void Initialise(const G4Track &aT)
const G4ParticleDefinition * GetDefinition() const
G4LorentzRotation & GetTrafoToLab()
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const
G4DynamicParticle * GetParticle()
G4double GetWeight() const
G4int GetCreatorModelType() const
G4double GetTime() const
void Report(std::ostream &aS) const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
virtual std::pair< G4double, G4double > GetEnergyMomentumCheckLevels() const
const G4String & GetModelName() const
void DeRegister(G4HadronicProcess *)
void RegisterParticle(G4HadronicProcess *, const G4ParticleDefinition *)
static G4HadronicProcessStore * Instance()
void RegisterInteraction(G4HadronicProcess *, G4HadronicInteraction *)
void Register(G4HadronicProcess *)
void PrintInfo(const G4ParticleDefinition *)
void FillResult(G4HadFinalState *aR, const G4Track &aT)
G4HadProjectile thePro
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
void ProcessDescription(std::ostream &outFile) const override
void BiasCrossSectionByFactor(G4double aScale)
G4HadFinalState * CheckResult(const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
G4HadronicInteraction * GetHadronicInteraction() const
G4ParticleChange * theTotalResult
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
G4double GetElementCrossSection(const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList()
void PreparePhysicsTable(const G4ParticleDefinition &) override
G4HadronicProcess(const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
G4HadronicInteraction * ChooseHadronicInteraction(const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement)
~G4HadronicProcess() override
void BuildPhysicsTable(const G4ParticleDefinition &) override
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
G4HadronicInteraction * GetHadronicModel(const G4String &)
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
void DumpPhysicsTable(const G4ParticleDefinition &p)
void MultiplyCrossSectionBy(G4double factor)
void RegisterMe(G4HadronicInteraction *a)
static G4KaonZeroLong * Definition()
static G4KaonZeroShort * Definition()
static G4KaonZero * Definition()
Definition: G4KaonZero.cc:52
const G4String & GetName() const
Definition: G4Material.hh:175
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
void AddSecondary(G4Track *aSecondary)
G4double GetEnergy() const
const G4ThreeVector * GetMomentumDirection() const
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual void Initialize(const G4Track &)
G4ProcessManager * GetProcessManager() const
G4double GetPDGWidth() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4int GetIndex(const G4String &)
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
Definition: G4Step.hh:62
std::size_t first(char) const
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
G4ThreeVector GetMomentum() 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 ProposeTrackStatus(G4TrackStatus status)
void SetSecondaryWeightByProcess(G4bool)
void ProposeWeight(G4double finalWeight)
G4int GetNumberOfSecondaries() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4Track * GetSecondary(G4int anIndex) const
G4TrackStatus GetTrackStatus() const
const G4String & GetName() const
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:424
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4double GetTotalNumberOfInteractionLengthTraversed() const
Definition: G4VProcess.hh:437
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62