Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPEnAngCorrelation Class Reference

#include <G4ParticleHPEnAngCorrelation.hh>

Public Member Functions

 G4ParticleHPEnAngCorrelation ()
 
 G4ParticleHPEnAngCorrelation (G4ParticleDefinition *proj)
 
 ~G4ParticleHPEnAngCorrelation ()
 
void Init (std::istream &aDataFile)
 
G4ReactionProductSampleOne (G4double anEnergy)
 
G4ReactionProductVectorSample (G4double anEnergy)
 
void SetTarget (G4ReactionProduct &aTarget)
 
void SetProjectileRP (G4ReactionProduct &aIncidentPart)
 
G4bool InCharge ()
 
G4double GetTargetMass ()
 
G4double GetTotalMeanEnergy ()
 

Detailed Description

Definition at line 46 of file G4ParticleHPEnAngCorrelation.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPEnAngCorrelation() [1/2]

G4ParticleHPEnAngCorrelation::G4ParticleHPEnAngCorrelation ( )
inline

Definition at line 59 of file G4ParticleHPEnAngCorrelation.hh.

60 {
61 theProjectile = G4Neutron::Neutron();
62 theProducts = 0;
63 inCharge = false;
64 toBeCached val;
65 fCache.Put( val );
66 //theTotalMeanEnergy = -1.;
67 fCache.Get().theTotalMeanEnergy = -1.;
68 targetMass = 0.0;
69 frameFlag = 0;
70 nProducts = 0;
71 }
value_type & Get() const
Definition: G4Cache.hh:315
void Put(const value_type &val) const
Definition: G4Cache.hh:321
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103

◆ G4ParticleHPEnAngCorrelation() [2/2]

G4ParticleHPEnAngCorrelation::G4ParticleHPEnAngCorrelation ( G4ParticleDefinition proj)
inline

Definition at line 73 of file G4ParticleHPEnAngCorrelation.hh.

74 : theProjectile(proj)
75 {
76 theProducts = 0;
77 inCharge = false;
78 toBeCached val;
79 fCache.Put( val );
80 //theTotalMeanEnergy = -1.;
81 fCache.Get().theTotalMeanEnergy = -1.;
82 targetMass = 0.0;
83 frameFlag = 0;
84 nProducts = 0;
85 }

◆ ~G4ParticleHPEnAngCorrelation()

G4ParticleHPEnAngCorrelation::~G4ParticleHPEnAngCorrelation ( )
inline

Definition at line 87 of file G4ParticleHPEnAngCorrelation.hh.

88 {
89 if(theProducts!=0) delete [] theProducts;
90 }

Member Function Documentation

◆ GetTargetMass()

G4double G4ParticleHPEnAngCorrelation::GetTargetMass ( )
inline

Definition at line 126 of file G4ParticleHPEnAngCorrelation.hh.

127 {
128 return targetMass;
129 }

Referenced by G4ParticleHPInelasticBaseFS::BaseApply().

◆ GetTotalMeanEnergy()

G4double G4ParticleHPEnAngCorrelation::GetTotalMeanEnergy ( )
inline

Definition at line 131 of file G4ParticleHPEnAngCorrelation.hh.

132 {
133 return fCache.Get().theTotalMeanEnergy; // cached in 'sample' call
134 }

Referenced by G4ParticleHPInelasticBaseFS::BaseApply().

◆ InCharge()

G4bool G4ParticleHPEnAngCorrelation::InCharge ( )
inline

Definition at line 121 of file G4ParticleHPEnAngCorrelation.hh.

122 {
123 return inCharge;
124 }

◆ Init()

void G4ParticleHPEnAngCorrelation::Init ( std::istream &  aDataFile)
inline

Definition at line 92 of file G4ParticleHPEnAngCorrelation.hh.

93 {
94 inCharge = true;
95 aDataFile>>targetMass>>frameFlag>>nProducts;
96 theProducts = new G4ParticleHPProduct[nProducts];
97 for(G4int i=0; i<nProducts; i++)
98 {
99 theProducts[i].Init(aDataFile,theProjectile);
100 }
101 }
int G4int
Definition: G4Types.hh:85
void Init(std::istream &aDataFile, G4ParticleDefinition *projectile)

Referenced by G4ParticleHPCaptureFS::Init(), G4ParticleHPInelasticCompFS::Init(), and G4ParticleHPInelasticBaseFS::Init().

◆ Sample()

G4ReactionProductVector * G4ParticleHPEnAngCorrelation::Sample ( G4double  anEnergy)

Definition at line 77 of file G4ParticleHPEnAngCorrelation.cc.

