Geant4 10.7.0
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 bAdjustFinalState = true;
72 }
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 74 of file G4ParticleHPEnAngCorrelation.hh.

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

◆ ~G4ParticleHPEnAngCorrelation()

G4ParticleHPEnAngCorrelation::~G4ParticleHPEnAngCorrelation ( )
inline

Definition at line 89 of file G4ParticleHPEnAngCorrelation.hh.

90 {
91 if(theProducts!=0) delete [] theProducts;
92 }

Member Function Documentation

◆ GetTargetMass()

G4double G4ParticleHPEnAngCorrelation::GetTargetMass ( )
inline

Definition at line 137 of file G4ParticleHPEnAngCorrelation.hh.

138 {
139 return targetMass;
140 }

Referenced by G4ParticleHPInelasticBaseFS::BaseApply().

◆ GetTotalMeanEnergy()

G4double G4ParticleHPEnAngCorrelation::GetTotalMeanEnergy ( )
inline

Definition at line 142 of file G4ParticleHPEnAngCorrelation.hh.

143 {
144 return fCache.Get().theTotalMeanEnergy; // cached in 'sample' call
145 }

Referenced by G4ParticleHPInelasticBaseFS::BaseApply().

◆ InCharge()

G4bool G4ParticleHPEnAngCorrelation::InCharge ( )
inline

Definition at line 132 of file G4ParticleHPEnAngCorrelation.hh.

133 {
134 return inCharge;
135 }

◆ Init()

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

Definition at line 94 of file G4ParticleHPEnAngCorrelation.hh.

95 {
96 bAdjustFinalState = true;
97 if ( G4ParticleHPManager::GetInstance()->GetDoNotAdjustFinalState() ) bAdjustFinalState = false;
98
99// T.K. Comment out following line to keep the condition at the
100// validation efforts comparing NeutronHP and PartileHP for neutrons (2015 Sep.)
101//#ifdef PHP_AS_HP
102// bAdjustFinalState = false;
103//#endif
104
105 inCharge = true;
106 aDataFile>>targetMass>>frameFlag>>nProducts;
107 theProducts = new G4ParticleHPProduct[nProducts];
108 for(G4int i=0; i<nProducts; i++)
109 {
110 theProducts[i].Init(aDataFile,theProjectile);
111 }
112 }
int G4int
Definition: G4Types.hh:85
static G4ParticleHPManager * GetInstance()
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( getenv("G4PHPTEST") ) G4cout << " EnAnG energy sampled " << it->operator[](0)->GetKineticEnergy() << " aMeanEnergy " << aMeanEnergy << G4endl; // GDEB
124 //if(aMeanEnergy>0)
125 //151120 TK Modified for solving reproducibility problem
126 //This change may have side effect.
127 if(aMeanEnergy>=0)
128 {
129 fCache.Get().theTotalMeanEnergy += aMeanEnergy;
130 }
131 else
132 {
133 fCache.Get().theTotalMeanEnergy = anEnergy/nProducts+theProducts[i].GetQValue();
134 }
135 if(it!=0)
136 {
137 for(unsigned int ii=0; ii<it->size(); ii++)
138 {
139 //if(!std::getenv("G4PHP_NO_LORENTZ_BOOST")) {
140 G4LorentzVector pTmp1 (it->operator[](ii)->GetMomentum(),
141 it->operator[](ii)->GetTotalEnergy());
142 pTmp1 = toLab*pTmp1;
143 if( std::getenv("G4PHPTEST") ) G4cout << " G4particleHPEnAngCorrelation COS THETA " << std::cos(it->operator[](ii)->GetMomentum().theta()) << G4endl;
144 it->operator[](ii)->SetMomentum(pTmp1.vect());
145 it->operator[](ii)->SetTotalEnergy(pTmp1.e());
146 if( std::getenv("G4PHPTEST") ) G4cout << " G4particleHPEnAngCorrelation COS THETA after toLab " << std::cos(it->operator[](ii)->GetMomentum().theta()) << G4endl;
147
148 if(frameFlag==1) // target rest //TK 100413 should be LAB?
149 {
150 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*(*fCache.Get().theTarget)); //TK 100413 Is this really need?
151 }
152 else if(frameFlag==2 ) // CMS
153 {
154#ifdef G4PHPDEBUG
155 if( std::getenv("G4ParticleHPDebug") )
156 G4cout <<"G4ParticleHPEnAngCorrelation: before Lorentz boost "<<
157 it->at(ii)->GetKineticEnergy()<<" "<<
158 it->at(ii)->GetMomentum()<<G4endl;
159#endif
160 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theCMS);
161#ifdef G4PHPDEBUG
162 if( std::getenv("G4ParticleHPDebug") )
163 G4cout <<"G4ParticleHPEnAngCorrelation: after Lorentz boost "<<
164 it->at(ii)->GetKineticEnergy()<<" "<<
165 it->at(ii)->GetMomentum()<<G4endl;
166#endif
167 }
168 //TK120515 migrate frameFlag (MF6 LCT) = 3
169 else if(frameFlag==3) // CMS A<=4 other LAB
170 {
171 if ( theProducts[i].GetMassCode() > 2004.5 ) //Alpha AWP 3.96713
172 {
173 //LAB
174 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*(*fCache.Get().theTarget)); //TK 100413 Is this really need?
175#ifdef G4PHPDEBUG
176 if( std::getenv("G4ParticleHPDebug") )
177 G4cout <<"G4ParticleHPEnAngCorrelation: after Lorentz boost "<<
178 it->at(ii)->GetKineticEnergy()<<" "<<
179 it->at(ii)->GetMomentum()<<G4endl;
180#endif
181 }
182 else
183 {
184 //CMS
185 it->operator[](ii)->Lorentz(*(it->operator[](ii)), -1.*theCMS);
186#ifdef G4PHPDEBUG
187 if( std::getenv("G4ParticleHPDebug") )
188 G4cout <<"G4ParticleHPEnAngCorrelation: after Lorentz boost "<<
189 it->at(ii)->GetKineticEnergy()<<" "<<
190 it->at(ii)->GetMomentum()<<G4endl;
191#endif
192 }
193 }
194 else
195 {
196 throw G4HadronicException(__FILE__, __LINE__, "G4ParticleHPEnAngCorrelation::Sample: The frame of the finalstate is not specified");
197 }
198 if( std::getenv("G4PHPTEST") ) G4cout << frameFlag << " G4particleHPEnAngCorrelation COS THETA after Lorentz " << std::cos(it->operator[](ii)->GetMomentum().theta()) << G4endl;
199
200 // }//getenv("G4PHP_NO_LORENTZ_BOOST"))
201 // G4cout << ii << " EnAnG energy after boost " << it->operator[](ii)->GetKineticEnergy() << G4endl; //GDEB
202 result->push_back(it->operator[](ii));
203 }
204 delete it;
205 }
206 }
207
208 return result;
209}
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 125 of file G4ParticleHPEnAngCorrelation.hh.

126 {
127 fCache.Get().theProjectileRP = &aIncidentPart;
128 for (G4int i=0;i<nProducts;i++)
129 theProducts[i].SetProjectileRP(fCache.Get().theProjectileRP);
130 }
void SetProjectileRP(G4ReactionProduct &aIncidentPart)

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

◆ SetTarget()

void G4ParticleHPEnAngCorrelation::SetTarget ( G4ReactionProduct aTarget)
inline

Definition at line 118 of file G4ParticleHPEnAngCorrelation.hh.

119 {
120 fCache.Get().theTarget = &aTarget;
121 for (G4int i=0;i<nProducts;i++)
122 theProducts[i].SetTarget(fCache.Get().theTarget);
123 }
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: