Geant4-11
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"
70#include "G4Exp.hh"
71
72#include <typeinfo>
73#include <sstream>
74#include <iostream>
75
76// File-scope variable to capture environment variable at startup
77
78static const char* G4Hadronic_Random_File = std::getenv("G4HADRONIC_RANDOM_FILE");
79
81
83 G4ProcessType procType)
84 : G4VDiscreteProcess(processName, procType)
85{
86 SetProcessSubType(fHadronInelastic); // Default unless subclass changes
88}
89
91
93 G4HadronicProcessType aHadSubType)
94 : G4VDiscreteProcess(processName, fHadronic)
95{
96 SetProcessSubType(aHadSubType);
98}
99
101{
103 delete theTotalResult;
105}
106
110 theInteraction = nullptr;
115 aScaleFactor = 1.0;
116 fWeight = 1.0;
117 nMatWarn = nKaonWarn = 0;
118 useIntegralXS = true;
120 nICelectrons = 0;
122 levelsSetByProcess = false;
123 epReportLevel = 0;
124 epCheckLevels.first = DBL_MAX;
125 epCheckLevels.second = DBL_MAX;
127}
128
130 if ( std::getenv("G4Hadronic_epReportLevel") ) {
131 epReportLevel = std::strtol(std::getenv("G4Hadronic_epReportLevel"),0,10);
132 }
133 if ( std::getenv("G4Hadronic_epCheckRelativeLevel") ) {
134 epCheckLevels.first = std::strtod(std::getenv("G4Hadronic_epCheckRelativeLevel"),0);
135 }
136 if ( std::getenv("G4Hadronic_epCheckAbsoluteLevel") ) {
137 epCheckLevels.second = std::strtod(std::getenv("G4Hadronic_epCheckAbsoluteLevel"),0);
138 }
139}
140
142{
143 if(!a) { return; }
146}
147
150 const G4Element * elm,
151 const G4Material* mat)
152{
153 if(!mat)
154 {
155 static const G4int nmax = 5;
156 if(nMatWarn < nmax) {
157 ++nMatWarn;
159 ed << "Cannot compute Element x-section for " << GetProcessName()
160 << " because no material defined \n"
161 << " Please, specify material pointer or define simple material"
162 << " for Z= " << elm->GetZasInt();
163 G4Exception("G4HadronicProcess::GetElementCrossSection", "had066",
164 JustWarning, ed);
165 }
166 }
167 return
169}
170
172{
173 if(std::getenv("G4HadronicProcess_debug")) {
175 }
177}
178
180{
184}
185
188{
189 //G4cout << "GetMeanFreePath " << aTrack.GetDefinition()->GetParticleName()
190 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
194 //G4cout << " xsection= " << theLastCrossSection << G4endl;
195 return res;
196}
197
200{
201 //G4cout << "PostStepDoIt " << aTrack.GetDefinition()->GetParticleName()
202 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
203 // if primary is not Alive then do nothing
205 theTotalResult->Initialize(aTrack);
206 fWeight = aTrack.GetWeight();
208 if(aTrack.GetTrackStatus() != fAlive) { return theTotalResult; }
209
210 // Find cross section at end of step and check if <= 0
211 //
212 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
213 const G4Material* aMaterial = aTrack.GetMaterial();
214
215 // check only for charged particles
216 if(aParticle->GetDefinition()->GetPDGCharge() != 0.0) {
218 theCrossSectionDataStore->ComputeCrossSection(aParticle,aMaterial);
219 if(xs <= 0.0 || xs < theLastCrossSection*G4UniformRand()) {
220 // No interaction
221 return theTotalResult;
222 }
223 }
224
225 const G4Element* anElement =
227
228 // Next check for illegal track status
229 //
230 if (aTrack.GetTrackStatus() != fAlive &&
231 aTrack.GetTrackStatus() != fSuspend) {
232 if (aTrack.GetTrackStatus() == fStopAndKill ||
236 ed << "G4HadronicProcess: track in unusable state - "
237 << aTrack.GetTrackStatus() << G4endl;
238 ed << "G4HadronicProcess: returning unchanged track " << G4endl;
239 DumpState(aTrack,"PostStepDoIt",ed);
240 G4Exception("G4HadronicProcess::PostStepDoIt", "had004", JustWarning, ed);
241 }
242 // No warning for fStopButAlive which is a legal status here
243 return theTotalResult;
244 }
245
246 // Initialize the hadronic projectile from the track
247 thePro.Initialise(aTrack);
248
250 aMaterial, anElement);
251 if(!theInteraction) {
253 ed << "Target element "<<anElement->GetName()<<" Z= "
254 << targetNucleus.GetZ_asInt() << " A= "
256 DumpState(aTrack,"ChooseHadronicInteraction",ed);
257 ed << " No HadronicInteraction found out" << G4endl;
258 G4Exception("G4HadronicProcess::PostStepDoIt", "had005", FatalException, ed);
259 return theTotalResult;
260 }
261
262 G4HadFinalState* result = nullptr;
263 G4int reentryCount = 0;
264 /*
265 G4cout << "### " << aParticle->GetDefinition()->GetParticleName()
266 << " Ekin(MeV)= " << aParticle->GetKineticEnergy()
267 << " Z= " << targetNucleus.GetZ_asInt()
268 << " A= " << targetNucleus.GetA_asInt()
269 << " by " << theInteraction->GetModelName()
270 << G4endl;
271 */
272 do
273 {
274 try
275 {
276 // Save random engine if requested for debugging
279 }
280 // Call the interaction
282 ++reentryCount;
283 }
284 catch(G4HadronicException & aR)
285 {
287 aR.Report(ed);
288 ed << "Call for " << theInteraction->GetModelName() << G4endl;
289 ed << "Target element "<<anElement->GetName()<<" Z= "
291 << " A= " << targetNucleus.GetA_asInt() << G4endl;
292 DumpState(aTrack,"ApplyYourself",ed);
293 ed << " ApplyYourself failed" << G4endl;
294 G4Exception("G4HadronicProcess::PostStepDoIt", "had006", FatalException,
295 ed);
296 }
297
298 // Check the result for catastrophic energy non-conservation
299 result = CheckResult(thePro, targetNucleus, result);
300
301 if(reentryCount>100) {
303 ed << "Call for " << theInteraction->GetModelName() << G4endl;
304 ed << "Target element "<<anElement->GetName()<<" Z= "
306 << " A= " << targetNucleus.GetA_asInt() << G4endl;
307 DumpState(aTrack,"ApplyYourself",ed);
308 ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
309 G4Exception("G4HadronicProcess::PostStepDoIt", "had006", FatalException,
310 ed);
311 }
312 }
313 while(!result); /* Loop checking, 30-Oct-2015, G.Folger */
314
315 // Check whether kaon0 or anti_kaon0 are present between the secondaries:
316 // if this is the case, transform them into either kaon0S or kaon0L,
317 // with equal, 50% probability, keeping their dynamical masses (and
318 // the other kinematical properties).
319 // When this happens - very rarely - a "JustWarning" exception is thrown.
320 G4int nSec = result->GetNumberOfSecondaries();
321 if ( nSec > 0 ) {
322 for ( G4int i = 0; i < nSec; ++i ) {
323 G4DynamicParticle* dynamicParticle = result->GetSecondary(i)->GetParticle();
324 const G4ParticleDefinition* particleDefinition =
325 dynamicParticle->GetParticleDefinition();
326 if ( particleDefinition == G4KaonZero::Definition() ||
327 particleDefinition == G4AntiKaonZero::Definition() ) {
328 G4ParticleDefinition* newPart;
329 if( G4UniformRand() > 0.5 ) { newPart = G4KaonZeroShort::Definition(); }
330 else { newPart = G4KaonZeroLong::Definition(); }
331 dynamicParticle->SetDefinition( newPart );
332 if(nKaonWarn < 5) {
333 ++nKaonWarn;
335 ed << " Hadronic model " << theInteraction->GetModelName() << G4endl;
336 ed << " created " << particleDefinition->GetParticleName() << G4endl;
337 ed << " -> forced to be " << newPart->GetParticleName() << G4endl;
338 G4Exception( "G4HadronicProcess::PostStepDoIt", "had007", JustWarning, ed );
339 }
340 }
341 }
342 }
343
345
347
348 FillResult(result, aTrack);
349
350 if (epReportLevel != 0) {
352 }
353 //G4cout << "PostStepDoIt done nICelectrons= " << nICelectrons << G4endl;
354 return theTotalResult;
355}
356
357void G4HadronicProcess::ProcessDescription(std::ostream& outFile) const
358{
359 outFile << "The description for this process has not been written yet.\n";
360}
361
362
364{
366 G4double biasedProbability = 1.-G4Exp(-nLTraversed);
367 G4double realProbability = 1-G4Exp(-nLTraversed/aScaleFactor);
368 G4double result = (biasedProbability-realProbability)/biasedProbability;
369 return result;
370}
371
373{
375 G4double result =
376 1./aScaleFactor*G4Exp(-nLTraversed/aScaleFactor*(1-1./aScaleFactor));
377 return result;
378}
379
380void
382{
384 const G4ThreeVector& dir = aT.GetMomentumDirection();
385
386 G4double efinal = std::max(aR->GetEnergyChange(), 0.0);
387
388 // check status of primary
389 if(aR->GetStatusChange() == stopAndKill) {
392
393 // check its final energy
394 } else if(0.0 == efinal) {
397 ->GetAtRestProcessVector()->size() > 0)
400
401 // primary is not killed apply rotation and Lorentz transformation
402 } else {
404 G4ThreeVector newDir = aR->GetMomentumChange();
405 newDir.rotateUz(dir);
408 }
409 //G4cout << "FillResult: Efinal= " << efinal << " status= "
410 // << theTotalResult->GetTrackStatus()
411 // << " fKill= " << fStopAndKill << G4endl;
412
413 // check secondaries
414 nICelectrons = 0;
415 G4int nSec = aR->GetNumberOfSecondaries();
417 G4double time0 = aT.GetGlobalTime();
418
419 for (G4int i = 0; i < nSec; ++i) {
420 G4DynamicParticle* dynParticle = aR->GetSecondary(i)->GetParticle();
421
422 // apply rotation
423 G4ThreeVector newDir = dynParticle->GetMomentumDirection();
424 newDir.rotateUz(dir);
425 dynParticle->SetMomentumDirection(newDir);
426
427 // check if secondary is on the mass shell
428 const G4ParticleDefinition* part = dynParticle->GetDefinition();
429 G4double mass = part->GetPDGMass();
430 G4double dmass= dynParticle->GetMass();
431 const G4double delta_mass_lim = 1.0*CLHEP::keV;
432 const G4double delta_ekin = 0.001*CLHEP::eV;
433 if(std::abs(dmass - mass) > delta_mass_lim) {
434 G4double e = std::max(dynParticle->GetKineticEnergy() + dmass - mass, delta_ekin);
437 ed << "TrackID= "<< aT.GetTrackID()
438 << " " << aT.GetParticleDefinition()->GetParticleName()
439 << " Target Z= " << targetNucleus.GetZ_asInt() << " A= "
441 << " Ekin(GeV)= " << aT.GetKineticEnergy()/CLHEP::GeV
442 << "\n Secondary is out of mass shell: " << part->GetParticleName()
443 << " EkinNew(MeV)= " << e
444 << " DeltaMass(MeV)= " << dmass - mass << G4endl;
445 G4Exception("G4HadronicProcess::FillResults", "had012", JustWarning, ed);
446 }
447 dynParticle->SetKineticEnergy(e);
448 dynParticle->SetMass(mass);
449 }
450 G4int idModel = aR->GetSecondary(i)->GetCreatorModelID();
451 if(part->GetPDGEncoding() == 11) { ++nICelectrons; }
452
453 // time of interaction starts from zero + global time
454 G4double time = std::max(aR->GetSecondary(i)->GetTime(), 0.0) + time0;
455
456 G4Track* track = new G4Track(dynParticle, time, aT.GetPosition());
457 track->SetCreatorModelID(idModel);
458 G4double newWeight = fWeight*aR->GetSecondary(i)->GetWeight();
459 track->SetWeight(newWeight);
463 G4double e = dynParticle->GetKineticEnergy();
464 if (e == 0.0) {
466 DumpState(aT,"Secondary has zero energy",ed);
467 ed << "Secondary " << part->GetParticleName()
468 << G4endl;
469 G4Exception("G4HadronicProcess::FillResults", "had011",
470 JustWarning,ed);
471 }
472 }
473 }
474 aR->Clear();
475 // G4cout << "FillResults done nICe= " << nICelectrons << G4endl;
476}
477
479{
481}
482
484{
485 if (aScale <= 0.0) {
487 ed << " Wrong biasing factor " << aScale << " for " << GetProcessName();
488 G4Exception("G4HadronicProcess::BiasCrossSectionByFactor", "had010",
489 JustWarning, ed, "Cross-section bias is ignored");
490 } else {
491 aScaleFactor = aScale;
492 }
493}
494
496 const G4Nucleus &aNucleus,
497 G4HadFinalState * result)
498{
499 // check for catastrophic energy non-conservation
500 // to re-sample the interaction
501
503 G4double nuclearMass(0);
504 if (theModel) {
505
506 // Compute final-state total energy
507 G4double finalE(0.);
508 G4int nSec = result->GetNumberOfSecondaries();
509
510 nuclearMass = G4NucleiProperties::GetNuclearMass(aNucleus.GetA_asInt(),
511 aNucleus.GetZ_asInt());
512 if (result->GetStatusChange() != stopAndKill) {
513 // Interaction didn't complete, returned "do nothing" state
514 // and reset nucleus or the primary survived the interaction
515 // (e.g. electro-nuclear ) => keep nucleus
516 finalE=result->GetLocalEnergyDeposit() +
517 aPro.GetDefinition()->GetPDGMass() + result->GetEnergyChange();
518 if( nSec == 0 ){
519 // Since there are no secondaries, there is no recoil nucleus.
520 // To check energy balance we must neglect the initial nucleus too.
521 nuclearMass=0.0;
522 }
523 }
524 for (G4int i = 0; i < nSec; i++) {
525 G4DynamicParticle *pdyn=result->GetSecondary(i)->GetParticle();
526 finalE += pdyn->GetTotalEnergy();
527 G4double mass_pdg=pdyn->GetDefinition()->GetPDGMass();
528 G4double mass_dyn=pdyn->GetMass();
529 if ( std::abs(mass_pdg - mass_dyn) > 0.1*mass_pdg + 1.*MeV ) {
530 // If it is shortlived, then a difference less than 3 times the width is acceptable
531 if ( pdyn->GetDefinition()->IsShortLived() &&
532 std::abs(mass_pdg - mass_dyn) < 3.0*pdyn->GetDefinition()->GetPDGWidth() ) {
533 continue;
534 }
535 result->Clear();
536 result = nullptr;
538 desc << "Warning: Secondary with off-shell dynamic mass detected: "
539 << G4endl
540 << " " << pdyn->GetDefinition()->GetParticleName()
541 << ", PDG mass: " << mass_pdg << ", dynamic mass: "
542 << mass_dyn << G4endl
543 << (epReportLevel<0 ? "abort the event"
544 : "re-sample the interaction") << G4endl
545 << " Process / Model: " << GetProcessName()<< " / "
546 << theModel->GetModelName() << G4endl
547 << " Primary: " << aPro.GetDefinition()->GetParticleName()
548 << " (" << aPro.GetDefinition()->GetPDGEncoding() << "), "
549 << " E= " << aPro.Get4Momentum().e()
550 << ", target nucleus (" << aNucleus.GetZ_asInt() << ", "
551 << aNucleus.GetA_asInt() << ")" << G4endl;
552 G4Exception("G4HadronicProcess:CheckResult()", "had012",
554 // must return here.....
555 return result;
556 }
557 }
558 G4double deltaE= nuclearMass + aPro.GetTotalEnergy() - finalE;
559
560 std::pair<G4double, G4double> checkLevels =
561 theModel->GetFatalEnergyCheckLevels(); // (relative, absolute)
562 if (std::abs(deltaE) > checkLevels.second &&
563 std::abs(deltaE) > checkLevels.first*aPro.GetKineticEnergy()){
564 // do not delete result, this is a pointer to a data member;
565 result->Clear();
566 result = nullptr;
568 desc << "Warning: Bad energy non-conservation detected, will "
569 << (epReportLevel<0 ? "abort the event"
570 : "re-sample the interaction") << G4endl
571 << " Process / Model: " << GetProcessName()<< " / "
572 << theModel->GetModelName() << G4endl
573 << " Primary: " << aPro.GetDefinition()->GetParticleName()
574 << " (" << aPro.GetDefinition()->GetPDGEncoding() << "), "
575 << " E= " << aPro.Get4Momentum().e()
576 << ", target nucleus (" << aNucleus.GetZ_asInt() << ", "
577 << aNucleus.GetA_asInt() << ")" << G4endl
578 << " E(initial - final) = " << deltaE << " MeV." << G4endl;
579 G4Exception("G4HadronicProcess:CheckResult()", "had012",
581 }
582 }
583 return result;
584}
585
586void
588 const G4Nucleus& aNucleus)
589{
590 G4int target_A=aNucleus.GetA_asInt();
591 G4int target_Z=aNucleus.GetZ_asInt();
592 G4double targetMass = G4NucleiProperties::GetNuclearMass(target_A,target_Z);
593 G4LorentzVector target4mom(0, 0, 0, targetMass
595
596 G4LorentzVector projectile4mom = aTrack.GetDynamicParticle()->Get4Momentum();
597 G4int track_A = aTrack.GetDefinition()->GetBaryonNumber();
598 G4int track_Z = G4lrint(aTrack.GetDefinition()->GetPDGCharge());
599
600 G4int initial_A = target_A + track_A;
601 G4int initial_Z = target_Z + track_Z - nICelectrons;
602
603 G4LorentzVector initial4mom = projectile4mom + target4mom;
604
605 // Compute final-state momentum for scattering and "do nothing" results
606 G4LorentzVector final4mom;
607 G4int final_A(0), final_Z(0);
608
610 if (theTotalResult->GetTrackStatus() != fStopAndKill) { // If it is Alive
611 // Either interaction didn't complete, returned "do nothing" state
612 // or the primary survived the interaction (e.g. electro-nucleus )
613
614 // Interaction didn't complete, returned "do nothing" state
615 // - or suppressed recoil (e.g. Neutron elastic )
616 final4mom = initial4mom;
617 final_A = initial_A;
618 final_Z = initial_Z;
619 if (nSec > 0) {
620 // The primary remains in final state (e.g. electro-nucleus )
621 // Use the final energy / momentum
624 G4double mass = aTrack.GetDefinition()->GetPDGMass();
625 G4double ptot = std::sqrt(ekin*(ekin + 2*mass));
626 final4mom.set(ptot*v.x(), ptot*v.y(), ptot*v.z(), mass + ekin);
627 final_A = track_A;
628 final_Z = track_Z;
629 // Expect that the target nucleus will have interacted,
630 // and its products, including recoil, will be included in secondaries.
631 }
632 }
633 if( nSec > 0 ) {
634 G4Track* sec;
635
636 for (G4int i = 0; i < nSec; i++) {
638 final4mom += sec->GetDynamicParticle()->Get4Momentum();
639 final_A += sec->GetDefinition()->GetBaryonNumber();
640 final_Z += G4lrint(sec->GetDefinition()->GetPDGCharge());
641 }
642 }
643
644 // Get level-checking information (used to cut-off relative checks)
645 G4String processName = GetProcessName();
647 G4String modelName("none");
648 if (theModel) modelName = theModel->GetModelName();
649 std::pair<G4double, G4double> checkLevels = epCheckLevels;
650 if (!levelsSetByProcess) {
651 if (theModel) checkLevels = theModel->GetEnergyMomentumCheckLevels();
652 checkLevels.first= std::min(checkLevels.first, epCheckLevels.first);
653 checkLevels.second=std::min(checkLevels.second, epCheckLevels.second);
654 }
655
656 // Compute absolute total-energy difference, and relative kinetic-energy
657 G4bool checkRelative = (aTrack.GetKineticEnergy() > checkLevels.second);
658
659 G4LorentzVector diff = initial4mom - final4mom;
660 G4double absolute = diff.e();
661 G4double relative = checkRelative ? absolute/aTrack.GetKineticEnergy() : 0.;
662
663 G4double absolute_mom = diff.vect().mag();
664 G4double relative_mom = checkRelative ? absolute_mom/aTrack.GetMomentum().mag() : 0.;
665
666 // Evaluate relative and absolute conservation
667 G4bool relPass = true;
668 G4String relResult = "pass";
669 if ( std::abs(relative) > checkLevels.first
670 || std::abs(relative_mom) > checkLevels.first) {
671 relPass = false;
672 relResult = checkRelative ? "fail" : "N/A";
673 }
674
675 G4bool absPass = true;
676 G4String absResult = "pass";
677 if ( std::abs(absolute) > checkLevels.second
678 || std::abs(absolute_mom) > checkLevels.second ) {
679 absPass = false ;
680 absResult = "fail";
681 }
682
683 G4bool chargePass = true;
684 G4String chargeResult = "pass";
685 if ( (initial_A-final_A)!=0
686 || (initial_Z-final_Z)!=0 ) {
687 chargePass = checkLevels.second < DBL_MAX ? false : true;
688 chargeResult = "fail";
689 }
690
691 G4bool conservationPass = (relPass || absPass) && chargePass;
692
693 std::stringstream Myout;
694 G4bool Myout_notempty(false);
695 // Options for level of reporting detail:
696 // 0. off
697 // 1. report only when E/p not conserved
698 // 2. report regardless of E/p conservation
699 // 3. report only when E/p not conserved, with model names, process names, and limits
700 // 4. report regardless of E/p conservation, with model names, process names, and limits
701 // negative -1.., as above, but send output to stderr
702
703 if( std::abs(epReportLevel) == 4
704 || ( std::abs(epReportLevel) == 3 && ! conservationPass ) ){
705 Myout << " Process: " << processName << " , Model: " << modelName << G4endl;
706 Myout << " Primary: " << aTrack.GetParticleDefinition()->GetParticleName()
707 << " (" << aTrack.GetParticleDefinition()->GetPDGEncoding() << "),"
708 << " E= " << aTrack.GetDynamicParticle()->Get4Momentum().e()
709 << ", target nucleus (" << aNucleus.GetZ_asInt() << ","
710 << aNucleus.GetA_asInt() << ")" << G4endl;
711 Myout_notempty=true;
712 }
713 if ( std::abs(epReportLevel) == 4
714 || std::abs(epReportLevel) == 2
715 || ! conservationPass ){
716
717 Myout << " "<< relResult <<" relative, limit " << checkLevels.first << ", values E/T(0) = "
718 << relative << " p/p(0)= " << relative_mom << G4endl;
719 Myout << " "<< absResult << " absolute, limit (MeV) " << checkLevels.second/MeV << ", values E / p (MeV) = "
720 << absolute/MeV << " / " << absolute_mom/MeV << " 3mom: " << (diff.vect())*1./MeV << G4endl;
721 Myout << " "<< chargeResult << " charge/baryon number balance " << (initial_Z-final_Z) << " / " << (initial_A-final_A) << " "<< G4endl;
722 Myout_notempty=true;
723
724 }
725 Myout.flush();
726 if ( Myout_notempty ) {
727 if (epReportLevel > 0) G4cout << Myout.str()<< G4endl;
728 else if (epReportLevel < 0) G4cerr << Myout.str()<< G4endl;
729 }
730}
731
733 const G4String& method,
735{
736 ed << "Unrecoverable error in the method " << method << " of "
737 << GetProcessName() << G4endl;
738 ed << "TrackID= "<< aTrack.GetTrackID() << " ParentID= "
739 << aTrack.GetParentID()
740 << " " << aTrack.GetParticleDefinition()->GetParticleName()
741 << G4endl;
742 ed << "Ekin(GeV)= " << aTrack.GetKineticEnergy()/CLHEP::GeV
743 << "; direction= " << aTrack.GetMomentumDirection() << G4endl;
744 ed << "Position(mm)= " << aTrack.GetPosition()/CLHEP::mm << ";";
745
746 if (aTrack.GetMaterial()) {
747 ed << " material " << aTrack.GetMaterial()->GetName();
748 }
749 ed << G4endl;
750
751 if (aTrack.GetVolume()) {
752 ed << "PhysicalVolume <" << aTrack.GetVolume()->GetName()
753 << ">" << G4endl;
754 }
755}
756
758{
760}
761
763{
765}
766
767std::vector<G4HadronicInteraction*>&
769{
771}
772
775{
776 std::vector<G4HadronicInteraction*>& list
778 for (size_t li=0; li<list.size(); li++) {
779 if (list[li]->GetModelName() == modelName) return list[li];
780 }
781 return nullptr;
782}
@ 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
static const char * G4Hadronic_Random_File
G4ProcessType
@ fHadronic
static constexpr double MeV
Definition: G4SIunits.hh:200
@ 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:127
G4int GetZasInt() const
Definition: G4Element.hh:132
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
G4double GetTime() const
G4int GetCreatorModelID() 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
G4CrossSectionDataStore * theCrossSectionDataStore
G4double theInitialNumberOfInteractionLength
void CheckEnergyMomentumConservation(const G4Track &, const G4Nucleus &)
G4HadronicInteraction * GetHadronicModel(const G4String &)
G4EnergyRangeManager theEnergyRangeManager
G4double XBiasSecondaryWeight()
G4HadronicProcessStore * theProcessStore
G4HadronicInteraction * theInteraction
G4double XBiasSurvivalProbability()
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
void GetEnergyMomentumCheckEnvvars()
void DumpPhysicsTable(const G4ParticleDefinition &p)
void MultiplyCrossSectionBy(G4double factor)
std::pair< G4double, G4double > epCheckLevels
void RegisterMe(G4HadronicInteraction *a)
static G4KaonZeroLong * Definition()
static G4KaonZeroShort * Definition()
static G4KaonZero * Definition()
Definition: G4KaonZero.cc:52
const G4String & GetName() const
Definition: G4Material.hh:173
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
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
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
G4ThreeVector GetMomentum() const
G4Material * GetMaterial() const
G4ParticleDefinition * GetDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetCreatorModelID(const G4int id)
G4int GetParentID() const
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
static constexpr double mm
Definition: SystemOfUnits.h:96
static constexpr double electron_mass_c2
static constexpr double GeV
static constexpr double keV
static constexpr double eV
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62