Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LowEGammaNuclearModel.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// Physics model class G4LowEGammaNuclearModel
27// Created: 15 May 2019
28// Author V.Ivanchenko
29//
30//
31
33#include "G4SystemOfUnits.hh"
35#include "G4Fragment.hh"
36#include "G4NucleiProperties.hh"
37#include "G4DynamicParticle.hh"
38#include "G4HadSecondary.hh"
39#include "G4ReactionProduct.hh"
43#include "G4PreCompoundModel.hh"
44#include "G4ThreeVector.hh"
46
48 : G4HadronicInteraction("GammaNPreco"),lab4mom(0.,0.,0.,0.), secID(-1)
49{
51 SetMinEnergy( 0.0*CLHEP::GeV );
53
54 // reuse existing pre-compound model
57 fPreco = static_cast<G4PreCompoundModel*>(p);
58 if(!fPreco) { fPreco = new G4PreCompoundModel(); }
59}
60
62{}
63
65{}
66
68 const G4HadProjectile& aTrack, G4Nucleus& theNucleus)
69{
71
72 G4int A = theNucleus.GetA_asInt();
73 G4int Z = theNucleus.GetZ_asInt();
74
75 // Create initial state
76 lab4mom.set(0.,0.,0.,G4NucleiProperties::GetNuclearMass(A, Z));
77 lab4mom += aTrack.Get4Momentum();
78
79 G4Fragment frag(A, Z, lab4mom);
80
81 frag.SetCreatorModelID(secID);
82
83 if (verboseLevel > 1) {
84 G4cout << "G4LowEGammaNuclearModel::ApplyYourself initial G4Fragmet:"
85 << G4endl;
86 G4cout << frag << G4endl;
87 }
88 G4ReactionProductVector* res = fPreco->DeExcite(frag);
89
90 // secondaries produced
91 if(res) {
92
94 std::size_t nsec = res->size();
95 if (verboseLevel > 1) {
96 G4cout << "G4LowEGammaNuclearModel: " << nsec << " secondaries" << G4endl;
97 }
98 for(std::size_t i=0; i<nsec; ++i) {
99 if((*res)[i]) {
100 G4double ekin = (*res)[i]->GetKineticEnergy();
101 G4ThreeVector dir(0.,0.,1.);
102 if(ekin > 0.0) { dir = (*res)[i]->GetMomentum().unit(); }
103 G4HadSecondary* news = new G4HadSecondary(
104 new G4DynamicParticle((*res)[i]->GetDefinition(), dir, ekin));
105 news->SetTime((*res)[i]->GetTOF());
106 news->SetCreatorModelID(secID);
108 if (verboseLevel > 1) {
109 G4cout << i << ". " << (*res)[i]->GetDefinition()->GetParticleName()
110 << " Ekin(MeV)= " << ekin/MeV
111 << " dir: " << dir << G4endl;
112 }
113 delete (*res)[i];
114 delete news;
115 }
116 }
117 delete res;
118 }
119 return &theParticleChange;
120}
121
@ stopAndKill
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
void set(double x, double y, double z, double t)
void SetCreatorModelID(G4int value)
Definition: G4Fragment.hh:433
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
const G4LorentzVector & Get4Momentum() const
void SetTime(G4double aT)
void SetCreatorModelID(G4int id)
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
void SetMinEnergy(G4double anEnergy)
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) final
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
static G4int GetModelID(const G4int modelIndex)
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment) final