78{
80 G4int i;
82 G4ReactionProduct theCMS;
84
85 if(frameFlag==2
86 || frameFlag==3) // Added for particle HP
87 {
88 // simplify and double check @
89 G4ThreeVector the3IncidentPart = fCache.Get().theProjectileRP->GetMomentum(); //theProjectileRP has value in LAB
90 G4double nEnergy = fCache.Get().theProjectileRP->GetTotalEnergy();
91 G4ThreeVector the3Target = fCache.Get().theTarget->GetMomentum(); //theTarget has value in LAB
92 G4double tEnergy = fCache.Get().theTarget->GetTotalEnergy();
93 G4double totE = nEnergy+tEnergy;
94 G4ThreeVector the3CMS = the3Target+the3IncidentPart;
95 theCMS.SetMomentum(the3CMS);
96 G4double cmsMom = std::sqrt(the3CMS*the3CMS);
97 G4double sqrts = std::sqrt((totE-cmsMom)*(totE+cmsMom));
98 theCMS.SetMass(sqrts);
99 theCMS.SetTotalEnergy(totE);
100 G4ReactionProduct aIncidentPart;
101 aIncidentPart.Lorentz(*fCache.Get().theProjectileRP, theCMS);
102 //TKDB 100413
103 //ENDF-6 Formats Manual ENDF-102
104 //CHAPTER 6. FILE 6: PRODUCT ENERGY-ANGLE DISTRIBUTIONS
105 //LCT Reference system for secondary energy and angle (incident energy is always given in the LAB system)
106 //anEnergy = aIncidentPart.GetKineticEnergy();
107 anEnergy = fCache.Get().theProjectileRP->GetKineticEnergy(); //should be same argumment of "anEnergy"
108
109 G4LorentzVector Ptmp (aIncidentPart.GetMomentum(), aIncidentPart.GetTotalEnergy());
110
111 toZ.rotateZ(-1*Ptmp.phi());
112 toZ.rotateY(-1*Ptmp.theta());
113 }
114 fCache.Get().theTotalMeanEnergy=0;
115 G4LorentzRotation toLab(toZ.inverse()); //toLab only change axis NOT to LAB system
116
117 for(i=0; i<nProducts; i++)
118 {
119 G4int nPart = theProducts[i].GetMultiplicity(anEnergy);
120 //- if( nParticles[i] == 0 ) continue;
121 it = theProducts[i].Sample(anEnergy,nPart);
122 G4double aMeanEnergy = theProducts[i].MeanEnergyOfThisInteraction();
123 //if(aMeanEnergy>0)
124 //151120 TK Modified for solving reproducibility problem
125 //This change may have side effect.
126 if(aMeanEnergy>=0)
127 {
128 fCache.Get().theTotalMeanEnergy += aMeanEnergy;
129 }
130 else
131 {
132 fCache.Get().theTotalMeanEnergy = anEnergy/nProducts+theProducts[i].GetQValue();
133 }
134 if(it!=0)
135 {
136 for(unsigned int ii=0; ii<it->size(); ii++)
137 {
138 //if(!std::getenv("G4PHP_NO_LORENTZ_BOOST")) {
139 G4LorentzVector pTmp1 (it->operator[](ii)->GetMomentum(),
140 it->operator[](ii)->GetTotalEnergy());
141 pTmp1 = toLab*pTmp1;
142 it->operator[](ii)->SetMomentum(pTmp1.vect());
143 it->operator[](ii)->SetTotalEnergy(pTmp1.e());
144
145 if(frameFlag==1) // target rest //TK 100413 should be LAB?
146 {
147 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*(*fCache.Get().theTarget)); //TK 100413 Is this really need?
148 }
149 else if(frameFlag==2 ) // CMS
150 {
151#ifdef G4PHPDEBUG
152 if( std::getenv("G4ParticleHPDebug") )
153 G4cout <<"G4ParticleHPEnAngCorrelation: before Lorentz boost "<<
154 it->at(ii)->GetKineticEnergy()<<" "<<
155 it->at(ii)->GetMomentum()<<G4endl;
156#endif
157 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theCMS);
158#ifdef G4PHPDEBUG
159 if( std::getenv("G4ParticleHPDebug") )
160 G4cout <<"G4ParticleHPEnAngCorrelation: after Lorentz boost "<<
161 it->at(ii)->GetKineticEnergy()<<" "<<
162 it->at(ii)->GetMomentum()<<G4endl;
163#endif
164 }
165 //TK120515 migrate frameFlag (MF6 LCT) = 3
166 else if(frameFlag==3) // CMS A<=4 other LAB
167 {
168 if ( theProducts[i].GetMassCode() > 2004.5 ) //Alpha AWP 3.96713
169 {
170 //LAB
171 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*(*fCache.Get().theTarget)); //TK 100413 Is this really need?
172#ifdef G4PHPDEBUG
173 if( std::getenv("G4ParticleHPDebug") )
174 G4cout <<"G4ParticleHPEnAngCorrelation: after Lorentz boost "<<
175 it->at(ii)->GetKineticEnergy()<<" "<<
176 it->at(ii)->GetMomentum()<<G4endl;
177#endif
178 }
179 else
180 {
181 //CMS
182 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theCMS);
183#ifdef G4PHPDEBUG
184 if( std::getenv("G4ParticleHPDebug") )
185 G4cout <<"G4ParticleHPEnAngCorrelation: after Lorentz boost "<<
186 it->at(ii)->GetKineticEnergy()<<" "<<
187 it->at(ii)->GetMomentum()<<G4endl;
188#endif
189 }
190 }
191 else
192 {
193 throw G4HadronicException(__FILE__, __LINE__, "G4ParticleHPEnAngCorrelation::Sample: The frame of the finalstate is not specified");
194 }
195
196 // }//getenv("G4PHP_NO_LORENTZ_BOOST"))
197 // G4cout << ii << " EnAnG energy after boost " << it->operator[](ii)->GetKineticEnergy() << G4endl; //GDEB
198 result->push_back(it->operator[](ii));
199 }
200 delete it;
201 }
202 }
203
204 return result;
205}
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
HepLorentzRotation & rotateY(double delta)
HepLorentzRotation & rotateZ(double delta)
HepLorentzRotation inverse() const
G4ReactionProductVector * Sample(G4double anEnergy, G4int nParticles)
G4double MeanEnergyOfThisInteraction()
G4int GetMultiplicity(G4double anEnergy)
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
void SetMass(const G4double mas)

Referenced by G4ParticleHPCaptureFS::ApplyYourself(), G4ParticleHPInelasticBaseFS::BaseApply(), and G4ParticleHPInelasticCompFS::CompositeApply().

◆ SampleOne()

G4ReactionProduct * G4ParticleHPEnAngCorrelation::SampleOne ( G4double  anEnergy)

Definition at line 42 of file G4ParticleHPEnAngCorrelation.cc.

43{
45
46 // do we have an appropriate distribution
47 if(nProducts!=1) throw G4HadronicException(__FILE__, __LINE__, "More than one product in SampleOne");
48
49 // get the result
51 G4int i=0;
52
53 G4int icounter=0;
54 G4int icounter_max=1024;
55 while(temp == 0) {
56 icounter++;
57 if ( icounter > icounter_max ) {
58 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
59 break;
60 }
61 temp = theProducts[i++].Sample(anEnergy,1);
62 }
63
64 // is the multiplicity correct
65 if(temp->size()!=1) throw G4HadronicException(__FILE__, __LINE__, "SampleOne: Yield not correct");
66
67 // fill result
68 result = temp->operator[](0);
69
70 // some garbage collection
71 delete temp;
72
73 // return result
74 return result;
75}

◆ SetProjectileRP()

void G4ParticleHPEnAngCorrelation::SetProjectileRP ( G4ReactionProduct aIncidentPart)
inline

Definition at line 114 of file G4ParticleHPEnAngCorrelation.hh.

115 {
116 fCache.Get().theProjectileRP = &aIncidentPart;
117 for (G4int i=0;i<nProducts;i++)
118 theProducts[i].SetProjectileRP(fCache.Get().theProjectileRP);
119 }
void SetProjectileRP(G4ReactionProduct &aIncidentPart)

Referenced by G4ParticleHPCaptureFS::ApplyYourself(), G4ParticleHPInelasticBaseFS::BaseApply(), G4ParticleHPInelasticCompFS::InitDistributionInitialState(), and SetProjectileRP().

◆ SetTarget()

void G4ParticleHPEnAngCorrelation::SetTarget ( G4ReactionProduct aTarget)
inline

Definition at line 107 of file G4ParticleHPEnAngCorrelation.hh.

108 {
109 fCache.Get().theTarget = &aTarget;
110 for (G4int i=0;i<nProducts;i++)
111 theProducts[i].SetTarget(fCache.Get().theTarget);
112 }
void SetTarget(G4ReactionProduct &aTarget)

Referenced by G4ParticleHPCaptureFS::ApplyYourself(), G4ParticleHPInelasticBaseFS::BaseApply(), G4ParticleHPInelasticCompFS::InitDistributionInitialState(), and SetTarget().


The documentation for this class was generated from the following files: