Geant4 11.1.1
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
68
69#include "G4NistManager.hh"
71#include "G4HadXSHelper.hh"
72#include "G4Threading.hh"
73#include "G4Exp.hh"
74
75#include <typeinfo>
76#include <sstream>
77#include <iostream>
78
79constexpr G4double lambdaFactor = 0.8;
81
82// File-scope variable to capture environment variable at startup
83static const char* G4Hadronic_Random_File = std::getenv("G4HADRONIC_RANDOM_FILE");
84
85//////////////////////////////////////////////////////////////////
86
88 G4ProcessType procType)
89 : G4VDiscreteProcess(processName, procType)
90{
91 SetProcessSubType(fHadronInelastic); // Default unless subclass changes
92 InitialiseLocal();
93}
94
96 G4HadronicProcessType aHadSubType)
97 : G4VDiscreteProcess(processName, fHadronic)
98{
99 SetProcessSubType(aHadSubType);
100 InitialiseLocal();
101}
102
104{
105 theProcessStore->DeRegister(this);
106 delete theTotalResult;
108 if(isMaster) {
109 delete fXSpeaks;
110 delete theEnergyOfCrossSectionMax;
111 }
112}
113
114void G4HadronicProcess::InitialiseLocal() {
118 theProcessStore = G4HadronicProcessStore::Instance();
119 theProcessStore->Register(this);
120 minKinEnergy = 1*CLHEP::MeV;
121 epCheckLevels.first = DBL_MAX;
122 epCheckLevels.second = DBL_MAX;
123 GetEnergyMomentumCheckEnvvars();
124 unitVector.set(0.0, 0.0, 0.1);
125 if(G4Threading::IsWorkerThread()) { isMaster = false; }
126}
127
128void G4HadronicProcess::GetEnergyMomentumCheckEnvvars() {
129 if ( std::getenv("G4Hadronic_epReportLevel") ) {
130 epReportLevel = std::strtol(std::getenv("G4Hadronic_epReportLevel"),0,10);
131 }
132 if ( std::getenv("G4Hadronic_epCheckRelativeLevel") ) {
133 epCheckLevels.first = std::strtod(std::getenv("G4Hadronic_epCheckRelativeLevel"),0);
134 }
135 if ( std::getenv("G4Hadronic_epCheckAbsoluteLevel") ) {
136 epCheckLevels.second = std::strtod(std::getenv("G4Hadronic_epCheckAbsoluteLevel"),0);
137 }
138}
139
141{
142 if(nullptr == a) { return; }
143 theEnergyRangeManager.RegisterMe( a );
145}
146
149 const G4Element * elm,
150 const G4Material* mat)
151{
152 if(nullptr == mat)
153 {
154 static const G4int nmax = 5;
155 if(nMatWarn < nmax) {
156 ++nMatWarn;
158 ed << "Cannot compute Element x-section for " << GetProcessName()
159 << " because no material defined \n"
160 << " Please, specify material pointer or define simple material"
161 << " for Z= " << elm->GetZasInt();
162 G4Exception("G4HadronicProcess::GetElementCrossSection", "had066",
163 JustWarning, ed);
164 }
165 }
166 return theCrossSectionDataStore->GetCrossSection(dp, elm, mat);
167}
168
170{
171 if(std::getenv("G4HadronicProcess_debug")) {
172 G4HadronicProcess_debug_flag = true;
173 }
174 if(nullptr == firstParticle) { firstParticle = &p; }
175 theProcessStore->RegisterParticle(this, &p);
176}
177
179{
180 if(firstParticle != &p) { return; }
181
183 theEnergyRangeManager.BuildPhysicsTable(p);
185
186 G4int subtype = GetProcessSubType();
187 if(useIntegralXS) {
188 if(subtype == fHadronInelastic) {
189 useIntegralXS = param->EnableIntegralInelasticXS();
190 } else if(subtype == fHadronElastic) {
191 useIntegralXS = param->EnableIntegralElasticXS();
192 }
193 }
195
196 // check particle for integral method
197 if(isMaster) {
198 G4double charge = p.GetPDGCharge()/eplus;
199 G4bool isLepton = (p.GetLeptonNumber() != 0);
200 G4bool ok = (p.GetAtomicNumber() != 0 || p.GetPDGMass() < GeV);
201
202 // select cross section shape
203 if(charge != 0.0 && useIntegralXS && !isLepton && ok) {
204 G4double tmax = param->GetMaxEnergy();
205 fXSType = (charge > 0.0) ? fHadIncreasing : fHadDecreasing;
206 currentParticle = firstParticle;
207 // initialisation in the master thread
208 G4int pdg = p.GetPDGEncoding();
209 if(std::abs(pdg) == 211) {
211 } else if(pdg == 321) {
213 } else if(pdg == 2212) {
215 }
216 delete theEnergyOfCrossSectionMax;
217 theEnergyOfCrossSectionMax = nullptr;
218 if(fXSType == fHadTwoPeaks) {
219 delete fXSpeaks;
220 fXSpeaks =
221 G4HadXSHelper::FillPeaksStructure(this, &p, minKinEnergy, tmax);
222 if(nullptr == fXSpeaks) {
224 }
225 }
226 if(fXSType == fHadOnePeak) {
227 theEnergyOfCrossSectionMax =
228 G4HadXSHelper::FindCrossSectionMax(this, &p, minKinEnergy, tmax);
229 if(nullptr == theEnergyOfCrossSectionMax) {
231 }
232 }
233 }
234 } else {
235 if(nullptr == masterProcess) {
236 masterProcess =
237 dynamic_cast<const G4HadronicProcess*>(GetMasterProcess());
238 }
239 if(nullptr == masterProcess) {
240 G4cout << "G4HadronicProcess::BuildPhysicsTable: for "
241 << GetProcessName() << " and " << p.GetParticleName()
242 << " fail due to undefined pointer to the master process"
243 << G4endl;
244 } else {
245 // initialisation in worker threads
246 fXSType = masterProcess->CrossSectionType();
247 fXSpeaks = masterProcess->TwoPeaksXS();
248 theEnergyOfCrossSectionMax =
249 masterProcess->EnergyOfCrossSectionMax();
250 }
251 }
252 if(isMaster && 1 < param->GetVerboseLevel()) {
253 G4cout << "G4HadronicProcess::BuildPhysicsTable: for "
254 << GetProcessName() << " and " << p.GetParticleName()
255 << " typeXS=" << fXSType << G4endl;
256 }
258}
259
261{
262 currentMat = nullptr;
263 currentParticle = track->GetDefinition();
264 fDynParticle = track->GetDynamicParticle();
266}
267
269 const G4Track& track,
270 G4double previousStepSize,
272{
274
275 const G4Material* mat = track.GetMaterial();
276 if(mat != currentMat) {
277 currentMat = mat;
279 matIdx = (G4int)track.GetMaterial()->GetIndex();
280 }
281 /*
282 G4cout << GetProcessName() << " E=" << track.GetKineticEnergy()
283 << " " << currentParticle->GetParticleName()
284 << " lastxs=" << theLastCrossSection
285 << " lastmfp=" << theMFP << G4endl;
286 */
287 UpdateCrossSectionAndMFP(track.GetKineticEnergy());
288 /*
289 G4cout << " xs=" << theLastCrossSection
290 << " mfp=" << theMFP << " nleft=" << theNumberOfInteractionLengthLeft
291 << G4endl;
292 */
293 // zero cross section
294 if(theLastCrossSection <= 0.0) {
297 return DBL_MAX;
298 }
299
300 // non-zero cross section
304 } else {
306 previousStepSize/currentInteractionLength;
309 }
312}
313
315 const G4Track &aTrack, G4double,
317{
320 return (xs > 0.0) ? 1.0/xs : DBL_MAX;
321}
322
325{
327
328 //G4cout << "PostStepDoIt " << aTrack.GetDefinition()->GetParticleName()
329 // << " Ekin= " << aTrack.GetKineticEnergy() << G4endl;
330 // if primary is not Alive then do nothing
332 theTotalResult->Initialize(aTrack);
333 fWeight = aTrack.GetWeight();
335 if(aTrack.GetTrackStatus() != fAlive) { return theTotalResult; }
336
337 // Find cross section at end of step and check if <= 0
338 //
339 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
340 const G4Material* aMaterial = aTrack.GetMaterial();
341
342 // check only for charged particles
343 if(fXSType != fHadNoIntegral) {
346 theCrossSectionDataStore->ComputeCrossSection(aParticle,aMaterial);
347 //G4cout << "xs=" << xs << " xs0=" << theLastCrossSection
348 // << " " << aMaterial->GetName() << G4endl;
350 // No interaction
351 return theTotalResult;
352 }
353 }
354
355 const G4Element* anElement =
356 theCrossSectionDataStore->SampleZandA(aParticle,aMaterial,targetNucleus);
357
358 // Next check for illegal track status
359 //
360 if (aTrack.GetTrackStatus() != fAlive &&
361 aTrack.GetTrackStatus() != fSuspend) {
362 if (aTrack.GetTrackStatus() == fStopAndKill ||
366 ed << "G4HadronicProcess: track in unusable state - "
367 << aTrack.GetTrackStatus() << G4endl;
368 ed << "G4HadronicProcess: returning unchanged track " << G4endl;
369 DumpState(aTrack,"PostStepDoIt",ed);
370 G4Exception("G4HadronicProcess::PostStepDoIt", "had004", JustWarning, ed);
371 }
372 // No warning for fStopButAlive which is a legal status here
373 return theTotalResult;
374 }
375
376 // Initialize the hadronic projectile from the track
377 thePro.Initialise(aTrack);
378
379 theInteraction = ChooseHadronicInteraction(thePro, targetNucleus,
380 aMaterial, anElement);
381 if(nullptr == theInteraction) {
383 ed << "Target element "<<anElement->GetName()<<" Z= "
384 << targetNucleus.GetZ_asInt() << " A= "
385 << targetNucleus.GetA_asInt() << G4endl;
386 DumpState(aTrack,"ChooseHadronicInteraction",ed);
387 ed << " No HadronicInteraction found out" << G4endl;
388 G4Exception("G4HadronicProcess::PostStepDoIt", "had005",
389 FatalException, ed);
390 return theTotalResult;
391 }
392
393 G4HadFinalState* result = nullptr;
394 G4int reentryCount = 0;
395 /*
396 G4cout << "### " << aParticle->GetDefinition()->GetParticleName()
397 << " Ekin(MeV)= " << aParticle->GetKineticEnergy()
398 << " Z= " << targetNucleus.GetZ_asInt()
399 << " A= " << targetNucleus.GetA_asInt()
400 << " by " << theInteraction->GetModelName()
401 << G4endl;
402 */
403 do
404 {
405 try
406 {
407 // Save random engine if requested for debugging
408 if (G4Hadronic_Random_File) {
409 CLHEP::HepRandom::saveEngineStatus(G4Hadronic_Random_File);
410 }
411 // Call the interaction
412 result = theInteraction->ApplyYourself( thePro, targetNucleus);
413 ++reentryCount;
414 }
415 catch(G4HadronicException & aR)
416 {
418 aR.Report(ed);
419 ed << "Call for " << theInteraction->GetModelName() << G4endl;
420 ed << "Target element "<<anElement->GetName()<<" Z= "
421 << targetNucleus.GetZ_asInt()
422 << " A= " << targetNucleus.GetA_asInt() << G4endl;
423 DumpState(aTrack,"ApplyYourself",ed);
424 ed << " ApplyYourself failed" << G4endl;
425 G4Exception("G4HadronicProcess::PostStepDoIt", "had006", FatalException,
426 ed);
427 }
428
429 // Check the result for catastrophic energy non-conservation
430 result = CheckResult(thePro, targetNucleus, result);
431
432 if(reentryCount>100) {
434 ed << "Call for " << theInteraction->GetModelName() << G4endl;
435 ed << "Target element "<<anElement->GetName()<<" Z= "
436 << targetNucleus.GetZ_asInt()
437 << " A= " << targetNucleus.GetA_asInt() << G4endl;
438 DumpState(aTrack,"ApplyYourself",ed);
439 ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
440 G4Exception("G4HadronicProcess::PostStepDoIt", "had006", FatalException,
441 ed);
442 }
443 }
444 while(!result); /* Loop checking, 30-Oct-2015, G.Folger */
445
446 // Check whether kaon0 or anti_kaon0 are present between the secondaries:
447 // if this is the case, transform them into either kaon0S or kaon0L,
448 // with equal, 50% probability, keeping their dynamical masses (and
449 // the other kinematical properties).
450 // When this happens - very rarely - a "JustWarning" exception is thrown.
451 G4int nSec = (G4int)result->GetNumberOfSecondaries();
452 if ( nSec > 0 ) {
453 for ( G4int i = 0; i < nSec; ++i ) {
454 auto dynamicParticle = result->GetSecondary(i)->GetParticle();
455 auto part = dynamicParticle->GetParticleDefinition();
456 if ( part == G4KaonZero::Definition() ||
457 part == G4AntiKaonZero::Definition() ) {
458 G4ParticleDefinition* newPart;
459 if( G4UniformRand() > 0.5 ) { newPart = G4KaonZeroShort::Definition(); }
460 else { newPart = G4KaonZeroLong::Definition(); }
461 dynamicParticle->SetDefinition( newPart );
462 if(nKaonWarn < 5) {
463 ++nKaonWarn;
465 ed << " Hadronic model " << theInteraction->GetModelName() << G4endl;
466 ed << " created " << part->GetParticleName() << G4endl;
467 ed << " -> forced to be " << newPart->GetParticleName() << G4endl;
468 G4Exception( "G4HadronicProcess::PostStepDoIt", "had007", JustWarning, ed );
469 }
470 }
471 }
472 }
473
475 FillResult(result, aTrack);
476
477 if (epReportLevel != 0) {
478 CheckEnergyMomentumConservation(aTrack, targetNucleus);
479 }
480 //G4cout << "PostStepDoIt done nICelectrons= " << nICelectrons << G4endl;
481 return theTotalResult;
482}
483
484void G4HadronicProcess::ProcessDescription(std::ostream& outFile) const
485{
486 outFile << "The description for this process has not been written yet.\n";
487}
488
489G4double G4HadronicProcess::XBiasSurvivalProbability()
490{
492 G4double biasedProbability = 1.-G4Exp(-nLTraversed);
493 G4double realProbability = 1-G4Exp(-nLTraversed/aScaleFactor);
494 G4double result = (biasedProbability-realProbability)/biasedProbability;
495 return result;
496}
497
498G4double G4HadronicProcess::XBiasSecondaryWeight()
499{
501 G4double result =
502 1./aScaleFactor*G4Exp(-nLTraversed/aScaleFactor*(1-1./aScaleFactor));
503 return result;
504}
505
506void
508{
510 const G4ThreeVector& dir = aT.GetMomentumDirection();
511
512 G4double efinal = std::max(aR->GetEnergyChange(), 0.0);
513
514 // check status of primary
515 if(aR->GetStatusChange() == stopAndKill) {
518
519 // check its final energy
520 } else if(0.0 == efinal) {
523 ->GetAtRestProcessVector()->size() > 0)
526
527 // primary is not killed apply rotation and Lorentz transformation
528 } else {
530 G4ThreeVector newDir = aR->GetMomentumChange();
531 newDir.rotateUz(dir);
534 }
535 //G4cout << "FillResult: Efinal= " << efinal << " status= "
536 // << theTotalResult->GetTrackStatus()
537 // << " fKill= " << fStopAndKill << G4endl;
538
539 // check secondaries
540 nICelectrons = 0;
541 G4int nSec = (G4int)aR->GetNumberOfSecondaries();
543 G4double time0 = aT.GetGlobalTime();
544
545 for (G4int i = 0; i < nSec; ++i) {
546 G4DynamicParticle* dynParticle = aR->GetSecondary(i)->GetParticle();
547
548 // apply rotation
549 G4ThreeVector newDir = dynParticle->GetMomentumDirection();
550 newDir.rotateUz(dir);
551 dynParticle->SetMomentumDirection(newDir);
552
553 // check if secondary is on the mass shell
554 const G4ParticleDefinition* part = dynParticle->GetDefinition();
555 G4double mass = part->GetPDGMass();
556 G4double dmass= dynParticle->GetMass();
557 const G4double delta_mass_lim = 1.0*CLHEP::keV;
558 const G4double delta_ekin = 0.001*CLHEP::eV;
559 if(std::abs(dmass - mass) > delta_mass_lim) {
560 G4double e =
561 std::max(dynParticle->GetKineticEnergy() + dmass - mass, delta_ekin);
562 if(G4HadronicProcess_debug_flag) {
564 ed << "TrackID= "<< aT.GetTrackID()
565 << " " << aT.GetParticleDefinition()->GetParticleName()
566 << " Target Z= " << targetNucleus.GetZ_asInt() << " A= "
567 << targetNucleus.GetA_asInt()
568 << " Ekin(GeV)= " << aT.GetKineticEnergy()/CLHEP::GeV
569 << "\n Secondary is out of mass shell: " << part->GetParticleName()
570 << " EkinNew(MeV)= " << e
571 << " DeltaMass(MeV)= " << dmass - mass << G4endl;
572 G4Exception("G4HadronicProcess::FillResults", "had012", JustWarning, ed);
573 }
574 dynParticle->SetKineticEnergy(e);
575 dynParticle->SetMass(mass);
576 }
577 G4int idModel = aR->GetSecondary(i)->GetCreatorModelID();
578 if(part->GetPDGEncoding() == 11) { ++nICelectrons; }
579
580 // time of interaction starts from zero + global time
581 G4double time = std::max(aR->GetSecondary(i)->GetTime(), 0.0) + time0;
582
583 G4Track* track = new G4Track(dynParticle, time, aT.GetPosition());
584 track->SetCreatorModelID(idModel);
587 G4double newWeight = fWeight*aR->GetSecondary(i)->GetWeight();
588 track->SetWeight(newWeight);
591 if (G4HadronicProcess_debug_flag) {
592 G4double e = dynParticle->GetKineticEnergy();
593 if (e == 0.0) {
595 DumpState(aT,"Secondary has zero energy",ed);
596 ed << "Secondary " << part->GetParticleName()
597 << G4endl;
598 G4Exception("G4HadronicProcess::FillResults", "had011",
599 JustWarning,ed);
600 }
601 }
602 }
603 aR->Clear();
604 // G4cout << "FillResults done nICe= " << nICelectrons << G4endl;
605}
606
608{
610}
611
613{
614 if (aScale <= 0.0) {
616 ed << " Wrong biasing factor " << aScale << " for " << GetProcessName();
617 G4Exception("G4HadronicProcess::BiasCrossSectionByFactor", "had010",
618 JustWarning, ed, "Cross-section bias is ignored");
619 } else {
620 aScaleFactor = aScale;
621 }
622}
623
625 const G4Nucleus &aNucleus,
626 G4HadFinalState * result)
627{
628 // check for catastrophic energy non-conservation
629 // to re-sample the interaction
631 G4double nuclearMass(0);
632 if (nullptr != theModel) {
633
634 // Compute final-state total energy
635 G4double finalE(0.);
636 G4int nSec = (G4int)result->GetNumberOfSecondaries();
637
638 nuclearMass = G4NucleiProperties::GetNuclearMass(aNucleus.GetA_asInt(),
639 aNucleus.GetZ_asInt());
640 if (result->GetStatusChange() != stopAndKill) {
641 // Interaction didn't complete, returned "do nothing" state
642 // and reset nucleus or the primary survived the interaction
643 // (e.g. electro-nuclear ) => keep nucleus
644 finalE=result->GetLocalEnergyDeposit() +
645 aPro.GetDefinition()->GetPDGMass() + result->GetEnergyChange();
646 if( nSec == 0 ){
647 // Since there are no secondaries, there is no recoil nucleus.
648 // To check energy balance we must neglect the initial nucleus too.
649 nuclearMass=0.0;
650 }
651 }
652 for (G4int i = 0; i < nSec; ++i) {
653 G4DynamicParticle *pdyn=result->GetSecondary(i)->GetParticle();
654 finalE += pdyn->GetTotalEnergy();
655 G4double mass_pdg=pdyn->GetDefinition()->GetPDGMass();
656 G4double mass_dyn=pdyn->GetMass();
657 if ( std::abs(mass_pdg - mass_dyn) > 0.1*mass_pdg + 1.*MeV ) {
658 // If it is shortlived, then a difference less than 3 times the width is acceptable
659 if ( pdyn->GetDefinition()->IsShortLived() &&
660 std::abs(mass_pdg - mass_dyn) < 3.0*pdyn->GetDefinition()->GetPDGWidth() ) {
661 continue;
662 }
663 result->Clear();
664 result = nullptr;
666 desc << "Warning: Secondary with off-shell dynamic mass detected: "
667 << G4endl
668 << " " << pdyn->GetDefinition()->GetParticleName()
669 << ", PDG mass: " << mass_pdg << ", dynamic mass: "
670 << mass_dyn << G4endl
671 << (epReportLevel<0 ? "abort the event"
672 : "re-sample the interaction") << G4endl
673 << " Process / Model: " << GetProcessName()<< " / "
674 << theModel->GetModelName() << G4endl
675 << " Primary: " << aPro.GetDefinition()->GetParticleName()
676 << " (" << aPro.GetDefinition()->GetPDGEncoding() << "), "
677 << " E= " << aPro.Get4Momentum().e()
678 << ", target nucleus (" << aNucleus.GetZ_asInt() << ", "
679 << aNucleus.GetA_asInt() << ")" << G4endl;
680 G4Exception("G4HadronicProcess:CheckResult()", "had012",
682 // must return here.....
683 return result;
684 }
685 }
686 G4double deltaE= nuclearMass + aPro.GetTotalEnergy() - finalE;
687
688 std::pair<G4double, G4double> checkLevels =
689 theModel->GetFatalEnergyCheckLevels(); // (relative, absolute)
690 if (std::abs(deltaE) > checkLevels.second &&
691 std::abs(deltaE) > checkLevels.first*aPro.GetKineticEnergy()){
692 // do not delete result, this is a pointer to a data member;
693 result->Clear();
694 result = nullptr;
696 desc << "Warning: Bad energy non-conservation detected, will "
697 << (epReportLevel<0 ? "abort the event"
698 : "re-sample the interaction") << G4endl
699 << " Process / Model: " << GetProcessName()<< " / "
700 << theModel->GetModelName() << G4endl
701 << " Primary: " << aPro.GetDefinition()->GetParticleName()
702 << " (" << aPro.GetDefinition()->GetPDGEncoding() << "), "
703 << " E= " << aPro.Get4Momentum().e()
704 << ", target nucleus (" << aNucleus.GetZ_asInt() << ", "
705 << aNucleus.GetA_asInt() << ")" << G4endl
706 << " E(initial - final) = " << deltaE << " MeV." << G4endl;
707 G4Exception("G4HadronicProcess:CheckResult()", "had012",
709 }
710 }
711 return result;
712}
713
714void
716 const G4Nucleus& aNucleus)
717{
718 G4int target_A=aNucleus.GetA_asInt();
719 G4int target_Z=aNucleus.GetZ_asInt();
720 G4double targetMass = G4NucleiProperties::GetNuclearMass(target_A,target_Z);
721 G4LorentzVector target4mom(0, 0, 0, targetMass
722 + nICelectrons*CLHEP::electron_mass_c2);
723
724 G4LorentzVector projectile4mom = aTrack.GetDynamicParticle()->Get4Momentum();
725 G4int track_A = aTrack.GetDefinition()->GetBaryonNumber();
726 G4int track_Z = G4lrint(aTrack.GetDefinition()->GetPDGCharge());
727
728 G4int initial_A = target_A + track_A;
729 G4int initial_Z = target_Z + track_Z - nICelectrons;
730
731 G4LorentzVector initial4mom = projectile4mom + target4mom;
732
733 // Compute final-state momentum for scattering and "do nothing" results
734 G4LorentzVector final4mom;
735 G4int final_A(0), final_Z(0);
736
738 if (theTotalResult->GetTrackStatus() != fStopAndKill) { // If it is Alive
739 // Either interaction didn't complete, returned "do nothing" state
740 // or the primary survived the interaction (e.g. electro-nucleus )
741
742 // Interaction didn't complete, returned "do nothing" state
743 // - or suppressed recoil (e.g. Neutron elastic )
744 final4mom = initial4mom;
745 final_A = initial_A;
746 final_Z = initial_Z;
747 if (nSec > 0) {
748 // The primary remains in final state (e.g. electro-nucleus )
749 // Use the final energy / momentum
752 G4double mass = aTrack.GetDefinition()->GetPDGMass();
753 G4double ptot = std::sqrt(ekin*(ekin + 2*mass));
754 final4mom.set(ptot*v.x(), ptot*v.y(), ptot*v.z(), mass + ekin);
755 final_A = track_A;
756 final_Z = track_Z;
757 // Expect that the target nucleus will have interacted,
758 // and its products, including recoil, will be included in secondaries.
759 }
760 }
761 if( nSec > 0 ) {
762 G4Track* sec;
763
764 for (G4int i = 0; i < nSec; i++) {
766 final4mom += sec->GetDynamicParticle()->Get4Momentum();
767 final_A += sec->GetDefinition()->GetBaryonNumber();
768 final_Z += G4lrint(sec->GetDefinition()->GetPDGCharge());
769 }
770 }
771
772 // Get level-checking information (used to cut-off relative checks)
773 G4String processName = GetProcessName();
775 G4String modelName("none");
776 if (theModel) modelName = theModel->GetModelName();
777 std::pair<G4double, G4double> checkLevels = epCheckLevels;
778 if (!levelsSetByProcess) {
779 if (theModel) checkLevels = theModel->GetEnergyMomentumCheckLevels();
780 checkLevels.first= std::min(checkLevels.first, epCheckLevels.first);
781 checkLevels.second=std::min(checkLevels.second, epCheckLevels.second);
782 }
783
784 // Compute absolute total-energy difference, and relative kinetic-energy
785 G4bool checkRelative = (aTrack.GetKineticEnergy() > checkLevels.second);
786
787 G4LorentzVector diff = initial4mom - final4mom;
788 G4double absolute = diff.e();
789 G4double relative = checkRelative ? absolute/aTrack.GetKineticEnergy() : 0.;
790
791 G4double absolute_mom = diff.vect().mag();
792 G4double relative_mom = checkRelative ? absolute_mom/aTrack.GetMomentum().mag() : 0.;
793
794 // Evaluate relative and absolute conservation
795 G4bool relPass = true;
796 G4String relResult = "pass";
797 if ( std::abs(relative) > checkLevels.first
798 || std::abs(relative_mom) > checkLevels.first) {
799 relPass = false;
800 relResult = checkRelative ? "fail" : "N/A";
801 }
802
803 G4bool absPass = true;
804 G4String absResult = "pass";
805 if ( std::abs(absolute) > checkLevels.second
806 || std::abs(absolute_mom) > checkLevels.second ) {
807 absPass = false ;
808 absResult = "fail";
809 }
810
811 G4bool chargePass = true;
812 G4String chargeResult = "pass";
813 if ( (initial_A-final_A)!=0
814 || (initial_Z-final_Z)!=0 ) {
815 chargePass = checkLevels.second < DBL_MAX ? false : true;
816 chargeResult = "fail";
817 }
818
819 G4bool conservationPass = (relPass || absPass) && chargePass;
820
821 std::stringstream Myout;
822 G4bool Myout_notempty(false);
823 // Options for level of reporting detail:
824 // 0. off
825 // 1. report only when E/p not conserved
826 // 2. report regardless of E/p conservation
827 // 3. report only when E/p not conserved, with model names, process names, and limits
828 // 4. report regardless of E/p conservation, with model names, process names, and limits
829 // negative -1.., as above, but send output to stderr
830
831 if( std::abs(epReportLevel) == 4
832 || ( std::abs(epReportLevel) == 3 && ! conservationPass ) ){
833 Myout << " Process: " << processName << " , Model: " << modelName << G4endl;
834 Myout << " Primary: " << aTrack.GetParticleDefinition()->GetParticleName()
835 << " (" << aTrack.GetParticleDefinition()->GetPDGEncoding() << "),"
836 << " E= " << aTrack.GetDynamicParticle()->Get4Momentum().e()
837 << ", target nucleus (" << aNucleus.GetZ_asInt() << ","
838 << aNucleus.GetA_asInt() << ")" << G4endl;
839 Myout_notempty=true;
840 }
841 if ( std::abs(epReportLevel) == 4
842 || std::abs(epReportLevel) == 2
843 || ! conservationPass ){
844
845 Myout << " "<< relResult <<" relative, limit " << checkLevels.first << ", values E/T(0) = "
846 << relative << " p/p(0)= " << relative_mom << G4endl;
847 Myout << " "<< absResult << " absolute, limit (MeV) " << checkLevels.second/MeV << ", values E / p (MeV) = "
848 << absolute/MeV << " / " << absolute_mom/MeV << " 3mom: " << (diff.vect())*1./MeV << G4endl;
849 Myout << " "<< chargeResult << " charge/baryon number balance " << (initial_Z-final_Z) << " / " << (initial_A-final_A) << " "<< G4endl;
850 Myout_notempty=true;
851
852 }
853 Myout.flush();
854 if ( Myout_notempty ) {
855 if (epReportLevel > 0) G4cout << Myout.str()<< G4endl;
856 else if (epReportLevel < 0) G4cerr << Myout.str()<< G4endl;
857 }
858}
859
861 const G4String& method,
863{
864 ed << "Unrecoverable error in the method " << method << " of "
865 << GetProcessName() << G4endl;
866 ed << "TrackID= "<< aTrack.GetTrackID() << " ParentID= "
867 << aTrack.GetParentID()
868 << " " << aTrack.GetParticleDefinition()->GetParticleName()
869 << G4endl;
870 ed << "Ekin(GeV)= " << aTrack.GetKineticEnergy()/CLHEP::GeV
871 << "; direction= " << aTrack.GetMomentumDirection() << G4endl;
872 ed << "Position(mm)= " << aTrack.GetPosition()/CLHEP::mm << ";";
873
874 if (aTrack.GetMaterial()) {
875 ed << " material " << aTrack.GetMaterial()->GetName();
876 }
877 ed << G4endl;
878
879 if (aTrack.GetVolume()) {
880 ed << "PhysicalVolume <" << aTrack.GetVolume()->GetName()
881 << ">" << G4endl;
882 }
883}
884
886{
888}
889
891{
893}
894
895std::vector<G4HadronicInteraction*>&
897{
898 return theEnergyRangeManager.GetHadronicInteractionList();
899}
900
903{
904 std::vector<G4HadronicInteraction*>& list
905 = theEnergyRangeManager.GetHadronicInteractionList();
906 for (auto & mod : list) {
907 if (mod->GetModelName() == modelName) return mod;
908 }
909 return nullptr;
910}
911
914 const G4Material* mat,
915 const G4double kinEnergy)
916{
917 auto dp = new G4DynamicParticle(part, unitVector, kinEnergy);
919 delete dp;
920 return xs;
921}
922
923void G4HadronicProcess::RecomputeXSandMFP(const G4double kinEnergy)
924{
925 auto dp = new G4DynamicParticle(currentParticle, unitVector, kinEnergy);
928 theMFP = (theLastCrossSection > 0.0) ? 1.0/theLastCrossSection : DBL_MAX;
929 delete dp;
930}
931
932void G4HadronicProcess::UpdateCrossSectionAndMFP(const G4double e)
933{
934 if(fXSType == fHadNoIntegral) {
935 DefineXSandMFP();
936
937 } else if(fXSType == fHadIncreasing) {
939 mfpKinEnergy = e;
940 ComputeXSandMFP();
941 }
942
943 } else if(fXSType == fHadDecreasing) {
944 if(e < mfpKinEnergy && mfpKinEnergy > minKinEnergy) {
945 G4double e1 = std::max(e*lambdaFactor, minKinEnergy);
946 mfpKinEnergy = e1;
947 RecomputeXSandMFP(e1);
948 }
949
950 } else if(fXSType == fHadOnePeak) {
951 G4double epeak = (*theEnergyOfCrossSectionMax)[matIdx];
952 if(e <= epeak) {
954 mfpKinEnergy = e;
955 ComputeXSandMFP();
956 }
957 } else if(e < mfpKinEnergy) {
958 G4double e1 = std::max(epeak, e*lambdaFactor);
959 mfpKinEnergy = e1;
960 RecomputeXSandMFP(e1);
961 }
962
963 } else if(fXSType == fHadTwoPeaks) {
964 G4TwoPeaksHadXS* xs = (*fXSpeaks)[matIdx];
965 const G4double e1peak = xs->e1peak;
966
967 // below the 1st peak
968 if(e <= e1peak) {
970 mfpKinEnergy = e;
971 ComputeXSandMFP();
972 }
973 return;
974 }
975 const G4double e1deep = xs->e1deep;
976 // above the 1st peak, below the deep
977 if(e <= e1deep) {
978 if(mfpKinEnergy >= e1deep || e <= mfpKinEnergy) {
979 const G4double e1 = std::max(e1peak, e*lambdaFactor);
980 mfpKinEnergy = e1;
981 RecomputeXSandMFP(e1);
982 }
983 return;
984 }
985 const G4double e2peak = xs->e2peak;
986 // above the deep, below 2nd peak
987 if(e <= e2peak) {
989 mfpKinEnergy = e;
990 ComputeXSandMFP();
991 }
992 return;
993 }
994 const G4double e2deep = xs->e2deep;
995 // above the 2nd peak, below the deep
996 if(e <= e2deep) {
997 if(mfpKinEnergy >= e2deep || e <= mfpKinEnergy) {
998 const G4double e1 = std::max(e2peak, e*lambdaFactor);
999 mfpKinEnergy = e1;
1000 RecomputeXSandMFP(e1);
1001 }
1002 return;
1003 }
1004 const G4double e3peak = xs->e3peak;
1005 // above the deep, below 3d peak
1006 if(e <= e3peak) {
1008 mfpKinEnergy = e;
1009 ComputeXSandMFP();
1010 }
1011 return;
1012 }
1013 // above 3d peak
1014 if(e <= mfpKinEnergy) {
1015 const G4double e1 = std::max(e3peak, e*lambdaFactor);
1016 mfpKinEnergy = e1;
1017 RecomputeXSandMFP(e1);
1018 }
1019
1020 } else {
1021 DefineXSandMFP();
1022 }
1023}
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
@ FatalException
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
G4ForceCondition
@ NotForced
@ stopAndKill
@ fHadTwoPeaks
Definition: G4HadXSTypes.hh:49
@ fHadIncreasing
Definition: G4HadXSTypes.hh:46
@ fHadDecreasing
Definition: G4HadXSTypes.hh:47
@ fHadNoIntegral
Definition: G4HadXSTypes.hh:45
@ fHadOnePeak
Definition: G4HadXSTypes.hh:48
G4HadronicProcessType
@ fHadronElastic
@ fHadronInelastic
constexpr G4double lambdaFactor
constexpr G4double invLambdaFactor
G4double G4Log(G4double x)
Definition: G4Log.hh:227
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
void set(double x, double y, double z)
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 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()
G4int GetParentResonanceID() const
G4double GetWeight() const
const G4ParticleDefinition * GetParentResonanceDef() const
G4double GetTime() const
G4int GetCreatorModelID() const
static std::vector< G4TwoPeaksHadXS * > * FillPeaksStructure(G4HadronicProcess *, const G4ParticleDefinition *, const G4double tmin, const G4double tmax)
static std::vector< G4double > * FindCrossSectionMax(G4HadronicProcess *, const G4ParticleDefinition *, const G4double tmin, const G4double tmax)
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
G4bool EnableIntegralInelasticXS() const
static G4HadronicParameters * Instance()
G4bool EnableIntegralElasticXS() const
G4double GetMaxEnergy() 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
G4double ComputeCrossSection(const G4ParticleDefinition *, const G4Material *, const G4double kinEnergy)
void BiasCrossSectionByFactor(G4double aScale)
void StartTracking(G4Track *track) override
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)
std::vector< G4TwoPeaksHadXS * > * TwoPeaksXS() const
G4HadXSType CrossSectionType() const
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
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)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double, G4ForceCondition *) override
void RegisterMe(G4HadronicInteraction *a)
std::vector< G4double > * EnergyOfCrossSectionMax() const
static G4KaonZeroLong * Definition()
static G4KaonZeroShort * Definition()
static G4KaonZero * Definition()
Definition: G4KaonZero.cc:52
const G4String & GetName() const
Definition: G4Material.hh:172
size_t GetIndex() const
Definition: G4Material.hh:255
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
void Initialize(const G4Track &) override
const G4ThreeVector * GetMomentumDirection() const
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4ProcessManager * GetProcessManager() const
G4int GetAtomicNumber() 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)
void SetParentResonanceID(const G4int parentID)
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 SetParentResonanceDef(const G4ParticleDefinition *parent)
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
G4double currentInteractionLength
Definition: G4VProcess.hh:339
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:342
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:522
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:335
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4double GetTotalNumberOfInteractionLengthTraversed() const
Definition: G4VProcess.hh:441
G4int GetProcessSubType() const
Definition: G4VProcess.hh:404
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
G4bool IsWorkerThread()
Definition: G4Threading.cc:123
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62