Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChiralInvariantPhaseSpace.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// Modified:
28// 24.08.10 A. Dotti ([email protected]) handle exceptions
29// thrown by Q4QEnvironment::Fragment retying interaction
30// 17.06.10 A. Dotti ([email protected]) handle case in which
31// Q4QEnvironment returns a 90000000 fragment (see code comments)
32
33// Created:
34// 16.01.08 V.Ivanchenko use initialization similar to other CHIPS models
35//
36
38#include "G4SystemOfUnits.hh"
39#include "G4ParticleTable.hh"
40#include "G4LorentzVector.hh"
41#include "G4DynamicParticle.hh"
42#include "G4IonTable.hh"
43#include "G4Neutron.hh"
44#include "G4Gamma.hh"
45
47{}
48
50{}
51
53 const G4HadProjectile& aTrack,
54 G4Nucleus& aTargetNucleus,
55 G4HadFinalState * aChange)
56{
57 G4HadFinalState * aResult;
58 if(aChange != 0)
59 {
60 aResult = aChange;
61 }
62 else
63 {
64 aResult = & theResult;
65 aResult->Clear();
67 }
68 //projectile properties needed in constructor of quasmon
69 G4LorentzVector proj4Mom;
70 proj4Mom = aTrack.Get4Momentum();
71 G4int projectilePDGCode = aTrack.GetDefinition()
73
74 //target properties needed in constructor of quasmon
75 G4int targetZ = aTargetNucleus.GetZ_asInt();
76 G4int targetA = aTargetNucleus.GetA_asInt();
77 G4int targetPDGCode = 90000000 + 1000*targetZ + (targetA-targetZ);
78 // NOT NECESSARY ______________
80 ->GetIonMass(targetZ, targetA);
81 G4LorentzVector targ4Mom(0.,0.,0.,targetMass);
82 // END OF NOT NECESSARY^^^^^^^^
83
84 // V.Ivanchenko set the same value as for other CHIPS models
85 G4int nop = 152;
86 //G4int nop = 164; // nuclear clusters up to A=21
87 G4double fractionOfSingleQuasiFreeNucleons = 0.4;
88 G4double fractionOfPairedQuasiFreeNucleons = 0.0;
89 if(targetA>27) fractionOfPairedQuasiFreeNucleons = 0.04;
90 G4double clusteringCoefficient = 4.;
91 G4double temperature = 180.;
92 G4double halfTheStrangenessOfSee = 0.1; // = s/d = s/u
93 G4double etaToEtaPrime = 0.3;
94
95 // construct and fragment the quasmon
96 //G4QCHIPSWorld aWorld(nop); // Create CHIPS World of nop particles
97 G4QCHIPSWorld::Get()->GetParticles(nop); // Create CHIPS World of nop particles
98 G4QNucleus::SetParameters(fractionOfSingleQuasiFreeNucleons,
99 fractionOfPairedQuasiFreeNucleons,
100 clusteringCoefficient);
101 G4Quasmon::SetParameters(temperature,
102 halfTheStrangenessOfSee,
103 etaToEtaPrime);
104 // G4QEnvironment::SetParameters(solidAngle);
105// G4cout << "Input info "<< projectilePDGCode << " "
106// << targetPDGCode <<" "
107// << 1./MeV*proj4Mom<<" "
108// << 1./MeV*targ4Mom << " "
109// << nop << G4endl;
110 G4QHadronVector projHV;
111 G4QHadron* iH = new G4QHadron(projectilePDGCode, 1./MeV*proj4Mom);
112 projHV.push_back(iH);
113
114 //AND
115 //A. Dotti 24 Aug. : Trying to handle situation when G4QEnvironment::Fragment throws an exception
116 // Seen by ATLAS for gamma+Nuclear interaction for [email protected] on Al
117 // The poor-man solution here is to re-try interaction if a G4QException is catched
118 // Warning: G4QExcpetion does NOT inherit from base class G4HadException
119 G4QHadronVector* output=0;
120
121 bool retry=true;
122 int retryes=0;
123 int maxretries=10;
124
125 while ( retry && retryes < maxretries )
126 {
127 G4QEnvironment* pan= new G4QEnvironment(projHV, targetPDGCode);
128
129 try
130 {
131 ++retryes;
132 output = pan->Fragment();
133 retry=false;//If here, Fragment did not throw! (AND)
134 }
135 catch( ... )
136 {
137 G4cerr << "***WARNING*** Exception thrown passing through G4ChiralInvariantPhaseSpace "<<G4endl;
138 G4cerr << " targetPDGCode = "<< targetPDGCode <<G4endl;
139 G4cerr << " Dumping the information in the pojectile list"<<G4endl;
140 for(size_t i=0; i< projHV.size(); i++)
141 {
142 G4cerr <<" Incoming 4-momentum and PDG code of "<<i<<"'th hadron: "
143 <<" "<< projHV[i]->Get4Momentum()<<" "<<projHV[i]->GetPDGCode()<<G4endl;
144 }
145 G4cerr << "Retrying interaction "<<G4endl; //AND
146 //throw; //AND
147 }
148 delete pan;
149 } //AND
150 if ( retryes >= maxretries ) //AND
151 {
152 G4cerr << "***ERROR*** Maximum number of retries ("<<maxretries<<") reached for G4QEnvironment::Fragment(), exception is being thrown" << G4endl;
153 throw G4HadronicException(__FILE__,__LINE__,"G4ChiralInvariantPhaseSpace::ApplyYourself(...) - Maximum number of re-tries reached, abandon interaction");
154 }
155
156 // Fill the particle change.
157 G4DynamicParticle * theSec;
158#ifdef CHIPSdebug
159 G4cout << "G4ChiralInvariantPhaseSpace: NEW EVENT #ofHadrons="
160 <<output->size()<<G4endl;
161#endif
162
163
164 unsigned int particle;
165 for( particle = 0; particle < output->size(); particle++)
166 {
167 if(output->operator[](particle)->GetNFragments() != 0)
168 {
169 delete output->operator[](particle);
170 continue;
171 }
172 theSec = new G4DynamicParticle;
173 G4int pdgCode = output->operator[](particle)->GetPDGCode();
174#ifdef CHIPSdebug
175 G4cout << "G4ChiralInvariantPhaseSpace: h#"<<particle
176 <<", PDG="<<pdgCode<<G4endl;
177#endif
178 G4ParticleDefinition * theDefinition;
179 // Note that I still have to take care of strange nuclei
180 // For this I need the mass calculation, and a changed interface
181 // for ion-tablel ==> work for Hisaya @@@@@@@
182 // Then I can sort out the pdgCode. I also need a decau process
183 // for strange nuclei; may be another chips interface
184 if(pdgCode>90000000)
185 {
186 G4int aZ = (pdgCode-90000000)/1000;
187 G4int anN = pdgCode-90000000-1000*aZ;
188 theDefinition = G4ParticleTable::GetParticleTable()->FindIon(aZ,anN+aZ,0,aZ);
189 if(aZ == 0 && anN == 1) theDefinition = G4Neutron::Neutron();
190 }
191 //AND
192 else if ( pdgCode == 90000000 && output->operator[](particle)->Get4Momentum().m()<1*MeV )
193 {
194 //A. Dotti:
195 //We cure the case the model returns a (A,Z)=0,0 G4QHadron with a very small mass
196 //We convert this to a gamma. According to the author of the model this is the
197 //correct thing to do and it is done also in other parts of the CHIPS package
198 theDefinition = G4Gamma::Gamma();
199 }
200 //AND
201 else
202 {
203 theDefinition = G4ParticleTable::GetParticleTable()
204 ->FindParticle(output->operator[](particle)->GetPDGCode());
205 }
206
207 //AND
208 //A. Dotti: Handle problematic cases in which one of the products has not been recognized
209 // This should never happen but we want to add an extra-protection
210 if ( theDefinition == NULL )
211 {
212 //If we arrived here something bad really happened. We do not know how to handle the products of the interaction, we give up, resetting the
213 //result and keeping the primary alive.
214 G4cerr<<"**WARNING*** G4ChiralInvariantPhaseSpace::ApplyYourself(...) : G4QEnvironment::Fragment() returns an invalid fragment\n with fourMom(MeV)=";
215 G4cerr<<output->operator[](particle)->Get4Momentum()<<" and mass(MeV)="<<output->operator[](particle)->Get4Momentum().m();
216 G4cerr<<". Offending PDG is:"<<pdgCode<<" abandon interaction. \n Taget PDG was:"<<targetPDGCode<<" \n Dumping the information in the projectile list:\n";
217 for(size_t i=0; i< projHV.size(); i++)
218 {
219 G4cerr <<" Incoming 4-momentum and PDG code of "<<i<<"'th hadron: "
220 <<" "<< projHV[i]->Get4Momentum()<<" "<<projHV[i]->GetPDGCode()<<"\n";
221 }
222 G4cerr<<"\n Please report as bug \n***END OF MESSAGE***"<<G4endl;
223
224 for ( unsigned int cparticle=0 ; cparticle<output->size();++cparticle)
225 delete output->operator[](cparticle);
226 delete output;
227 std::for_each(projHV.begin(), projHV.end(), DeleteQHadron());
228 projHV.clear();
229
230 aResult->Clear();
231 aResult->SetStatusChange(isAlive);
232 aResult->SetEnergyChange(aTrack.GetKineticEnergy());
233 aResult->SetMomentumChange(aTrack.Get4Momentum().vect().unit());
234 return aResult;
235 }
236 //AND
237
238
239 theSec->SetDefinition(theDefinition);
240 theSec->SetMomentum(output->operator[](particle)->Get4Momentum().vect());
241 aResult->AddSecondary(theSec);
242 delete output->operator[](particle);
243 }
244 delete output;
245 std::for_each(projHV.begin(), projHV.end(), DeleteQHadron());
246 projHV.clear();
247 return aResult;
248}
249
@ isAlive
@ stopAndKill
std::vector< G4QHadron * > G4QHadronVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
G4DLLIMPORT std::ostream G4cout
Hep3Vector unit() const
Hep3Vector vect() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus, G4HadFinalState *aChange=0)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMomentum(const G4ThreeVector &momentum)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetIonMass(G4int Z, G4int A, G4int L=0) const
!! Only ground states are supported now
Definition: G4IonTable.cc:774
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4ParticleDefinition * FindIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4IonTable * GetIonTable()
static G4QCHIPSWorld * Get()
G4QParticleVector * GetParticles(G4int nOfParts=0)
G4QHadronVector * Fragment()
static void SetParameters(G4double fN=.1, G4double fD=.05, G4double cP=4., G4double mR=1., G4double nD=.8 *CLHEP::fermi)
Definition: G4QNucleus.cc:347
static void SetParameters(G4double temper=180., G4double ssin2g=.3, G4double etaetap=.3)
Definition: G4Quasmon.cc:192