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

#include <G4VLEPTSModel.hh>

+ Inheritance diagram for G4VLEPTSModel:

Public Member Functions

 G4VLEPTSModel (const G4String &processName)
 
 ~G4VLEPTSModel () override
 
void BuildPhysicsTable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Material *mate, const G4ParticleDefinition *aParticle, G4double kineticEnergy)
 
G4ThreeVector SampleNewDirection (const G4Material *aMaterial, G4ThreeVector Dir, G4double e, G4double el)
 
G4double SampleAngle (const G4Material *aMaterial, G4double e, G4double el)
 
G4ThreeVector SampleNewDirection (G4ThreeVector Dir, G4double ang)
 
G4VLEPTSModeloperator= (const G4VLEPTSModel &right)
 
 G4VLEPTSModel (const G4VLEPTSModel &)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Member Functions

void Init ()
 
G4bool ReadParam (G4String fileName, const G4Material *aMaterial)
 
virtual std::map< G4int, std::vector< G4double > > ReadIXS (G4String fileName, const G4Material *aMaterial)
 
G4double SampleEnergyLoss (const G4Material *aMaterial, G4double eMin, G4double eMax)
 
void BuildMeanFreePathTable (const G4Material *aMaterial, std::map< G4int, std::vector< G4double > > &integralXS)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4PhysicsTabletheMeanFreePathTable
 
G4double theLowestEnergyLimit
 
G4double theHighestEnergyLimit
 
G4int theNumbBinTable
 
std::map< const G4Material *, G4doubletheIonisPot
 
std::map< const G4Material *, G4doubletheIonisPotInt
 
std::map< const G4Material *, G4doubletheMolecularMass
 
std::map< const G4Material *, G4LEPTSDiffXS * > theDiffXS
 
std::map< const G4Material *, G4LEPTSDistribution * > theRMTDistr
 
std::map< const G4Material *, G4LEPTSElossDistr * > theElostDistr
 
std::map< const G4Material *, G4LEPTSDistribution * > theElostDistr2
 
std::map< const G4Material *, G4inttheNXSdat
 
std::map< const G4Material *, G4inttheNXSsub
 
G4bool isInitialised {false}
 
XSType theXSType
 
G4int verboseLevel
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 59 of file G4VLEPTSModel.hh.

Constructor & Destructor Documentation

◆ G4VLEPTSModel() [1/2]

G4VLEPTSModel::G4VLEPTSModel ( const G4String & processName)

Definition at line 31 of file G4VLEPTSModel.cc.

31 : G4VEmModel(modelName)
32{
34
36
37 verboseLevel = 0;
38
39 theLowestEnergyLimit = 0.5*CLHEP::eV;
40
41 theHighestEnergyLimit = 1.0*CLHEP::MeV;
42
44}
@ XSEnergy
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
G4PhysicsTable * theMeanFreePathTable
G4double theHighestEnergyLimit
G4double theLowestEnergyLimit

◆ ~G4VLEPTSModel()

G4VLEPTSModel::~G4VLEPTSModel ( )
override

Definition at line 48 of file G4VLEPTSModel.cc.

49{
50
51 if(theMeanFreePathTable != nullptr) {
54 }
55}
void clearAndDestroy()

◆ G4VLEPTSModel() [2/2]

G4VLEPTSModel::G4VLEPTSModel ( const G4VLEPTSModel & )

Member Function Documentation

◆ BuildMeanFreePathTable()

void G4VLEPTSModel::BuildMeanFreePathTable ( const G4Material * aMaterial,
std::map< G4int, std::vector< G4double > > & integralXS )
protected

Definition at line 181 of file G4VLEPTSModel.cc.

182{
183 G4double LowEdgeEnergy, fValue;
184
185 //BUILD MEAN FREE PATH TABLE FROM INTEGRAL CROSS SECTION
186 std::size_t matIdx = aMaterial->GetIndex();
188
189 for (G4int ii=0; ii < theNumbBinTable; ++ii) {
190 LowEdgeEnergy = ptrVector->Energy(ii);
191 if( verboseLevel >= 2 ) G4cout << GetName() << " " << ii << " Energy " << LowEdgeEnergy << " > " << theLowestEnergyLimit << " < " << theHighestEnergyLimit << G4endl;
192 //- fValue = ComputeMFP(LowEdgeEnergy, material, aParticleName);
193 fValue = 0.;
194 if( LowEdgeEnergy >= theLowestEnergyLimit &&
195 LowEdgeEnergy <= theHighestEnergyLimit) {
196 G4double NbOfMoleculesPerVolume = aMaterial->GetDensity()/theMolecularMass[aMaterial]*CLHEP::Avogadro;
197
198 G4double SIGMA = 0. ;
199 //- for ( std::size_t elm=0 ; elm < aMaterial->GetNumberOfElements() ; elm++ ) {
200 G4double crossSection = 0.;
201
202 G4double eVEnergy = LowEdgeEnergy/CLHEP::eV;
203
204 //- if( verboseLevel >= 2 ) G4cout << " eVEnergy " << eVEnergy << " LowEdgeE " << LowEdgeEnergy << " " << integralXS[theXSType][1] << G4endl;
205
206 if(eVEnergy < integralXS[0][1] ) {
207 crossSection = 0.;
208 } else {
209 G4int Bin = 0; // locate bin
210 G4double aa, bb;
211 for( G4int jj=1; jj<theNXSdat[aMaterial]; jj++) { // Extrapolate for E > Emax !!!
212 if( verboseLevel >= 3 ) G4cout << " GET BIN " << jj << " "<< eVEnergy << " > " << integralXS[0][jj] << G4endl;
213 if( eVEnergy > integralXS[0][jj]) {
214 Bin = jj;
215 } else {
216 break;
217 }
218 }
219 aa = integralXS[0][Bin];
220 bb = integralXS[0][Bin+1];
221 crossSection = (integralXS[theXSType][Bin] + (integralXS[theXSType][Bin+1]-integralXS[theXSType][Bin])/(bb-aa)*(eVEnergy-aa) ) * 1.e-16*CLHEP::cm2;
222
223 if( verboseLevel >= 3 ) G4cout << " crossSection " << crossSection << " " <<integralXS[theXSType][Bin] << " + " << (integralXS[theXSType][Bin+1]-integralXS[theXSType][Bin]) << " / " << (bb-aa) << " *" << (eVEnergy-aa) << " * " << 1.e-16*CLHEP::cm2 << G4endl;;
224
225 // SIGMA += NbOfAtomsPerVolume[elm] * crossSection;
226 SIGMA = NbOfMoleculesPerVolume * crossSection;
227 if( verboseLevel >= 2 ) G4cout << GetName() << " ADDING SIGMA " << SIGMA << " NAtoms " << NbOfMoleculesPerVolume
228 << " Bin " << Bin << " TOTAL " << aa << " " << bb
229 << " XS " << integralXS[theXSType][Bin] << " " << integralXS[theXSType][Bin+1] << G4endl;
230 }
231
232 //-}
233
234 fValue = SIGMA > DBL_MIN ? 1./SIGMA : DBL_MAX;
235 }
236
237 ptrVector->PutValue(ii, fValue);
238 if( verboseLevel >= 2 ) G4cout << GetName() << " BUILDXS " << ii << " : " << LowEdgeEnergy << " = " << fValue << G4endl;
239 }
240
241 theMeanFreePathTable->insertAt( matIdx , ptrVector ) ;
242}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetDensity() const
std::size_t GetIndex() const
void insertAt(std::size_t, G4PhysicsVector *)
const G4String & GetName() const
std::map< const G4Material *, G4double > theMolecularMass
std::map< const G4Material *, G4int > theNXSdat
#define DBL_MIN
Definition templates.hh:54
#define DBL_MAX
Definition templates.hh:62

Referenced by BuildPhysicsTable().

◆ BuildPhysicsTable()

void G4VLEPTSModel::BuildPhysicsTable ( const G4ParticleDefinition & aParticleType)

Definition at line 90 of file G4VLEPTSModel.cc.

91{
92 //CHECK IF PATH VARIABLE IS DEFINED
93 const char* path = G4FindDataDir("G4LEDATA");
94 if( path == nullptr ) {
95 G4Exception("G4VLEPTSModel",
96 "",
98 "variable G4LEDATA not defined");
99 }
100
101 // Build microscopic cross section table and mean free path table
102
103 G4String aParticleName = aParticleType.GetParticleName();
104
105 if (theMeanFreePathTable != nullptr) {
108 }
109
111
112 //LOOP TO MATERIALS IN GEOMETRY
113 const G4MaterialTable * materialTable = G4Material::GetMaterialTable() ;
114 std::vector<G4Material*>::const_iterator matite;
115 for( matite = materialTable->begin(); matite != materialTable->end(); matite++ ) {
116 const G4Material * aMaterial = (*matite);
117 G4String mateName = aMaterial->GetName();
118
119 //READ PARAMETERS FOR THIS MATERIAL
120 std::string dirName = std::string(path) + "/lepts/";
121 std::string fnParam = dirName + mateName + "." + aParticleName + ".param.dat";
122 std::string baseName = std::string(path) + "/lepts/" + mateName + "." + aParticleName;
123 G4bool bData = ReadParam( fnParam, aMaterial );
124 if( !bData ) continue; // MATERIAL NOT EXISTING, DO NOT READ OTHER FILES
125
126 //READ INTEGRAL CROSS SECTION FOR THIS MATERIAL
127 std::string fnIXS = baseName + ".IXS.dat";
128
129 std::map< G4int, std::vector<G4double> > integralXS = ReadIXS(fnIXS, aMaterial);
130 if( verboseLevel >= 2 ) G4cout << GetName() << " : " << theXSType << " " << mateName << " INTEGRALXS " << integralXS.size() << G4endl;
131
132 if( integralXS.empty() ) {
133 G4cerr << " Integral cross sections will be set to 0. for material " << mateName << G4endl;
135 ptrVector->PutValue(0, DBL_MAX);
136 ptrVector->PutValue(1, DBL_MAX);
137
138 std::size_t matIdx = aMaterial->GetIndex();
139 theMeanFreePathTable->insertAt( matIdx , ptrVector ) ;
140
141 } else {
142
143 if( verboseLevel >= 2 ) {
144 std::map< G4int, std::vector<G4double> >::const_iterator itei;
145 for( itei = integralXS.begin(); itei != integralXS.end(); itei++ ){
146 G4cout << GetName() << " : " << (*itei).first << " INTEGRALXS NDATA " << (*itei).second.size() << G4endl;
147 }
148 }
149
150 BuildMeanFreePathTable( aMaterial, integralXS );
151
152 std::string fnDXS = baseName + ".DXS.dat";
153 std::string fnRMT = baseName + ".RMT.dat";
154 std::string fnEloss = baseName + ".Eloss.dat";
155 std::string fnEloss2 = baseName + ".Eloss2.dat";
156
157 theDiffXS[aMaterial] = new G4LEPTSDiffXS(fnDXS);
158 if( !theDiffXS[aMaterial]->IsFileFound() ) {
159 G4Exception("G4VLEPTSModel::BuildPhysicsTable",
160 "",
162 (G4String("File not found :" + fnDXS).c_str()));
163 }
164
165 theRMTDistr[aMaterial] = new G4LEPTSDistribution();
166 theRMTDistr[aMaterial]->ReadFile(fnRMT);
167
168 theElostDistr[aMaterial] = new G4LEPTSElossDistr(fnEloss);
169 if( !theElostDistr[aMaterial]->IsFileFound() ) {
170 G4Exception("G4VLEPTSModel::BuildPhysicsTable",
171 "",
173 (G4String("File not found :" + fnEloss).c_str()));
174 }
175 }
176
177 }
178
179}
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::vector< G4Material * > G4MaterialTable
bool G4bool
Definition G4Types.hh:86
G4GLOB_DLL std::ostream G4cerr
static std::size_t GetNumberOfMaterials()
static G4MaterialTable * GetMaterialTable()
const G4String & GetName() const
const G4String & GetParticleName() const
G4bool ReadParam(G4String fileName, const G4Material *aMaterial)
std::map< const G4Material *, G4LEPTSDiffXS * > theDiffXS
std::map< const G4Material *, G4LEPTSElossDistr * > theElostDistr
virtual std::map< G4int, std::vector< G4double > > ReadIXS(G4String fileName, const G4Material *aMaterial)
std::map< const G4Material *, G4LEPTSDistribution * > theRMTDistr
void BuildMeanFreePathTable(const G4Material *aMaterial, std::map< G4int, std::vector< G4double > > &integralXS)

Referenced by G4LEPTSAttachmentModel::Initialise(), G4LEPTSDissociationModel::Initialise(), G4LEPTSElasticModel::Initialise(), G4LEPTSExcitationModel::Initialise(), G4LEPTSIonisationModel::Initialise(), G4LEPTSPositroniumModel::Initialise(), G4LEPTSRotExcitationModel::Initialise(), and G4LEPTSVibExcitationModel::Initialise().

◆ GetMeanFreePath()

G4double G4VLEPTSModel::GetMeanFreePath ( const G4Material * mate,
const G4ParticleDefinition * aParticle,
G4double kineticEnergy )

Definition at line 73 of file G4VLEPTSModel.cc.

76{
77 G4double MeanFreePath;
78
79 if( verboseLevel >= 3 ) G4cout << aMaterial->GetIndex() << " G4VLEPTSModel::GetMeanFreePath " << kineticEnergy << " > " << theHighestEnergyLimit << " < " << theLowestEnergyLimit << G4endl;
80 if (kineticEnergy > theHighestEnergyLimit || kineticEnergy < theLowestEnergyLimit)
81 MeanFreePath = DBL_MAX;
82 else
83 MeanFreePath = (*theMeanFreePathTable)(aMaterial->GetIndex())->Value(kineticEnergy);
84
85 return MeanFreePath;
86}

Referenced by G4LEPTSAttachmentModel::CrossSectionPerVolume(), G4LEPTSDissociationModel::CrossSectionPerVolume(), G4LEPTSElasticModel::CrossSectionPerVolume(), G4LEPTSExcitationModel::CrossSectionPerVolume(), G4LEPTSIonisationModel::CrossSectionPerVolume(), G4LEPTSPositroniumModel::CrossSectionPerVolume(), G4LEPTSRotExcitationModel::CrossSectionPerVolume(), and G4LEPTSVibExcitationModel::CrossSectionPerVolume().

◆ Init()

◆ operator=()

G4VLEPTSModel & G4VLEPTSModel::operator= ( const G4VLEPTSModel & right)

◆ ReadIXS()

std::map< G4int, std::vector< G4double > > G4VLEPTSModel::ReadIXS ( G4String fileName,
const G4Material * aMaterial )
protectedvirtual

Reimplemented in G4LEPTSExcitationModel.

Definition at line 364 of file G4VLEPTSModel.cc.

365{
366 std::map< G4int, std::vector<G4double> > integralXS; // process type - energy
367 //G4cout << "fnIXS (" << fnIXS << ")" << G4endl;
368
369 std::ifstream fin(fnIXS);
370 if (!fin.is_open()) {
371 G4Exception("G4VLEPTSModel::ReadIXS",
372 "",
374 (G4String("File not found: ")+ fnIXS).c_str());
375 return integralXS;
376 }
377
378 G4int nXSdat, nXSsub;
379 fin >> nXSdat >> nXSsub;
380 if( verboseLevel >= 1 ) G4cout << "Read IXS (" << fnIXS << ")\t nXSdat: " << nXSdat
381 << " nXSsub: " << nXSsub << G4endl;
382 theNXSdat[aMaterial] = nXSdat;
383 theNXSsub[aMaterial] = nXSsub;
384
385 G4double xsdat;
386 for (G4int ip=0; ip<=nXSsub; ip++) {
387 integralXS[ip].push_back(0.);
388 }
389 for (G4int ie=1; ie<=nXSdat; ie++) {
390 for (G4int ip=0; ip<=nXSsub; ip++) {
391 fin >> xsdat;
392 integralXS[ip].push_back(xsdat);
393 if( verboseLevel >= 3 ) G4cout << GetName() << " FILL IXS " << ip << " " << ie << " = " << integralXS[ip][ie] << " " << xsdat << G4endl;
394 // xsdat 1e-16*cm2
395 }
396 }
397 fin.close();
398
399 return integralXS;
400}
@ JustWarning
std::map< const G4Material *, G4int > theNXSsub

Referenced by BuildPhysicsTable(), and G4LEPTSExcitationModel::ReadIXS().

◆ ReadParam()

G4bool G4VLEPTSModel::ReadParam ( G4String fileName,
const G4Material * aMaterial )
protected

Definition at line 325 of file G4VLEPTSModel.cc.

326{
327 std::ifstream fin(fnParam);
328 if (!fin.is_open()) {
329 G4Exception("G4VLEPTSModel::ReadParam",
330 "",
332 (G4String("File not found: ")+ fnParam).c_str());
333 return false;
334 }
335
336 G4double IonisPot, IonisPotInt;
337
338 fin >> IonisPot >> IonisPotInt;
339 if( verboseLevel >= 1 ) G4cout << "Read param (" << fnParam << ")\t IonisPot: " << IonisPot
340 << " IonisPotInt: " << IonisPotInt << G4endl;
341
342 theIonisPot[aMaterial] = IonisPot * CLHEP::eV;
343 theIonisPotInt[aMaterial] = IonisPotInt * CLHEP::eV;
344
345 G4double MolecularMass = 0;
346 auto nelem = (G4int)aMaterial->GetNumberOfElements();
347 const G4int* atomsV = aMaterial->GetAtomsVector();
348 for( G4int ii = 0; ii < nelem; ++ii ) {
349 MolecularMass += aMaterial->GetElement(ii)->GetA()*atomsV[ii]/CLHEP::g;
350 // G4cout << " MMASS1 " << mmass/CLHEP::g << " " << aMaterial->GetElement(ii)->GetName() << " " << aMaterial->GetElement(ii)->GetA()/CLHEP::g << G4endl;
351 }
352 // G4cout << " MMASS " << MolecularMass << " " << MolecularMass*CLHEP::g << " ME " << mmass << " " << mmass/CLHEP::g << G4endl;
353 theMolecularMass[aMaterial] = MolecularMass* CLHEP::g/CLHEP::mole;
354 // theMolecularMass[aMaterial] = aMaterial->GetMassOfMolecule()*CLHEP::Avogadro; // Material mixtures do not calculate molecular mass
355
356 if( verboseLevel >= 1) G4cout << " IonisPot: " << IonisPot/CLHEP::eV << " eV "
357 << " IonisPotInt: " << IonisPotInt/CLHEP::eV << " eV"
358 << " MolecularMass " << MolecularMass/(CLHEP::g/CLHEP::mole) << " g/mole" << G4endl;
359
360 return true;
361}
G4double GetA() const
Definition G4Element.hh:127
const G4Element * GetElement(G4int iel) const
const G4int * GetAtomsVector() const
std::size_t GetNumberOfElements() const
std::map< const G4Material *, G4double > theIonisPotInt
std::map< const G4Material *, G4double > theIonisPot

Referenced by BuildPhysicsTable().

◆ SampleAngle()

G4double G4VLEPTSModel::SampleAngle ( const G4Material * aMaterial,
G4double e,
G4double el )

Definition at line 246 of file G4VLEPTSModel.cc.

247{
248 G4double x;
249
250 if( e < 10001) {
251 x = theDiffXS[aMaterial]->SampleAngleMT(e, el);
252 }
253 else {
254 G4double Ei = e; //incidente
255 G4double Ed = e -el; //dispersado
256
257 G4double Pi = std::sqrt( std::pow( (Ei/27.2/137),2) +2*Ei/27.2); //incidente
258 G4double Pd = std::sqrt( std::pow( (Ed/27.2/137),2) +2*Ed/27.2); //dispersado
259
260 G4double Kmin = Pi - Pd;
261 G4double Kmax = Pi + Pd;
262
263 G4double KR = theRMTDistr[aMaterial]->Sample(Kmin, Kmax); //sorteo mom. transf.
264
265 G4double co = (Pi*Pi + Pd*Pd - KR*KR) / (2*Pi*Pd); //cos ang. disp.
266 if( co > 1. ) co = 1.;
267 x = std::acos(co); //*360/twopi; //ang. dispers.
268 }
269 return(x);
270}

Referenced by SampleNewDirection(), and G4LEPTSElasticModel::SampleSecondaries().

◆ SampleEnergyLoss()

G4double G4VLEPTSModel::SampleEnergyLoss ( const G4Material * aMaterial,
G4double eMin,
G4double eMax )
protected

Definition at line 311 of file G4VLEPTSModel.cc.

312{
313 G4double el;
314 el = theElostDistr[aMaterial]->Sample(eMin/CLHEP::eV, eMax/CLHEP::eV)*CLHEP::eV;
315
316#ifdef DEBUG_LEPTS
317 if( verboseLevel >= 2 ) G4cout << aMaterial->GetName() <<"SampleEnergyLoss/eV " << el/CLHEP::eV << " eMax/eV " << eMax/CLHEP::eV << " "
318 << " " << GetName() << G4endl;
319#endif
320 return el;
321}

Referenced by G4LEPTSDissociationModel::SampleSecondaries(), G4LEPTSExcitationModel::SampleSecondaries(), and G4LEPTSIonisationModel::SampleSecondaries().

◆ SampleNewDirection() [1/2]

G4ThreeVector G4VLEPTSModel::SampleNewDirection ( const G4Material * aMaterial,
G4ThreeVector Dir,
G4double e,
G4double el )

Definition at line 273 of file G4VLEPTSModel.cc.

273 {
274
275 G4double x = SampleAngle(aMaterial, e, el);
276
277 G4double cosTeta = std::cos(x); //*twopi/360.0);
278 G4double sinTeta = std::sqrt(1.0-cosTeta*cosTeta);
279 G4double Phi = CLHEP::twopi * G4UniformRand() ;
280 G4double dirx = sinTeta*std::cos(Phi) , diry = sinTeta*std::sin(Phi) , dirz = cosTeta ;
281
282 G4ThreeVector P1Dir(dirx, diry, dirz);
283#ifdef DEBUG_LEPTS
284 if( verboseLevel >= 2 ) G4cout << " G4VLEPTSModel::SampleNewDirection " <<P1Dir << G4endl;
285#endif
286 P1Dir.rotateUz(P0Dir);
287#ifdef DEBUG_LEPTS
288 if( verboseLevel >= 2 ) G4cout << " G4VLEPTSModel::SampleNewDirection rotated " <<P1Dir << " " << P0Dir << G4endl;
289#endif
290
291 return(P1Dir);
292}
#define G4UniformRand()
Definition Randomize.hh:52
G4double SampleAngle(const G4Material *aMaterial, G4double e, G4double el)

Referenced by G4LEPTSDissociationModel::SampleSecondaries(), G4LEPTSElasticModel::SampleSecondaries(), G4LEPTSExcitationModel::SampleSecondaries(), G4LEPTSIonisationModel::SampleSecondaries(), G4LEPTSPositroniumModel::SampleSecondaries(), G4LEPTSRotExcitationModel::SampleSecondaries(), and G4LEPTSVibExcitationModel::SampleSecondaries().

◆ SampleNewDirection() [2/2]

G4ThreeVector G4VLEPTSModel::SampleNewDirection ( G4ThreeVector Dir,
G4double ang )

Definition at line 296 of file G4VLEPTSModel.cc.

297{
298 G4double cosTeta = std::cos(x); //*twopi/360.0);
299 G4double sinTeta = std::sqrt(1.0-cosTeta*cosTeta);
300 G4double Phi = CLHEP::twopi * G4UniformRand() ;
301 G4double dirx = sinTeta*std::cos(Phi) , diry = sinTeta*std::sin(Phi) , dirz = cosTeta ;
302
303 G4ThreeVector P1Dir( dirx,diry,dirz );
304 P1Dir.rotateUz(P0Dir);
305
306 return(P1Dir);
307}

Member Data Documentation

◆ isInitialised

G4bool G4VLEPTSModel::isInitialised {false}
protected

Definition at line 106 of file G4VLEPTSModel.hh.

106{false};

◆ theDiffXS

std::map<const G4Material*, G4LEPTSDiffXS*> G4VLEPTSModel::theDiffXS
protected

Definition at line 97 of file G4VLEPTSModel.hh.

Referenced by BuildPhysicsTable(), and SampleAngle().

◆ theElostDistr

std::map<const G4Material*, G4LEPTSElossDistr*> G4VLEPTSModel::theElostDistr
protected

Definition at line 100 of file G4VLEPTSModel.hh.

Referenced by BuildPhysicsTable(), and SampleEnergyLoss().

◆ theElostDistr2

std::map<const G4Material*, G4LEPTSDistribution*> G4VLEPTSModel::theElostDistr2
protected

Definition at line 101 of file G4VLEPTSModel.hh.

◆ theHighestEnergyLimit

G4double G4VLEPTSModel::theHighestEnergyLimit
protected

◆ theIonisPot

std::map<const G4Material*, G4double > G4VLEPTSModel::theIonisPot
protected

◆ theIonisPotInt

std::map<const G4Material*, G4double > G4VLEPTSModel::theIonisPotInt
protected

Definition at line 94 of file G4VLEPTSModel.hh.

Referenced by ReadParam(), and G4LEPTSIonisationModel::SampleSecondaries().

◆ theLowestEnergyLimit

◆ theMeanFreePathTable

G4PhysicsTable* G4VLEPTSModel::theMeanFreePathTable
protected

◆ theMolecularMass

std::map<const G4Material*, G4double > G4VLEPTSModel::theMolecularMass
protected

◆ theNumbBinTable

G4int G4VLEPTSModel::theNumbBinTable
protected

Definition at line 91 of file G4VLEPTSModel.hh.

Referenced by BuildMeanFreePathTable(), G4VLEPTSModel(), and Init().

◆ theNXSdat

std::map<const G4Material*, G4int> G4VLEPTSModel::theNXSdat
protected

◆ theNXSsub

std::map<const G4Material*, G4int> G4VLEPTSModel::theNXSsub
protected

Definition at line 104 of file G4VLEPTSModel.hh.

Referenced by ReadIXS().

◆ theRMTDistr

std::map<const G4Material*, G4LEPTSDistribution*> G4VLEPTSModel::theRMTDistr
protected

Definition at line 98 of file G4VLEPTSModel.hh.

Referenced by BuildPhysicsTable(), and SampleAngle().

◆ theXSType

◆ verboseLevel


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