Geant4 11.2.2
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 (const G4ParticleDefinition *proj=nullptr)
 
 ~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 GetNumberOfProducts ()
 
G4double GetTotalMeanEnergy ()
 
 G4ParticleHPEnAngCorrelation (G4ParticleHPEnAngCorrelation &)=delete
 
G4ParticleHPEnAngCorrelationoperator= (const G4ParticleHPEnAngCorrelation &right)=delete
 

Detailed Description

Definition at line 46 of file G4ParticleHPEnAngCorrelation.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPEnAngCorrelation() [1/2]

G4ParticleHPEnAngCorrelation::G4ParticleHPEnAngCorrelation ( const G4ParticleDefinition * proj = nullptr)
explicit

Definition at line 50 of file G4ParticleHPEnAngCorrelation.cc.

51{
52 theProjectile = (nullptr == proj) ? G4Neutron::Neutron() : proj;
53 toBeCached val;
54 fCache.Put( val );
55}
void Put(const value_type &val) const
Definition G4Cache.hh:321
static G4Neutron * Neutron()
Definition G4Neutron.cc:101

◆ ~G4ParticleHPEnAngCorrelation()

G4ParticleHPEnAngCorrelation::~G4ParticleHPEnAngCorrelation ( )

Definition at line 57 of file G4ParticleHPEnAngCorrelation.cc.

58{
59 delete[] theProducts;
60}

◆ G4ParticleHPEnAngCorrelation() [2/2]

G4ParticleHPEnAngCorrelation::G4ParticleHPEnAngCorrelation ( G4ParticleHPEnAngCorrelation & )
delete

Member Function Documentation

◆ GetNumberOfProducts()

G4double G4ParticleHPEnAngCorrelation::GetNumberOfProducts ( )
inline

Definition at line 84 of file G4ParticleHPEnAngCorrelation.hh.

84{ return nProducts; }

◆ GetTargetMass()

G4double G4ParticleHPEnAngCorrelation::GetTargetMass ( )
inline

Definition at line 82 of file G4ParticleHPEnAngCorrelation.hh.

82{ return targetMass; }

Referenced by G4ParticleHPInelasticBaseFS::BaseApply().

◆ GetTotalMeanEnergy()

G4double G4ParticleHPEnAngCorrelation::GetTotalMeanEnergy ( )
inline

Definition at line 86 of file G4ParticleHPEnAngCorrelation.hh.

87 {
88 return fCache.Get().theTotalMeanEnergy; // cached in 'sample' call
89 }
value_type & Get() const
Definition G4Cache.hh:315

Referenced by G4ParticleHPInelasticBaseFS::BaseApply().

◆ InCharge()

G4bool G4ParticleHPEnAngCorrelation::InCharge ( )
inline

Definition at line 80 of file G4ParticleHPEnAngCorrelation.hh.

80{ return inCharge; }

◆ Init()

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

Definition at line 62 of file G4ParticleHPEnAngCorrelation.cc.

63{
64 inCharge = true;
65 aDataFile >> targetMass >> frameFlag >> nProducts;
66 const std::size_t psize = nProducts > 0 ? nProducts : 1;
67 theProducts = new G4ParticleHPProduct[psize];
68 for (std::size_t i = 0; i < psize; ++i) {
69 theProducts[i].Init(aDataFile, theProjectile);
70 }
71}
void Init(std::istream &aDataFile, const G4ParticleDefinition *projectile)

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

◆ operator=()

G4ParticleHPEnAngCorrelation & G4ParticleHPEnAngCorrelation::operator= ( const G4ParticleHPEnAngCorrelation & right)
delete

◆ Sample()

G4ReactionProductVector * G4ParticleHPEnAngCorrelation::Sample ( G4double anEnergy)

Definition at line 116 of file G4ParticleHPEnAngCorrelation.cc.

117{
118 auto result = new G4ReactionProductVector;
119 G4int i;
121 G4ReactionProduct theCMS;
123 /*
124 G4cout << "G4ParticleHPEnAngCorrelation::Sample for E=" << anEnergy
125 << " flag=" << frameFlag << " nProducts=" << nProducts <<G4endl;
126 */
127 if (frameFlag == 2 || frameFlag == 3) // Added for particle HP
128 {
129 G4ThreeVector the3IncidentPart = fCache.Get().theProjectileRP->GetMomentum();
130 // theProjectileRP has value in LAB
131 G4double nEnergy = fCache.Get().theProjectileRP->GetTotalEnergy();
132 G4ThreeVector the3Target = fCache.Get().theTarget->GetMomentum();
133 // theTarget has value in LAB
134 G4double tEnergy = fCache.Get().theTarget->GetTotalEnergy();
135 G4double totE = nEnergy + tEnergy;
136 G4ThreeVector the3CMS = the3Target + the3IncidentPart;
137 theCMS.SetMomentum(the3CMS);
138 G4double cmsMom = std::sqrt(the3CMS * the3CMS);
139 G4double sqrts = std::sqrt((totE - cmsMom) * (totE + cmsMom));
140 theCMS.SetMass(sqrts);
141 theCMS.SetTotalEnergy(totE);
142 G4ReactionProduct aIncidentPart;
143 aIncidentPart.Lorentz(*fCache.Get().theProjectileRP, theCMS);
144 // TKDB 100413
145 // ENDF-6 Formats Manual ENDF-102
146 // CHAPTER 6. FILE 6: PRODUCT ENERGY-ANGLE DISTRIBUTIONS
147 // LCT Reference system for secondary energy and angle
148 // incident energy is always given in the LAB system
149 anEnergy = fCache.Get().theProjectileRP->GetKineticEnergy();
150 // should be same argumment of "anEnergy"
151
152 G4LorentzVector Ptmp(aIncidentPart.GetMomentum(),
153 aIncidentPart.GetTotalEnergy());
154
155 toZ.rotateZ(-1 * Ptmp.phi());
156 toZ.rotateY(-1 * Ptmp.theta());
157 }
158 fCache.Get().theTotalMeanEnergy = 0;
159 G4LorentzRotation toLab(toZ.inverse());
160 // toLab only change axis NOT to LAB system
161
162 for (i = 0; i < nProducts; ++i) {
163 G4int nPart = theProducts[i].GetMultiplicity(anEnergy);
164 it = theProducts[i].Sample(anEnergy, nPart);
165 G4double aMeanEnergy = theProducts[i].MeanEnergyOfThisInteraction();
166
167 // 151120 TK Modified for solving reproducibility problem
168 // This change may have side effect.
169 if (aMeanEnergy >= 0) {
170 fCache.Get().theTotalMeanEnergy += aMeanEnergy;
171 }
172 else {
173 fCache.Get().theTotalMeanEnergy = anEnergy/nProducts + theProducts[i].GetQValue();
174 }
175 if (it != nullptr) {
176 for (unsigned int ii = 0; ii < it->size(); ++ii) {
177 G4LorentzVector pTmp1(it->operator[](ii)->GetMomentum(),
178 it->operator[](ii)->GetTotalEnergy());
179 pTmp1 = toLab * pTmp1;
180 it->operator[](ii)->SetMomentum(pTmp1.vect());
181 it->operator[](ii)->SetTotalEnergy(pTmp1.e());
182
183 if (frameFlag == 1) // target rest //TK 100413 should be LAB?
184 {
185 it->operator[](ii)->Lorentz(
186 *(it->operator[](ii)),
187 -1. * (*fCache.Get().theTarget)); // TK 100413 Is this really need?
188 }
189 else if (frameFlag == 2) // CMS
190 {
191#ifdef G4PHPDEBUG
192 if (G4ParticleHPManager::GetInstance()->GetDEBUG())
193 G4cout << "G4ParticleHPEnAngCorrelation: before Lorentz boost "
194 << it->at(ii)->GetKineticEnergy() << " "
195 << it->at(ii)->GetMomentum() << G4endl;
196#endif
197 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1. * theCMS);
198#ifdef G4PHPDEBUG
199 if (G4ParticleHPManager::GetInstance()->GetDEBUG())
200 G4cout << "G4ParticleHPEnAngCorrelation: after Lorentz boost "
201 << it->at(ii)->GetKineticEnergy() << " "
202 << it->at(ii)->GetMomentum() << G4endl;
203#endif
204 }
205 // TK120515 migrate frameFlag (MF6 LCT) = 3
206 else if (frameFlag == 3) // CMS A<=4 other LAB
207 {
208 if (theProducts[i].GetMassCode() > 2004.5) // Alpha AWP 3.96713
209 {
210 // LAB
211 it->operator[](ii)->Lorentz(
212 *(it->operator[](ii)),
213 -1. * (*fCache.Get().theTarget)); // TK 100413 Is this really need?
214#ifdef G4PHPDEBUG
215 if (G4ParticleHPManager::GetInstance()->GetDEBUG())
216 G4cout << "G4ParticleHPEnAngCorrelation: after Lorentz boost "
217 << it->at(ii)->GetKineticEnergy() << " "
218 << it->at(ii)->GetMomentum() << G4endl;
219#endif
220 }
221 else {
222 // CMS
223 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1. * theCMS);
224#ifdef G4PHPDEBUG
225 if (G4ParticleHPManager::GetInstance()->GetDEBUG())
226 G4cout << "G4ParticleHPEnAngCorrelation: after Lorentz boost "
227 << it->at(ii)->GetKineticEnergy() << " "
228 << it->at(ii)->GetMomentum() << G4endl;
229#endif
230 }
231 }
232 else
233 {
234 throw G4HadronicException(__FILE__, __LINE__,
235 "Sample: The frame of the final state is not specified");
236 }
237 result->push_back(it->operator[](ii));
238 }
239 delete it;
240 }
241 }
242 return result;
243}
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
HepLorentzRotation & rotateY(double delta)
HepLorentzRotation & rotateZ(double delta)
HepLorentzRotation inverse() const
static G4ParticleHPManager * GetInstance()
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 G4NeutronHPCaptureFS::ApplyYourself(), G4ParticleHPInelasticBaseFS::BaseApply(), and G4ParticleHPInelasticCompFS::CompositeApply().

◆ SampleOne()

G4ReactionProduct * G4ParticleHPEnAngCorrelation::SampleOne ( G4double anEnergy)

Definition at line 73 of file G4ParticleHPEnAngCorrelation.cc.

74{
75 auto result = new G4ReactionProduct;
76
77 // do we have an appropriate distribution
78 if (nProducts != 1)
79 throw G4HadronicException(__FILE__, __LINE__, "More than one product in SampleOne");
80
81 // get the result
82 G4ReactionProductVector* temp = nullptr;
83 G4int i = 0;
84
85 G4int icounter = 0;
86 G4int icounter_max = 1024;
87 while (temp == nullptr)
88 {
89 ++icounter;
90 if (icounter > icounter_max)
91 {
92 G4cout << "Loop-counter exceeded the threshold value at "
93 << __LINE__ << "th line of "
94 << __FILE__ << "." << G4endl;
95 break;
96 }
97 temp = theProducts[++i].Sample(anEnergy, 1);
98 }
99 if (nullptr != temp)
100 {
101 if (temp->size() == 1)
102 {
103 result = (*temp)[0];
104 }
105 else
106 {
107 for (auto const& ptr : *temp) { delete ptr; }
108 throw G4HadronicException(__FILE__, __LINE__,
109 "SampleOne: Yield not correct");
110 }
111 }
112 delete temp;
113 return result;
114}

◆ SetProjectileRP()

void G4ParticleHPEnAngCorrelation::SetProjectileRP ( G4ReactionProduct & aIncidentPart)
inline

Definition at line 73 of file G4ParticleHPEnAngCorrelation.hh.

74 {
75 fCache.Get().theProjectileRP = &aIncidentPart;
76 for (G4int i = 0; i < nProducts; ++i)
77 theProducts[i].SetProjectileRP(&aIncidentPart);
78 }
void SetProjectileRP(G4ReactionProduct &aIncidentPart)

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

◆ SetTarget()

void G4ParticleHPEnAngCorrelation::SetTarget ( G4ReactionProduct & aTarget)
inline

Definition at line 66 of file G4ParticleHPEnAngCorrelation.hh.

67 {
68 fCache.Get().theTarget = &aTarget;
69 for (G4int i = 0; i < nProducts; ++i)
70 theProducts[i].SetTarget(&aTarget);
71 }
void SetTarget(G4ReactionProduct &aTarget)

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


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