Geant4 10.7.0
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 ()
 
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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 LPMFlag () 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 SetLPMFlag (G4bool val)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (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
 
XSType theXSType
 
G4int verboseLevel
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

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),isInitialised(false)
32{
34
36
37 verboseLevel = 0;
38}
G4int theNumbBinTable
G4PhysicsTable * theMeanFreePathTable
G4bool isInitialised

◆ ~G4VLEPTSModel()

G4VLEPTSModel::~G4VLEPTSModel ( )

Definition at line 42 of file G4VLEPTSModel.cc.

43{
44
48 }
49}
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 178 of file G4VLEPTSModel.cc.

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

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

70{
71
72 G4double MeanFreePath;
73 G4bool isOutRange ;
74
75 if( verboseLevel >= 3 ) G4cout << aMaterial->GetIndex() << " G4VLEPTSModel::GetMeanFreePath " << kineticEnergy << " > " << theHighestEnergyLimit << " < " << theLowestEnergyLimit << G4endl;
76 if (kineticEnergy > theHighestEnergyLimit || kineticEnergy < theLowestEnergyLimit)
77 MeanFreePath = DBL_MAX;
78 else
79 MeanFreePath = (*theMeanFreePathTable)(aMaterial->GetIndex())->
80 GetValue(kineticEnergy, isOutRange);
81
82 return MeanFreePath;
83}

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 361 of file G4VLEPTSModel.cc.

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

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

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

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

◆ SampleEnergyLoss()

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

Definition at line 308 of file G4VLEPTSModel.cc.

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

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 270 of file G4VLEPTSModel.cc.

270 {
271
272 G4double x = SampleAngle(aMaterial, e, el);
273
274 G4double cosTeta = std::cos(x); //*twopi/360.0);
275 G4double sinTeta = std::sqrt(1.0-cosTeta*cosTeta);
276 G4double Phi = CLHEP::twopi * G4UniformRand() ;
277 G4double dirx = sinTeta*std::cos(Phi) , diry = sinTeta*std::sin(Phi) , dirz = cosTeta ;
278
279 G4ThreeVector P1Dir(dirx, diry, dirz);
280#ifdef DEBUG_LEPTS
281 if( verboseLevel >= 2 ) G4cout << " G4VLEPTSModel::SampleNewDirection " <<P1Dir << G4endl;
282#endif
283 P1Dir.rotateUz(P0Dir);
284#ifdef DEBUG_LEPTS
285 if( verboseLevel >= 2 ) G4cout << " G4VLEPTSModel::SampleNewDirection rotated " <<P1Dir << " " << P0Dir << G4endl;
286#endif
287
288 return(P1Dir);
289}
#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 293 of file G4VLEPTSModel.cc.

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

Member Data Documentation

◆ isInitialised

G4bool G4VLEPTSModel::isInitialised
protected

Definition at line 106 of file G4VLEPTSModel.hh.

◆ 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

Definition at line 90 of file G4VLEPTSModel.hh.

Referenced by BuildMeanFreePathTable(), BuildPhysicsTable(), GetMeanFreePath(), and Init().

◆ 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: