Geant4 11.1.1
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 *, 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 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 *)
 
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 = 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
 
size_t currentCoupleIndex = 0
 
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),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 175 of file G4VLEPTSModel.cc.

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

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

70{
71 G4double MeanFreePath;
72
73 if( verboseLevel >= 3 ) G4cout << aMaterial->GetIndex() << " G4VLEPTSModel::GetMeanFreePath " << kineticEnergy << " > " << theHighestEnergyLimit << " < " << theLowestEnergyLimit << G4endl;
74 if (kineticEnergy > theHighestEnergyLimit || kineticEnergy < theLowestEnergyLimit)
75 MeanFreePath = DBL_MAX;
76 else
77 MeanFreePath = (*theMeanFreePathTable)(aMaterial->GetIndex())->Value(kineticEnergy);
78
79 return MeanFreePath;
80}
virtual G4double Value(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
Definition: G4VEmModel.cc:349

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

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

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

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

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

◆ SampleEnergyLoss()

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

Definition at line 305 of file G4VLEPTSModel.cc.

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

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

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

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

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: