Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4AdjointhIonisationModel Class Reference

#include <G4AdjointhIonisationModel.hh>

Inheritance diagram for G4AdjointhIonisationModel:
G4VEmAdjointModel

Public Member Functions

G4double AdjointCrossSection (const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool isScatProjToProj) override
 
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerAtomForScatProj (G4double kinEnergyProd, G4double Z, G4double A=0., G4int nbin_pro_decade=10)
 
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerAtomForSecond (G4double kinEnergyProd, G4double Z, G4double A=0., G4int nbin_pro_decade=10)
 
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerVolumeForScatProj (G4Material *aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade=10)
 
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerVolumeForSecond (G4Material *aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade=10)
 
void DefineCurrentMaterial (const G4MaterialCutsCouple *couple)
 
void DefineDirectEMModel (G4VEmModel *aModel)
 
virtual G4double DiffCrossSectionPerAtomPrimToScatPrim (G4double kinEnergyProj, G4double kinEnergyScatProj, G4double Z, G4double A=0.)
 
G4double DiffCrossSectionPerAtomPrimToSecond (G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.) override
 
virtual G4double DiffCrossSectionPerVolumePrimToScatPrim (const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyScatProj)
 
virtual G4double DiffCrossSectionPerVolumePrimToSecond (const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
 
 G4AdjointhIonisationModel (G4AdjointhIonisationModel &)=delete
 
 G4AdjointhIonisationModel (G4ParticleDefinition *pDef)
 
G4ParticleDefinitionGetAdjointEquivalentOfDirectPrimaryParticleDefinition ()
 
G4ParticleDefinitionGetAdjointEquivalentOfDirectSecondaryParticleDefinition ()
 
G4bool GetApplyCutInRange ()
 
G4double GetHighEnergyLimit ()
 
G4double GetLowEnergyLimit ()
 
G4String GetName ()
 
G4double GetSecondAdjEnergyMaxForProdToProj (G4double primAdjEnergy) override
 
G4double GetSecondAdjEnergyMaxForScatProjToProj (G4double primAdjEnergy) override
 
G4double GetSecondAdjEnergyMinForProdToProj (G4double primAdjEnergy) override
 
G4double GetSecondAdjEnergyMinForScatProjToProj (G4double primAdjEnergy, G4double tcut=0.) override
 
G4bool GetSecondPartOfSameType ()
 
G4bool GetUseMatrix ()
 
G4bool GetUseMatrixPerElement ()
 
G4bool GetUseOnlyOneMatrixForAllElements ()
 
G4AdjointhIonisationModeloperator= (const G4AdjointhIonisationModel &right)=delete
 
void RapidSampleSecondaries (const G4Track &aTrack, G4bool isScatProjToProj, G4ParticleChange *fParticleChange)
 
void SampleSecondaries (const G4Track &aTrack, G4bool isScatProjToProj, G4ParticleChange *fParticleChange) override
 
void SetAdditionalWeightCorrectionFactorForPostStepOutsideModel (G4double factor)
 
void SetAdjointEquivalentOfDirectPrimaryParticleDefinition (G4ParticleDefinition *aPart)
 
void SetAdjointEquivalentOfDirectSecondaryParticleDefinition (G4ParticleDefinition *aPart)
 
void SetApplyCutInRange (G4bool aBool)
 
void SetCorrectWeightForPostStepInModel (G4bool aBool)
 
virtual void SetCSBiasingFactor (G4double aVal)
 
void SetCSMatrices (std::vector< G4AdjointCSMatrix * > *Vec1CSMatrix, std::vector< G4AdjointCSMatrix * > *Vec2CSMatrix)
 
void SetHighEnergyLimit (G4double aVal)
 
void SetLowEnergyLimit (G4double aVal)
 
void SetSecondPartOfSameType (G4bool aBool)
 
void SetUseMatrix (G4bool aBool)
 
void SetUseMatrixPerElement (G4bool aBool)
 
void SetUseOnlyOneMatrixForAllElements (G4bool aBool)
 
 ~G4AdjointhIonisationModel () override
 

Protected Member Functions

virtual void CorrectPostStepWeight (G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool isScatProjToProj)
 
G4double DiffCrossSectionFunction1 (G4double kinEnergyProj)
 
G4double DiffCrossSectionFunction2 (G4double kinEnergyProj)
 
G4double SampleAdjSecEnergyFromCSMatrix (G4double prim_energy, G4bool isScatProjToProj)
 
G4double SampleAdjSecEnergyFromCSMatrix (size_t MatrixIndex, G4double prim_energy, G4bool isScatProjToProj)
 
virtual G4double SampleAdjSecEnergyFromDiffCrossSectionPerAtom (G4double prim_energy, G4bool isScatProjToProj)
 
void SelectCSMatrix (G4bool isScatProjToProj)
 

Protected Attributes

G4ParticleDefinitionfAdjEquivDirectPrimPart = nullptr
 
G4ParticleDefinitionfAdjEquivDirectSecondPart = nullptr
 
G4bool fApplyCutInRange = true
 
G4int fASelectedNucleus = 0
 
G4double fCsBiasingFactor = 1.
 
G4AdjointCSManagerfCSManager
 
std::vector< G4AdjointCSMatrix * > * fCSMatrixProdToProjBackScat = nullptr
 
std::vector< G4AdjointCSMatrix * > * fCSMatrixProjToProjBackScat = nullptr
 
size_t fCSMatrixUsed = 0
 
G4MaterialCutsCouplefCurrentCouple = nullptr
 
G4MaterialfCurrentMaterial = nullptr
 
G4VEmModelfDirectModel = nullptr
 
G4ParticleDefinitionfDirectPrimaryPart = nullptr
 
std::vector< G4doublefElementCSProdToProj
 
std::vector< G4doublefElementCSScatProjToProj
 
G4double fHighEnergyLimit = 0.
 
G4bool fInModelWeightCorr
 
G4double fKinEnergyProdForIntegration = 0.
 
G4double fKinEnergyScatProjForIntegration = 0.
 
G4double fLastAdjointCSForProdToProj = 0.
 
G4double fLastAdjointCSForScatProjToProj = 0.
 
G4double fLastCS = 0.
 
G4double fLowEnergyLimit = 0.
 
const G4String fName
 
G4bool fOneMatrixForAllElements = false
 
G4double fOutsideWeightFactor = 1.
 
G4double fPreStepEnergy = 0.
 
G4bool fSecondPartSameType = false
 
G4MaterialfSelectedMaterial = nullptr
 
G4double fTcutPrim = 0.
 
G4double fTcutSecond = 0.
 
G4bool fUseMatrix = false
 
G4bool fUseMatrixPerElement = false
 
G4int fZSelectedNucleus = 0
 

Private Member Functions

void DefineProjectileProperty ()
 

Private Attributes

G4VEmModelfBraggDirectEMModel
 
G4double fFormFact = 0.
 
G4double fMagMoment2 = 0.
 
G4double fMass = 0.
 
G4double fMassRatio = 0.
 
G4double fOneMinusRatio2 = 0.
 
G4double fOnePlusRatio2 = 0.
 
G4double fSpin = 0.
 

Detailed Description

Definition at line 47 of file G4AdjointhIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4AdjointhIonisationModel() [1/2]

G4AdjointhIonisationModel::G4AdjointhIonisationModel ( G4ParticleDefinition pDef)
explicit

Definition at line 42 of file G4AdjointhIonisationModel.cc.

43 : G4VEmAdjointModel("Adjoint_hIonisation")
44{
45 fUseMatrix = true;
47 fApplyCutInRange = true;
49 fSecondPartSameType = false;
50
51 // The direct EM Model is taken as BetheBloch. It is only used for the
52 // computation of the differential cross section.
53 // The Bragg model could be used as an alternative as it offers the same
54 // differential cross section
55
59 fDirectPrimaryPart = pDef;
60
61 if(pDef == G4Proton::Proton())
62 {
64 }
65
67}
static G4AdjointElectron * AdjointElectron()
static G4AdjointProton * AdjointProton()
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4ParticleDefinition * fAdjEquivDirectSecondPart
G4ParticleDefinition * fAdjEquivDirectPrimPart
G4VEmModel * fDirectModel
G4ParticleDefinition * fDirectPrimaryPart
G4VEmAdjointModel(const G4String &nam)

References G4AdjointElectron::AdjointElectron(), G4AdjointProton::AdjointProton(), DefineProjectileProperty(), G4VEmAdjointModel::fAdjEquivDirectPrimPart, G4VEmAdjointModel::fAdjEquivDirectSecondPart, G4VEmAdjointModel::fApplyCutInRange, fBraggDirectEMModel, G4VEmAdjointModel::fDirectModel, G4VEmAdjointModel::fDirectPrimaryPart, G4VEmAdjointModel::fOneMatrixForAllElements, G4VEmAdjointModel::fSecondPartSameType, G4VEmAdjointModel::fUseMatrix, G4VEmAdjointModel::fUseMatrixPerElement, and G4Proton::Proton().

◆ ~G4AdjointhIonisationModel()

G4AdjointhIonisationModel::~G4AdjointhIonisationModel ( )
override

Definition at line 70 of file G4AdjointhIonisationModel.cc.

70{}

◆ G4AdjointhIonisationModel() [2/2]

G4AdjointhIonisationModel::G4AdjointhIonisationModel ( G4AdjointhIonisationModel )
delete

Member Function Documentation

◆ AdjointCrossSection()

G4double G4AdjointhIonisationModel::AdjointCrossSection ( const G4MaterialCutsCouple aCouple,
G4double  primEnergy,
G4bool  isScatProjToProj 
)
overridevirtual

Reimplemented from G4VEmAdjointModel.

Definition at line 397 of file G4AdjointhIonisationModel.cc.

400{
401 if(fUseMatrix)
402 return G4VEmAdjointModel::AdjointCrossSection(aCouple, primEnergy,
403 isScatProjToProj);
404 DefineCurrentMaterial(aCouple);
405
406 G4double Cross =
408
409 if(!isScatProjToProj)
410 {
411 G4double Emax_proj = GetSecondAdjEnergyMaxForProdToProj(primEnergy);
412 G4double Emin_proj = GetSecondAdjEnergyMinForProdToProj(primEnergy);
413 if(Emax_proj > Emin_proj && primEnergy > fTcutSecond)
414 {
415 Cross *= (1. / Emin_proj - 1. / Emax_proj) / primEnergy;
416 }
417 else
418 Cross = 0.;
419 }
420 else
421 {
422 G4double Emax_proj = GetSecondAdjEnergyMaxForScatProjToProj(primEnergy);
423 G4double Emin_proj =
425 G4double diff1 = Emin_proj - primEnergy;
426 G4double diff2 = Emax_proj - primEnergy;
427 G4double t1 =
428 (1. / diff1 + 1. / Emin_proj - 1. / diff2 - 1. / Emax_proj) / primEnergy;
429 G4double t2 =
430 2. * std::log(Emax_proj / Emin_proj) / primEnergy / primEnergy;
431 Cross *= (t1 + t2);
432 }
433 fLastCS = Cross;
434 return Cross;
435}
double G4double
Definition: G4Types.hh:83
G4double GetSecondAdjEnergyMinForScatProjToProj(G4double primAdjEnergy, G4double tcut=0.) override
G4double GetSecondAdjEnergyMaxForScatProjToProj(G4double primAdjEnergy) override
G4double GetSecondAdjEnergyMinForProdToProj(G4double primAdjEnergy) override
G4double GetSecondAdjEnergyMaxForProdToProj(G4double primAdjEnergy) override
G4double GetElectronDensity() const
Definition: G4Material.hh:213
G4Material * fCurrentMaterial
void DefineCurrentMaterial(const G4MaterialCutsCouple *couple)
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool isScatProjToProj)
int twopi_mc2_rcl2
Definition: hepunit.py:293

References G4VEmAdjointModel::AdjointCrossSection(), G4VEmAdjointModel::DefineCurrentMaterial(), G4VEmAdjointModel::fCurrentMaterial, G4VEmAdjointModel::fLastCS, fMass, G4VEmAdjointModel::fTcutSecond, G4VEmAdjointModel::fUseMatrix, G4Material::GetElectronDensity(), GetSecondAdjEnergyMaxForProdToProj(), GetSecondAdjEnergyMaxForScatProjToProj(), GetSecondAdjEnergyMinForProdToProj(), GetSecondAdjEnergyMinForScatProjToProj(), and source.hepunit::twopi_mc2_rcl2.

◆ ComputeAdjointCrossSectionVectorPerAtomForScatProj()

std::vector< std::vector< G4double > * > G4VEmAdjointModel::ComputeAdjointCrossSectionVectorPerAtomForScatProj ( G4double  kinEnergyProd,
G4double  Z,
G4double  A = 0.,
G4int  nbin_pro_decade = 10 
)
inherited

Definition at line 252 of file G4VEmAdjointModel.cc.

255{
257 integral;
260 fKinEnergyScatProjForIntegration = kinEnergyScatProj;
261
262 // compute the vector of integrated cross sections
263 G4double minEProj = GetSecondAdjEnergyMinForScatProjToProj(kinEnergyScatProj);
264 G4double maxEProj = GetSecondAdjEnergyMaxForScatProjToProj(kinEnergyScatProj);
265 G4double dEmax = maxEProj - kinEnergyScatProj;
266 G4double dEmin = GetLowEnergyLimit();
267 G4double dE1 = dEmin;
268 G4double dE2 = dEmin;
269
270 std::vector<G4double>* log_ESec_vector = new std::vector<G4double>();
271 std::vector<G4double>* log_Prob_vector = new std::vector<G4double>();
272 log_ESec_vector->push_back(std::log(dEmin));
273 log_Prob_vector->push_back(-50.);
274 G4int nbins = std::max(G4int(std::log10(dEmax / dEmin)) * nbin_pro_decade, 5);
275 G4double fE = std::pow(dEmax / dEmin, 1. / nbins);
276
277 G4double int_cross_section = 0.;
278 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
279 while(dE1 < dEmax * 0.9999999999999)
280 {
281 dE2 = dE1 * fE;
282 int_cross_section +=
283 integral.Simpson(this, &G4VEmAdjointModel::DiffCrossSectionFunction2,
284 minEProj + dE1, std::min(minEProj + dE2, maxEProj), 5);
285 log_ESec_vector->push_back(std::log(std::min(dE2, maxEProj - minEProj)));
286 log_Prob_vector->push_back(std::log(int_cross_section));
287 dE1 = dE2;
288 }
289
290 std::vector<std::vector<G4double>*> res_mat;
291 if(int_cross_section > 0.)
292 {
293 res_mat.push_back(log_ESec_vector);
294 res_mat.push_back(log_Prob_vector);
295 }
296 else {
297 delete log_ESec_vector;
298 delete log_Prob_vector;
299 }
300
301 return res_mat;
302}
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4double fKinEnergyScatProjForIntegration
virtual G4double GetSecondAdjEnergyMinForScatProjToProj(G4double primAdjEnergy, G4double tcut=0.)
virtual G4double GetSecondAdjEnergyMaxForScatProjToProj(G4double primAdjEnergy)
G4double GetLowEnergyLimit()
G4double DiffCrossSectionFunction2(G4double kinEnergyProj)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
int G4lrint(double ad)
Definition: templates.hh:134

References A, G4VEmAdjointModel::DiffCrossSectionFunction2(), G4VEmAdjointModel::fASelectedNucleus, G4VEmAdjointModel::fKinEnergyScatProjForIntegration, G4VEmAdjointModel::fZSelectedNucleus, G4lrint(), G4VEmAdjointModel::G4VEmAdjointModel(), G4VEmAdjointModel::GetLowEnergyLimit(), G4VEmAdjointModel::GetSecondAdjEnergyMaxForScatProjToProj(), G4VEmAdjointModel::GetSecondAdjEnergyMinForScatProjToProj(), G4INCL::Math::max(), G4INCL::Math::min(), and Z.

Referenced by G4AdjointCSManager::BuildCrossSectionsModelAndElement().

◆ ComputeAdjointCrossSectionVectorPerAtomForSecond()

std::vector< std::vector< G4double > * > G4VEmAdjointModel::ComputeAdjointCrossSectionVectorPerAtomForSecond ( G4double  kinEnergyProd,
G4double  Z,
G4double  A = 0.,
G4int  nbin_pro_decade = 10 
)
inherited

Definition at line 198 of file G4VEmAdjointModel.cc.

201{
203 integral;
206 fKinEnergyProdForIntegration = kinEnergyProd;
207
208 // compute the vector of integrated cross sections
209 G4double minEProj = GetSecondAdjEnergyMinForProdToProj(kinEnergyProd);
210 G4double maxEProj = GetSecondAdjEnergyMaxForProdToProj(kinEnergyProd);
211 G4double E1 = minEProj;
212 std::vector<G4double>* log_ESec_vector = new std::vector<G4double>();
213 std::vector<G4double>* log_Prob_vector = new std::vector<G4double>();
214 log_ESec_vector->push_back(std::log(E1));
215 log_Prob_vector->push_back(-50.);
216
217 G4double E2 =
218 std::pow(10., G4double(G4int(std::log10(minEProj) * nbin_pro_decade) + 1) /
219 nbin_pro_decade);
220 G4double fE = std::pow(10., 1. / nbin_pro_decade);
221
222 if(std::pow(fE, 5.) > (maxEProj / minEProj))
223 fE = std::pow(maxEProj / minEProj, 0.2);
224
225 G4double int_cross_section = 0.;
226 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
227 while(E1 < maxEProj * 0.9999999)
228 {
229 int_cross_section +=
230 integral.Simpson(this, &G4VEmAdjointModel::DiffCrossSectionFunction1, E1,
231 std::min(E2, maxEProj * 0.99999999), 5);
232 log_ESec_vector->push_back(std::log(std::min(E2, maxEProj)));
233 log_Prob_vector->push_back(std::log(int_cross_section));
234 E1 = E2;
235 E2 *= fE;
236 }
237 std::vector<std::vector<G4double>*> res_mat;
238 if(int_cross_section > 0.)
239 {
240 res_mat.push_back(log_ESec_vector);
241 res_mat.push_back(log_Prob_vector);
242 }
243 else {
244 delete log_ESec_vector;
245 delete log_Prob_vector;
246 }
247 return res_mat;
248}
virtual G4double GetSecondAdjEnergyMaxForProdToProj(G4double primAdjEnergy)
G4double DiffCrossSectionFunction1(G4double kinEnergyProj)
G4double fKinEnergyProdForIntegration
virtual G4double GetSecondAdjEnergyMinForProdToProj(G4double primAdjEnergy)

References A, G4VEmAdjointModel::DiffCrossSectionFunction1(), G4VEmAdjointModel::fASelectedNucleus, G4VEmAdjointModel::fKinEnergyProdForIntegration, G4VEmAdjointModel::fZSelectedNucleus, G4lrint(), G4VEmAdjointModel::G4VEmAdjointModel(), G4VEmAdjointModel::GetSecondAdjEnergyMaxForProdToProj(), G4VEmAdjointModel::GetSecondAdjEnergyMinForProdToProj(), G4INCL::Math::min(), and Z.

Referenced by G4AdjointCSManager::BuildCrossSectionsModelAndElement().

◆ ComputeAdjointCrossSectionVectorPerVolumeForScatProj()

std::vector< std::vector< G4double > * > G4VEmAdjointModel::ComputeAdjointCrossSectionVectorPerVolumeForScatProj ( G4Material aMaterial,
G4double  kinEnergyProd,
G4int  nbin_pro_decade = 10 
)
inherited

Definition at line 360 of file G4VEmAdjointModel.cc.

363{
365 integral;
366 fSelectedMaterial = aMaterial;
367 fKinEnergyScatProjForIntegration = kinEnergyScatProj;
368
369 // compute the vector of integrated cross sections
370 G4double minEProj = GetSecondAdjEnergyMinForScatProjToProj(kinEnergyScatProj);
371 G4double maxEProj = GetSecondAdjEnergyMaxForScatProjToProj(kinEnergyScatProj);
372
373 G4double dEmax = maxEProj - kinEnergyScatProj;
374 G4double dEmin = GetLowEnergyLimit();
375 G4double dE1 = dEmin;
376 G4double dE2 = dEmin;
377
378 std::vector<G4double>* log_ESec_vector = new std::vector<G4double>();
379 std::vector<G4double>* log_Prob_vector = new std::vector<G4double>();
380 log_ESec_vector->push_back(std::log(dEmin));
381 log_Prob_vector->push_back(-50.);
382 G4int nbins = std::max(int(std::log10(dEmax / dEmin)) * nbin_pro_decade, 5);
383 G4double fE = std::pow(dEmax / dEmin, 1. / nbins);
384
385 G4double int_cross_section = 0.;
386 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
387 while(dE1 < dEmax * 0.9999999999999)
388 {
389 dE2 = dE1 * fE;
390 int_cross_section +=
391 integral.Simpson(this, &G4VEmAdjointModel::DiffCrossSectionFunction2,
392 minEProj + dE1, std::min(minEProj + dE2, maxEProj), 5);
393 log_ESec_vector->push_back(std::log(std::min(dE2, maxEProj - minEProj)));
394 log_Prob_vector->push_back(std::log(int_cross_section));
395 dE1 = dE2;
396 }
397
398 std::vector<std::vector<G4double>*> res_mat;
399 if(int_cross_section > 0.)
400 {
401 res_mat.push_back(log_ESec_vector);
402 res_mat.push_back(log_Prob_vector);
403 }
404 else {
405 delete log_ESec_vector;
406 delete log_Prob_vector;
407 }
408
409 return res_mat;
410}
G4Material * fSelectedMaterial

References G4VEmAdjointModel::DiffCrossSectionFunction2(), G4VEmAdjointModel::fKinEnergyScatProjForIntegration, G4VEmAdjointModel::fSelectedMaterial, G4VEmAdjointModel::G4VEmAdjointModel(), G4VEmAdjointModel::GetLowEnergyLimit(), G4VEmAdjointModel::GetSecondAdjEnergyMaxForScatProjToProj(), G4VEmAdjointModel::GetSecondAdjEnergyMinForScatProjToProj(), G4INCL::Math::max(), and G4INCL::Math::min().

Referenced by G4AdjointCSManager::BuildCrossSectionsModelAndMaterial().

◆ ComputeAdjointCrossSectionVectorPerVolumeForSecond()

std::vector< std::vector< G4double > * > G4VEmAdjointModel::ComputeAdjointCrossSectionVectorPerVolumeForSecond ( G4Material aMaterial,
G4double  kinEnergyProd,
G4int  nbin_pro_decade = 10 
)
inherited

Definition at line 306 of file G4VEmAdjointModel.cc.

309{
311 integral;
312 fSelectedMaterial = aMaterial;
313 fKinEnergyProdForIntegration = kinEnergyProd;
314
315 // compute the vector of integrated cross sections
316 G4double minEProj = GetSecondAdjEnergyMinForProdToProj(kinEnergyProd);
317 G4double maxEProj = GetSecondAdjEnergyMaxForProdToProj(kinEnergyProd);
318 G4double E1 = minEProj;
319 std::vector<G4double>* log_ESec_vector = new std::vector<G4double>();
320 std::vector<G4double>* log_Prob_vector = new std::vector<G4double>();
321 log_ESec_vector->push_back(std::log(E1));
322 log_Prob_vector->push_back(-50.);
323
324 G4double E2 =
325 std::pow(10., G4double(G4int(std::log10(minEProj) * nbin_pro_decade) + 1) /
326 nbin_pro_decade);
327 G4double fE = std::pow(10., 1. / nbin_pro_decade);
328
329 if(std::pow(fE, 5.) > (maxEProj / minEProj))
330 fE = std::pow(maxEProj / minEProj, 0.2);
331
332 G4double int_cross_section = 0.;
333 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
334 while(E1 < maxEProj * 0.9999999)
335 {
336 int_cross_section +=
337 integral.Simpson(this, &G4VEmAdjointModel::DiffCrossSectionFunction1, E1,
338 std::min(E2, maxEProj * 0.99999999), 5);
339 log_ESec_vector->push_back(std::log(std::min(E2, maxEProj)));
340 log_Prob_vector->push_back(std::log(int_cross_section));
341 E1 = E2;
342 E2 *= fE;
343 }
344 std::vector<std::vector<G4double>*> res_mat;
345
346 if(int_cross_section > 0.)
347 {
348 res_mat.push_back(log_ESec_vector);
349 res_mat.push_back(log_Prob_vector);
350 }
351 else {
352 delete log_ESec_vector;
353 delete log_Prob_vector;
354 }
355 return res_mat;
356}

References G4VEmAdjointModel::DiffCrossSectionFunction1(), G4VEmAdjointModel::fKinEnergyProdForIntegration, G4VEmAdjointModel::fSelectedMaterial, G4VEmAdjointModel::G4VEmAdjointModel(), G4VEmAdjointModel::GetSecondAdjEnergyMaxForProdToProj(), G4VEmAdjointModel::GetSecondAdjEnergyMinForProdToProj(), and G4INCL::Math::min().

Referenced by G4AdjointCSManager::BuildCrossSectionsModelAndMaterial().

◆ CorrectPostStepWeight()

void G4VEmAdjointModel::CorrectPostStepWeight ( G4ParticleChange fParticleChange,
G4double  old_weight,
G4double  adjointPrimKinEnergy,
G4double  projectileKinEnergy,
G4bool  isScatProjToProj 
)
protectedvirtualinherited

Reimplemented in G4AdjointIonIonisationModel, and G4AdjointPhotoElectricModel.

Definition at line 616 of file G4VEmAdjointModel.cc.

621{
622 G4double new_weight = old_weight;
623 G4double w_corr =
625
627 if(!isScatProjToProj)
629 if((adjointPrimKinEnergy - fPreStepEnergy) / fPreStepEnergy > 0.001)
630 {
631 G4double post_stepCS = AdjointCrossSection(
632 fCurrentCouple, adjointPrimKinEnergy, isScatProjToProj);
633 if(post_stepCS > 0. && fLastCS > 0.)
634 w_corr *= post_stepCS / fLastCS;
635 }
636
637 new_weight *= w_corr;
638 new_weight *= projectileKinEnergy / adjointPrimKinEnergy;
639 // This is needed due to the biasing of diff CS
640 // by the factor adjointPrimKinEnergy/projectileKinEnergy
641
642 fParticleChange->SetParentWeightByProcess(false);
643 fParticleChange->SetSecondaryWeightByProcess(false);
644 fParticleChange->ProposeParentWeight(new_weight);
645}
G4double GetPostStepWeightCorrection()
G4double fLastAdjointCSForScatProjToProj
G4MaterialCutsCouple * fCurrentCouple
G4AdjointCSManager * fCSManager
G4double fLastAdjointCSForProdToProj
void SetSecondaryWeightByProcess(G4bool)
void SetParentWeightByProcess(G4bool)
void ProposeParentWeight(G4double finalWeight)

References G4VEmAdjointModel::AdjointCrossSection(), G4VEmAdjointModel::fCsBiasingFactor, G4VEmAdjointModel::fCSManager, G4VEmAdjointModel::fCurrentCouple, G4VEmAdjointModel::fLastAdjointCSForProdToProj, G4VEmAdjointModel::fLastAdjointCSForScatProjToProj, G4VEmAdjointModel::fLastCS, G4VEmAdjointModel::fPreStepEnergy, G4AdjointCSManager::GetPostStepWeightCorrection(), G4VParticleChange::ProposeParentWeight(), G4VParticleChange::SetParentWeightByProcess(), and G4VParticleChange::SetSecondaryWeightByProcess().

Referenced by G4AdjointBremsstrahlungModel::SampleSecondaries(), G4AdjointComptonModel::SampleSecondaries(), G4AdjointeIonisationModel::SampleSecondaries(), and SampleSecondaries().

◆ DefineCurrentMaterial()

void G4VEmAdjointModel::DefineCurrentMaterial ( const G4MaterialCutsCouple couple)
inherited

Definition at line 684 of file G4VEmAdjointModel.cc.

686{
687 if(couple != fCurrentCouple)
688 {
689 fCurrentCouple = const_cast<G4MaterialCutsCouple*>(couple);
690 fCurrentMaterial = const_cast<G4Material*>(couple->GetMaterial());
691 size_t idx = 56;
692 fTcutSecond = 1.e-11;
694 {
696 idx = 0;
698 idx = 1;
700 idx = 2;
701 if(idx < 56)
702 {
703 const std::vector<G4double>* aVec =
705 idx);
706 fTcutSecond = (*aVec)[couple->GetIndex()];
707 }
708 }
709 }
710}
static G4AdjointGamma * AdjointGamma()
static G4AdjointPositron * AdjointPositron()
const G4Material * GetMaterial() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()

References G4AdjointElectron::AdjointElectron(), G4AdjointGamma::AdjointGamma(), G4AdjointPositron::AdjointPositron(), G4VEmAdjointModel::fAdjEquivDirectSecondPart, G4VEmAdjointModel::fCurrentCouple, G4VEmAdjointModel::fCurrentMaterial, G4VEmAdjointModel::fTcutSecond, G4ProductionCutsTable::GetEnergyCutsVector(), G4MaterialCutsCouple::GetIndex(), G4MaterialCutsCouple::GetMaterial(), and G4ProductionCutsTable::GetProductionCutsTable().

Referenced by G4VEmAdjointModel::AdjointCrossSection(), G4AdjointBremsstrahlungModel::AdjointCrossSection(), G4AdjointComptonModel::AdjointCrossSection(), AdjointCrossSection(), G4AdjointBremsstrahlungModel::RapidSampleSecondaries(), G4AdjointComptonModel::RapidSampleSecondaries(), RapidSampleSecondaries(), and G4AdjointBremsstrahlungModel::SampleSecondaries().

◆ DefineDirectEMModel()

void G4VEmAdjointModel::DefineDirectEMModel ( G4VEmModel aModel)
inlineinherited

Definition at line 160 of file G4VEmAdjointModel.hh.

160{ fDirectModel = aModel; }

References G4VEmAdjointModel::fDirectModel.

◆ DefineProjectileProperty()

void G4AdjointhIonisationModel::DefineProjectileProperty ( )
private

Definition at line 367 of file G4AdjointhIonisationModel.cc.

368{
369 // Slightly modified code taken from G4BetheBlochModel::SetParticle
371
375 fOnePlusRatio2 = (1. + fMassRatio) * (1. + fMassRatio);
376 fOneMinusRatio2 = (1. - fMassRatio) * (1. - fMassRatio);
378 (0.5 * eplus * hbar_Planck * c_squared);
379 fMagMoment2 = magmom * magmom - 1.0;
380 fFormFact = 0.0;
382 {
383 G4double x = 0.8426 * GeV;
384 if(fSpin == 0.0 && fMass < GeV)
385 {
386 x = 0.736 * GeV;
387 }
388 else if(fMass > GeV)
389 {
391 }
392 fFormFact = 2.0 * electron_mass_c2 / (x * x);
393 }
394}
static constexpr double eplus
Definition: G4SIunits.hh:184
static constexpr double GeV
Definition: G4SIunits.hh:203
G4double GetZ13(G4double Z) const
static G4NistManager * Instance()
G4double GetPDGMagneticMoment() const
const G4String & GetParticleName() const
string pname
Definition: eplot.py:33
float electron_mass_c2
Definition: hepunit.py:273
float proton_mass_c2
Definition: hepunit.py:274
float hbar_Planck
Definition: hepunit.py:263
float c_squared
Definition: hepunit.py:257

References source.hepunit::c_squared, source.hepunit::electron_mass_c2, eplus, G4VEmAdjointModel::fDirectPrimaryPart, fFormFact, fMagMoment2, fMass, fMassRatio, fOneMinusRatio2, fOnePlusRatio2, fSpin, G4ParticleDefinition::GetLeptonNumber(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMagneticMoment(), G4ParticleDefinition::GetPDGMass(), G4ParticleDefinition::GetPDGSpin(), G4NistManager::GetZ13(), GeV, source.hepunit::hbar_Planck, G4NistManager::Instance(), eplot::pname, and source.hepunit::proton_mass_c2.

Referenced by G4AdjointhIonisationModel().

◆ DiffCrossSectionFunction1()

G4double G4VEmAdjointModel::DiffCrossSectionFunction1 ( G4double  kinEnergyProj)
protectedinherited

◆ DiffCrossSectionFunction2()

G4double G4VEmAdjointModel::DiffCrossSectionFunction2 ( G4double  kinEnergyProj)
protectedinherited

Definition at line 176 of file G4VEmAdjointModel.cc.

177{
178 G4double bias_factor =
181 {
185 bias_factor;
186 }
187 else
188 {
190 fSelectedMaterial, kinEnergyProj,
192 bias_factor;
193 }
194}
virtual G4double DiffCrossSectionPerVolumePrimToScatPrim(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyScatProj)
virtual G4double DiffCrossSectionPerAtomPrimToScatPrim(G4double kinEnergyProj, G4double kinEnergyScatProj, G4double Z, G4double A=0.)

References G4VEmAdjointModel::DiffCrossSectionPerAtomPrimToScatPrim(), G4VEmAdjointModel::DiffCrossSectionPerVolumePrimToScatPrim(), G4VEmAdjointModel::fASelectedNucleus, G4VEmAdjointModel::fCsBiasingFactor, G4VEmAdjointModel::fKinEnergyScatProjForIntegration, G4VEmAdjointModel::fSelectedMaterial, G4VEmAdjointModel::fUseMatrixPerElement, and G4VEmAdjointModel::fZSelectedNucleus.

Referenced by G4VEmAdjointModel::ComputeAdjointCrossSectionVectorPerAtomForScatProj(), and G4VEmAdjointModel::ComputeAdjointCrossSectionVectorPerVolumeForScatProj().

◆ DiffCrossSectionPerAtomPrimToScatPrim()

G4double G4VEmAdjointModel::DiffCrossSectionPerAtomPrimToScatPrim ( G4double  kinEnergyProj,
G4double  kinEnergyScatProj,
G4double  Z,
G4double  A = 0. 
)
virtualinherited

Reimplemented in G4AdjointComptonModel.

Definition at line 105 of file G4VEmAdjointModel.cc.

107{
108 G4double kinEnergyProd = kinEnergyProj - kinEnergyScatProj;
109 G4double dSigmadEprod;
110 if(kinEnergyProd <= 0.)
111 dSigmadEprod = 0.;
112 else
113 dSigmadEprod =
114 DiffCrossSectionPerAtomPrimToSecond(kinEnergyProj, kinEnergyProd, Z, A);
115 return dSigmadEprod;
116}

References A, G4VEmAdjointModel::DiffCrossSectionPerAtomPrimToSecond(), and Z.

Referenced by G4VEmAdjointModel::DiffCrossSectionFunction2(), and G4VEmAdjointModel::SampleAdjSecEnergyFromDiffCrossSectionPerAtom().

◆ DiffCrossSectionPerAtomPrimToSecond()

G4double G4AdjointhIonisationModel::DiffCrossSectionPerAtomPrimToSecond ( G4double  kinEnergyProj,
G4double  kinEnergyProd,
G4double  Z,
G4double  A = 0. 
)
overridevirtual

Reimplemented from G4VEmAdjointModel.

Definition at line 280 of file G4AdjointhIonisationModel.cc.

282{ // Probably here the Bragg Model should be also used for
283 // kinEnergyProj/nuc < 2 MeV
284
285 G4double dSigmadEprod = 0.;
286 G4double Emax_proj = GetSecondAdjEnergyMaxForProdToProj(kinEnergyProd);
287 G4double Emin_proj = GetSecondAdjEnergyMinForProdToProj(kinEnergyProd);
288
289 // the produced particle should have a kinetic energy smaller than the
290 // projectile
291 if(kinEnergyProj > Emin_proj && kinEnergyProj <= Emax_proj)
292 {
293 G4double Tmax = kinEnergyProj;
294 G4double E1 = kinEnergyProd;
295 //1.0006 factor seems to give the best diff CS, important impact on proton correction factor
296 G4double E2 = kinEnergyProd *1.0006;
297 G4double sigma1, sigma2;
298 if(kinEnergyProj > 2. * MeV)
299 {
301 fDirectPrimaryPart, kinEnergyProj, Z, A, E1, 1.e20);
303 fDirectPrimaryPart, kinEnergyProj, Z, A, E2, 1.e20);
304 }
305 else
306 {
308 fDirectPrimaryPart, kinEnergyProj, Z, A, E1, 1.e20);
310 fDirectPrimaryPart, kinEnergyProj, Z, A, E2, 1.e20);
311 }
312
313 dSigmadEprod = (sigma1 - sigma2) / (E2 - E1);
314 if(dSigmadEprod > 1.)
315 {
316 G4cout << "sigma1 " << kinEnergyProj / MeV << '\t' << kinEnergyProd / MeV
317 << '\t' << sigma1 << G4endl;
318 G4cout << "sigma2 " << kinEnergyProj / MeV << '\t' << kinEnergyProd / MeV
319 << '\t' << sigma2 << G4endl;
320 G4cout << "dsigma " << kinEnergyProj / MeV << '\t' << kinEnergyProd / MeV
321 << '\t' << dSigmadEprod << G4endl;
322 }
323
324 // correction of differential cross section at high energy to correct for
325 // the suppression of particle at secondary at high energy used in the Bethe
326 // Bloch Model. This correction consists of multiplying by g the probability
327 // function used to test the rejection of a secondary. Source code taken
328 // from G4BetheBlochModel::SampleSecondaries
329 G4double deltaKinEnergy = kinEnergyProd;
330
331 // projectile formfactor - suppression of high energy
332 // delta-electron production at high energy
333 G4double x = fFormFact * deltaKinEnergy;
334 if(x > 1.e-6)
335 {
336 G4double totEnergy = kinEnergyProj + fMass;
337 G4double etot2 = totEnergy * totEnergy;
338 G4double beta2 = kinEnergyProj * (kinEnergyProj + 2.0 * fMass) / etot2;
339 G4double f = 1.0 - beta2 * deltaKinEnergy / Tmax;
340 G4double f1 = 0.0;
341 if(0.5 == fSpin)
342 {
343 f1 = 0.5 * deltaKinEnergy * deltaKinEnergy / etot2;
344 f += f1;
345 }
346 G4double x1 = 1.0 + x;
347 G4double gg = 1.0 / (x1 * x1);
348 if(0.5 == fSpin)
349 {
350 G4double x2 = 0.5 * electron_mass_c2 * deltaKinEnergy / (fMass * fMass);
351 gg *= (1.0 + fMagMoment2 * (x2 - f1 / f) / (1.0 + x2));
352 }
353 if(gg > 1.0)
354 {
355 G4cout << "### G4BetheBlochModel in Adjoint Sim WARNING: g= " << g
356 << G4endl;
357 gg = 1.;
358 }
359 dSigmadEprod *= gg;
360 }
361 }
362
363 return dSigmadEprod;
364}
static constexpr double g
Definition: G4SIunits.hh:168
static constexpr double MeV
Definition: G4SIunits.hh:200
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:341

References A, G4VEmModel::ComputeCrossSectionPerAtom(), source.hepunit::electron_mass_c2, fBraggDirectEMModel, G4VEmAdjointModel::fDirectModel, G4VEmAdjointModel::fDirectPrimaryPart, fFormFact, fMagMoment2, fMass, fSpin, g, G4cout, G4endl, GetSecondAdjEnergyMaxForProdToProj(), GetSecondAdjEnergyMinForProdToProj(), MeV, and Z.

Referenced by RapidSampleSecondaries().

◆ DiffCrossSectionPerVolumePrimToScatPrim()

G4double G4VEmAdjointModel::DiffCrossSectionPerVolumePrimToScatPrim ( const G4Material aMaterial,
G4double  kinEnergyProj,
G4double  kinEnergyScatProj 
)
virtualinherited

Definition at line 140 of file G4VEmAdjointModel.cc.

143{
144 G4double kinEnergyProd = kinEnergyProj - kinEnergyScatProj;
145 G4double dSigmadEprod;
146 if(kinEnergyProd <= 0.)
147 dSigmadEprod = 0.;
148 else
150 aMaterial, kinEnergyProj, kinEnergyProd);
151 return dSigmadEprod;
152}

References G4VEmAdjointModel::DiffCrossSectionPerVolumePrimToSecond().

Referenced by G4VEmAdjointModel::DiffCrossSectionFunction2(), and G4AdjointeIonisationModel::SampleSecondaries().

◆ DiffCrossSectionPerVolumePrimToSecond()

G4double G4VEmAdjointModel::DiffCrossSectionPerVolumePrimToSecond ( const G4Material aMaterial,
G4double  kinEnergyProj,
G4double  kinEnergyProd 
)
virtualinherited

Reimplemented in G4AdjointBremsstrahlungModel.

Definition at line 119 of file G4VEmAdjointModel.cc.

121{
122 G4double dSigmadEprod = 0.;
123 G4double Emax_proj = GetSecondAdjEnergyMaxForProdToProj(kinEnergyProd);
124 G4double Emin_proj = GetSecondAdjEnergyMinForProdToProj(kinEnergyProd);
125
126 if(kinEnergyProj > Emin_proj && kinEnergyProj <= Emax_proj)
127 {
128 G4double E1 = kinEnergyProd;
129 G4double E2 = kinEnergyProd * 1.0001;
131 aMaterial, fDirectPrimaryPart, kinEnergyProj, E1, 1.e20);
133 aMaterial, fDirectPrimaryPart, kinEnergyProj, E2, 1.e20);
134 dSigmadEprod = (sigma1 - sigma2) / (E2 - E1);
135 }
136 return dSigmadEprod;
137}
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:237

References G4VEmModel::CrossSectionPerVolume(), G4VEmAdjointModel::fDirectModel, G4VEmAdjointModel::fDirectPrimaryPart, G4VEmAdjointModel::GetSecondAdjEnergyMaxForProdToProj(), and G4VEmAdjointModel::GetSecondAdjEnergyMinForProdToProj().

Referenced by G4VEmAdjointModel::DiffCrossSectionFunction1(), G4VEmAdjointModel::DiffCrossSectionPerVolumePrimToScatPrim(), G4AdjointBremsstrahlungModel::DiffCrossSectionPerVolumePrimToSecond(), and G4AdjointeIonisationModel::SampleSecondaries().

◆ GetAdjointEquivalentOfDirectPrimaryParticleDefinition()

G4ParticleDefinition * G4VEmAdjointModel::GetAdjointEquivalentOfDirectPrimaryParticleDefinition ( )
inlineinherited

◆ GetAdjointEquivalentOfDirectSecondaryParticleDefinition()

G4ParticleDefinition * G4VEmAdjointModel::GetAdjointEquivalentOfDirectSecondaryParticleDefinition ( )
inlineinherited

◆ GetApplyCutInRange()

G4bool G4VEmAdjointModel::GetApplyCutInRange ( )
inlineinherited

◆ GetHighEnergyLimit()

G4double G4VEmAdjointModel::GetHighEnergyLimit ( )
inlineinherited

◆ GetLowEnergyLimit()

G4double G4VEmAdjointModel::GetLowEnergyLimit ( )
inlineinherited

◆ GetName()

G4String G4VEmAdjointModel::GetName ( )
inlineinherited

Definition at line 203 of file G4VEmAdjointModel.hh.

203{ return fName; }
const G4String fName

References G4VEmAdjointModel::fName.

◆ GetSecondAdjEnergyMaxForProdToProj()

G4double G4AdjointhIonisationModel::GetSecondAdjEnergyMaxForProdToProj ( G4double  primAdjEnergy)
overridevirtual

◆ GetSecondAdjEnergyMaxForScatProjToProj()

G4double G4AdjointhIonisationModel::GetSecondAdjEnergyMaxForScatProjToProj ( G4double  primAdjEnergy)
overridevirtual

Reimplemented from G4VEmAdjointModel.

Definition at line 438 of file G4AdjointhIonisationModel.cc.

440{
441 G4double Tmax = primAdjEnergy * fOnePlusRatio2 /
442 (fOneMinusRatio2 - 2. * fMassRatio * primAdjEnergy / fMass);
443 return Tmax;
444}

References fMass, fMassRatio, fOneMinusRatio2, and fOnePlusRatio2.

Referenced by AdjointCrossSection(), and RapidSampleSecondaries().

◆ GetSecondAdjEnergyMinForProdToProj()

G4double G4AdjointhIonisationModel::GetSecondAdjEnergyMinForProdToProj ( G4double  primAdjEnergy)
overridevirtual

Reimplemented from G4VEmAdjointModel.

Definition at line 460 of file G4AdjointhIonisationModel.cc.

462{
463 G4double Tmin =
464 (2. * primAdjEnergy - 4. * fMass +
465 std::sqrt(4. * primAdjEnergy * primAdjEnergy + 16. * fMass * fMass +
466 8. * primAdjEnergy * fMass * (1. / fMassRatio + fMassRatio))) /
467 4.;
468 return Tmin;
469}

References fMass, and fMassRatio.

Referenced by AdjointCrossSection(), DiffCrossSectionPerAtomPrimToSecond(), and RapidSampleSecondaries().

◆ GetSecondAdjEnergyMinForScatProjToProj()

G4double G4AdjointhIonisationModel::GetSecondAdjEnergyMinForScatProjToProj ( G4double  primAdjEnergy,
G4double  tcut = 0. 
)
overridevirtual

Reimplemented from G4VEmAdjointModel.

Definition at line 447 of file G4AdjointhIonisationModel.cc.

449{
450 return primAdjEnergy + tcut;
451}

Referenced by AdjointCrossSection(), and RapidSampleSecondaries().

◆ GetSecondPartOfSameType()

G4bool G4VEmAdjointModel::GetSecondPartOfSameType ( )
inlineinherited

◆ GetUseMatrix()

G4bool G4VEmAdjointModel::GetUseMatrix ( )
inlineinherited

Definition at line 192 of file G4VEmAdjointModel.hh.

192{ return fUseMatrix; }

References G4VEmAdjointModel::fUseMatrix.

Referenced by G4AdjointCSManager::ComputeAdjointCS().

◆ GetUseMatrixPerElement()

G4bool G4VEmAdjointModel::GetUseMatrixPerElement ( )
inlineinherited

◆ GetUseOnlyOneMatrixForAllElements()

G4bool G4VEmAdjointModel::GetUseOnlyOneMatrixForAllElements ( )
inlineinherited

◆ operator=()

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

◆ RapidSampleSecondaries()

void G4AdjointhIonisationModel::RapidSampleSecondaries ( const G4Track aTrack,
G4bool  isScatProjToProj,
G4ParticleChange fParticleChange 
)

Definition at line 143 of file G4AdjointhIonisationModel.cc.

146{
147 const G4DynamicParticle* theAdjointPrimary = aTrack.GetDynamicParticle();
149
150 G4double adjointPrimKinEnergy = theAdjointPrimary->GetKineticEnergy();
151 G4double adjointPrimP = theAdjointPrimary->GetTotalMomentum();
152
153 if(adjointPrimKinEnergy > GetHighEnergyLimit() * 0.999)
154 {
155 return;
156 }
157
158 G4double projectileKinEnergy = 0.;
159 G4double eEnergy = 0.;
160 G4double newCS =
162 if(!isScatProjToProj)
163 { // 1/E^2 distribution
164
165 eEnergy = adjointPrimKinEnergy;
166 G4double Emax = GetSecondAdjEnergyMaxForProdToProj(adjointPrimKinEnergy);
167 G4double Emin = GetSecondAdjEnergyMinForProdToProj(adjointPrimKinEnergy);
168 if(Emin >= Emax)
169 return;
170 G4double a = 1. / Emax;
171 G4double b = 1. / Emin;
172 newCS = newCS * (b - a) / eEnergy;
173
174 projectileKinEnergy = 1. / (b - (b - a) * G4UniformRand());
175 }
176 else
177 {
178 G4double Emax =
179 GetSecondAdjEnergyMaxForScatProjToProj(adjointPrimKinEnergy);
180 G4double Emin =
182 if(Emin >= Emax)
183 return;
184 G4double diff1 = Emin - adjointPrimKinEnergy;
185 G4double diff2 = Emax - adjointPrimKinEnergy;
186
187 G4double t1 = adjointPrimKinEnergy * (1. / diff1 - 1. / diff2);
188 G4double t2 = adjointPrimKinEnergy * (1. / Emin - 1. / Emax);
189 G4double t3 = 2. * std::log(Emax / Emin);
190 G4double sum_t = t1 + t2 + t3;
191 newCS = newCS * sum_t / adjointPrimKinEnergy / adjointPrimKinEnergy;
192 G4double t = G4UniformRand() * sum_t;
193 if(t <= t1)
194 {
195 G4double q = G4UniformRand() * t1 / adjointPrimKinEnergy;
196 projectileKinEnergy = adjointPrimKinEnergy + 1. / (1. / diff1 - q);
197 }
198 else if(t <= t2)
199 {
200 G4double q = G4UniformRand() * t2 / adjointPrimKinEnergy;
201 projectileKinEnergy = 1. / (1. / Emin - q);
202 }
203 else
204 {
205 projectileKinEnergy = Emin * std::pow(Emax / Emin, G4UniformRand());
206 }
207 eEnergy = projectileKinEnergy - adjointPrimKinEnergy;
208 }
209
210 G4double diffCS_perAtom_Used = twopi_mc2_rcl2 * fMass * adjointPrimKinEnergy /
211 projectileKinEnergy / projectileKinEnergy /
212 eEnergy / eEnergy;
213
214 // Weight correction
215 // First w_corr is set to the ratio between adjoint total CS and fwd total CS
216 G4double w_corr =
218
219 w_corr *= newCS / fLastCS;
220 // Then another correction is needed due to the fact that a biaised
221 // differential CS has been used rather than the one consistent with the
222 // direct model. Here we consider the true diffCS as the one obtained by the
223 // numerical differentiation over Tcut of the direct CS
224
225 G4double diffCS =
226 DiffCrossSectionPerAtomPrimToSecond(projectileKinEnergy, eEnergy, 1, 1);
227 w_corr *= diffCS / diffCS_perAtom_Used;
228
229 if (isScatProjToProj && fTcutSecond>0.005) w_corr=1.;
230
231 G4double new_weight = aTrack.GetWeight() * w_corr;
232 fParticleChange->SetParentWeightByProcess(false);
233 fParticleChange->SetSecondaryWeightByProcess(false);
234 fParticleChange->ProposeParentWeight(new_weight);
235
236 // Kinematic:
237 // we consider a two body elastic scattering for the forward processes where
238 // the projectile knocks on an e- at rest and gives it part of its energy
240 G4double projectileTotalEnergy = projectileM0 + projectileKinEnergy;
241 G4double projectileP2 =
242 projectileTotalEnergy * projectileTotalEnergy - projectileM0 * projectileM0;
243
244 // Companion
246 if(isScatProjToProj)
247 {
248 companionM0 = fAdjEquivDirectSecondPart->GetPDGMass();
249 }
250 G4double companionTotalEnergy =
251 companionM0 + projectileKinEnergy - adjointPrimKinEnergy;
252 G4double companionP2 =
253 companionTotalEnergy * companionTotalEnergy - companionM0 * companionM0;
254
255 // Projectile momentum
256 G4double P_parallel =
257 (adjointPrimP * adjointPrimP + projectileP2 - companionP2) /
258 (2. * adjointPrimP);
259 G4double P_perp = std::sqrt(projectileP2 - P_parallel * P_parallel);
260 G4ThreeVector dir_parallel = theAdjointPrimary->GetMomentumDirection();
261 G4double phi = G4UniformRand() * twopi;
262 G4ThreeVector projectileMomentum =
263 G4ThreeVector(P_perp * std::cos(phi), P_perp * std::sin(phi), P_parallel);
264 projectileMomentum.rotateUz(dir_parallel);
265
266 if(!isScatProjToProj)
267 { // kill the primary and add a secondary
268 fParticleChange->ProposeTrackStatus(fStopAndKill);
269 fParticleChange->AddSecondary(
270 new G4DynamicParticle(fAdjEquivDirectPrimPart, projectileMomentum));
271 }
272 else
273 {
274 fParticleChange->ProposeEnergy(projectileKinEnergy);
275 fParticleChange->ProposeMomentumDirection(projectileMomentum.unit());
276 }
277}
static const G4double Emax
static const G4double Emin
static constexpr double twopi
Definition: G4SIunits.hh:56
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
static G4AdjointCSManager * GetAdjointCSManager()
G4double DiffCrossSectionPerAtomPrimToSecond(G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.) override
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetTotalMomentum() const
void AddSecondary(G4Track *aSecondary)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double GetWeight() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void ProposeTrackStatus(G4TrackStatus status)

References G4ParticleChange::AddSecondary(), G4VEmAdjointModel::DefineCurrentMaterial(), DiffCrossSectionPerAtomPrimToSecond(), Emax, Emin, G4VEmAdjointModel::fAdjEquivDirectPrimPart, G4VEmAdjointModel::fAdjEquivDirectSecondPart, G4VEmAdjointModel::fCurrentMaterial, G4VEmAdjointModel::fLastCS, fMass, fStopAndKill, G4VEmAdjointModel::fTcutSecond, G4UniformRand, G4AdjointCSManager::GetAdjointCSManager(), G4Track::GetDynamicParticle(), G4Material::GetElectronDensity(), G4VEmAdjointModel::GetHighEnergyLimit(), G4DynamicParticle::GetKineticEnergy(), G4Track::GetMaterialCutsCouple(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetPDGMass(), G4AdjointCSManager::GetPostStepWeightCorrection(), GetSecondAdjEnergyMaxForProdToProj(), GetSecondAdjEnergyMaxForScatProjToProj(), GetSecondAdjEnergyMinForProdToProj(), GetSecondAdjEnergyMinForScatProjToProj(), G4DynamicParticle::GetTotalMomentum(), G4Track::GetWeight(), G4ParticleChange::ProposeEnergy(), G4ParticleChange::ProposeMomentumDirection(), G4VParticleChange::ProposeParentWeight(), G4VParticleChange::ProposeTrackStatus(), CLHEP::Hep3Vector::rotateUz(), G4VParticleChange::SetParentWeightByProcess(), G4VParticleChange::SetSecondaryWeightByProcess(), twopi, source.hepunit::twopi_mc2_rcl2, and CLHEP::Hep3Vector::unit().

Referenced by SampleSecondaries().

◆ SampleAdjSecEnergyFromCSMatrix() [1/2]

G4double G4VEmAdjointModel::SampleAdjSecEnergyFromCSMatrix ( G4double  prim_energy,
G4bool  isScatProjToProj 
)
protectedinherited

Definition at line 518 of file G4VEmAdjointModel.cc.

520{
521 SelectCSMatrix(isScatProjToProj);
523 isScatProjToProj);
524}
G4double SampleAdjSecEnergyFromCSMatrix(size_t MatrixIndex, G4double prim_energy, G4bool isScatProjToProj)
void SelectCSMatrix(G4bool isScatProjToProj)

References G4VEmAdjointModel::fCSMatrixUsed, G4VEmAdjointModel::SampleAdjSecEnergyFromCSMatrix(), and G4VEmAdjointModel::SelectCSMatrix().

◆ SampleAdjSecEnergyFromCSMatrix() [2/2]

G4double G4VEmAdjointModel::SampleAdjSecEnergyFromCSMatrix ( size_t  MatrixIndex,
G4double  prim_energy,
G4bool  isScatProjToProj 
)
protectedinherited

Definition at line 413 of file G4VEmAdjointModel.cc.

415{
416 G4AdjointCSMatrix* theMatrix = (*fCSMatrixProdToProjBackScat)[MatrixIndex];
417 if(isScatProjToProj)
418 theMatrix = (*fCSMatrixProjToProjBackScat)[MatrixIndex];
419 std::vector<G4double>* theLogPrimEnergyVector =
420 theMatrix->GetLogPrimEnergyVector();
421
422 if(theLogPrimEnergyVector->empty())
423 {
424 G4cout << "No data are contained in the given AdjointCSMatrix!" << G4endl;
425 G4cout << "The sampling procedure will be stopped." << G4endl;
426 return 0.;
427 }
428
430 G4double aLogPrimEnergy = std::log(aPrimEnergy);
431 size_t ind = theInterpolator->FindPositionForLogVector(
432 aLogPrimEnergy, *theLogPrimEnergyVector);
433
434 G4double aLogPrimEnergy1, aLogPrimEnergy2;
435 G4double aLogCS1, aLogCS2;
436 G4double log01, log02;
437 std::vector<double>* aLogSecondEnergyVector1 = nullptr;
438 std::vector<double>* aLogSecondEnergyVector2 = nullptr;
439 std::vector<double>* aLogProbVector1 = nullptr;
440 std::vector<double>* aLogProbVector2 = nullptr;
441 std::vector<size_t>* aLogProbVectorIndex1 = nullptr;
442 std::vector<size_t>* aLogProbVectorIndex2 = nullptr;
443
444 theMatrix->GetData(ind, aLogPrimEnergy1, aLogCS1, log01,
445 aLogSecondEnergyVector1, aLogProbVector1,
446 aLogProbVectorIndex1 );
447 theMatrix->GetData(ind + 1, aLogPrimEnergy2, aLogCS2, log02,
448 aLogSecondEnergyVector2, aLogProbVector2,
449 aLogProbVectorIndex2);
450
451 if (! (aLogProbVector1 && aLogProbVector2 &&
452 aLogSecondEnergyVector1 && aLogSecondEnergyVector2)){
453 return 0.;
454 }
455
456 G4double rand_var = G4UniformRand();
457 G4double log_rand_var = std::log(rand_var);
458 G4double log_Tcut = std::log(fTcutSecond);
459 G4double Esec = 0.;
460 G4double log_dE1, log_dE2;
461 G4double log_rand_var1, log_rand_var2;
462 G4double log_E1, log_E2;
463 log_rand_var1 = log_rand_var;
464 log_rand_var2 = log_rand_var;
465
466 G4double Emin = 0.;
467 G4double Emax = 0.;
468 if(theMatrix->IsScatProjToProj())
469 { // case where Tcut plays a role
472 G4double dE = 0.;
473 if(Emin < Emax)
474 {
476 {
477 if(fSecondPartSameType && fTcutSecond > aPrimEnergy)
478 return aPrimEnergy;
479
480 log_rand_var1 = log_rand_var +
481 theInterpolator->InterpolateForLogVector(
482 log_Tcut, *aLogSecondEnergyVector1, *aLogProbVector1);
483 log_rand_var2 = log_rand_var +
484 theInterpolator->InterpolateForLogVector(
485 log_Tcut, *aLogSecondEnergyVector2, *aLogProbVector2);
486 }
487 log_dE1 = theInterpolator->Interpolate(log_rand_var1, *aLogProbVector1,
488 *aLogSecondEnergyVector1, "Lin");
489 log_dE2 = theInterpolator->Interpolate(log_rand_var2, *aLogProbVector2,
490 *aLogSecondEnergyVector2, "Lin");
491 dE = std::exp(theInterpolator->LinearInterpolation(
492 aLogPrimEnergy, aLogPrimEnergy1, aLogPrimEnergy2, log_dE1, log_dE2));
493 }
494
495 Esec = aPrimEnergy + dE;
496 Esec = std::max(Esec, Emin);
497 Esec = std::min(Esec, Emax);
498 }
499 else
500 { // Tcut condition is already full-filled
501
502 log_E1 = theInterpolator->Interpolate(log_rand_var, *aLogProbVector1,
503 *aLogSecondEnergyVector1, "Lin");
504 log_E2 = theInterpolator->Interpolate(log_rand_var, *aLogProbVector2,
505 *aLogSecondEnergyVector2, "Lin");
506
507 Esec = std::exp(theInterpolator->LinearInterpolation(
508 aLogPrimEnergy, aLogPrimEnergy1, aLogPrimEnergy2, log_E1, log_E2));
511 Esec = std::max(Esec, Emin);
512 Esec = std::min(Esec, Emax);
513 }
514 return Esec;
515}
static const G4double dE
G4bool GetData(unsigned int i, G4double &aPrimEnergy, G4double &aCS, G4double &log0, std::vector< double > *&aLogSecondEnergyVector, std::vector< double > *&aLogProbVector, std::vector< size_t > *&aLogProbVectorIndex)
std::vector< double > * GetLogPrimEnergyVector()
G4double LinearInterpolation(G4double &x, G4double &x1, G4double &x2, G4double &y1, G4double &y2)
G4double Interpolate(G4double &x, std::vector< G4double > &x_vec, std::vector< G4double > &y_vec, G4String InterPolMethod="Log")
static G4AdjointInterpolator * GetInstance()
size_t FindPositionForLogVector(G4double &x, std::vector< G4double > &x_vec)
G4double InterpolateForLogVector(G4double &x, std::vector< G4double > &x_vec, std::vector< G4double > &y_vec)

References dE, Emax, Emin, G4VEmAdjointModel::fApplyCutInRange, G4AdjointInterpolator::FindPositionForLogVector(), G4VEmAdjointModel::fSecondPartSameType, G4VEmAdjointModel::fTcutSecond, G4cout, G4endl, G4UniformRand, G4AdjointCSMatrix::GetData(), G4AdjointInterpolator::GetInstance(), G4AdjointCSMatrix::GetLogPrimEnergyVector(), G4VEmAdjointModel::GetSecondAdjEnergyMaxForProdToProj(), G4VEmAdjointModel::GetSecondAdjEnergyMaxForScatProjToProj(), G4VEmAdjointModel::GetSecondAdjEnergyMinForProdToProj(), G4VEmAdjointModel::GetSecondAdjEnergyMinForScatProjToProj(), G4AdjointInterpolator::Interpolate(), G4AdjointInterpolator::InterpolateForLogVector(), G4AdjointCSMatrix::IsScatProjToProj(), G4AdjointInterpolator::LinearInterpolation(), G4INCL::Math::max(), and G4INCL::Math::min().

Referenced by G4VEmAdjointModel::SampleAdjSecEnergyFromCSMatrix(), G4AdjointBremsstrahlungModel::SampleSecondaries(), G4AdjointComptonModel::SampleSecondaries(), G4AdjointeIonisationModel::SampleSecondaries(), SampleSecondaries(), and G4AdjointIonIonisationModel::SampleSecondaries().

◆ SampleAdjSecEnergyFromDiffCrossSectionPerAtom()

G4double G4VEmAdjointModel::SampleAdjSecEnergyFromDiffCrossSectionPerAtom ( G4double  prim_energy,
G4bool  isScatProjToProj 
)
protectedvirtualinherited

Definition at line 557 of file G4VEmAdjointModel.cc.

559{
560 // here we try to use the rejection method
561 constexpr G4int iimax = 1000;
562 G4double E = 0.;
563 G4double x, xmin, greject;
564 if(isScatProjToProj)
565 {
567 G4double Emin = prim_energy + fTcutSecond;
568 xmin = Emin / Emax;
569 G4double grejmax =
570 DiffCrossSectionPerAtomPrimToScatPrim(Emin, prim_energy, 1) * prim_energy;
571
572 G4int ii = 0;
573 do
574 {
575 // q = G4UniformRand();
576 x = 1. / (G4UniformRand() * (1. / xmin - 1.) + 1.);
577 E = x * Emax;
578 greject =
579 DiffCrossSectionPerAtomPrimToScatPrim(E, prim_energy, 1) * prim_energy;
580 ++ii;
581 if(ii >= iimax)
582 {
583 break;
584 }
585 }
586 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
587 while(greject < G4UniformRand() * grejmax);
588 }
589 else
590 {
593 xmin = Emin / Emax;
594 G4double grejmax =
596 G4int ii = 0;
597 do
598 {
599 x = std::pow(xmin, G4UniformRand());
600 E = x * Emax;
601 greject = DiffCrossSectionPerAtomPrimToSecond(E, prim_energy, 1);
602 ++ii;
603 if(ii >= iimax)
604 {
605 break;
606 }
607 }
608 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
609 while(greject < G4UniformRand() * grejmax);
610 }
611
612 return E;
613}

References G4VEmAdjointModel::DiffCrossSectionPerAtomPrimToScatPrim(), G4VEmAdjointModel::DiffCrossSectionPerAtomPrimToSecond(), Emax, Emin, G4VEmAdjointModel::fTcutSecond, G4UniformRand, G4VEmAdjointModel::GetSecondAdjEnergyMaxForProdToProj(), G4VEmAdjointModel::GetSecondAdjEnergyMaxForScatProjToProj(), and G4VEmAdjointModel::GetSecondAdjEnergyMinForProdToProj().

◆ SampleSecondaries()

void G4AdjointhIonisationModel::SampleSecondaries ( const G4Track aTrack,
G4bool  isScatProjToProj,
G4ParticleChange fParticleChange 
)
overridevirtual

Implements G4VEmAdjointModel.

Definition at line 73 of file G4AdjointhIonisationModel.cc.

76{
77 if(!fUseMatrix)
78 return RapidSampleSecondaries(aTrack, isScatProjToProj, fParticleChange);
79
80 const G4DynamicParticle* theAdjointPrimary = aTrack.GetDynamicParticle();
81
82 // Elastic inverse scattering
83 G4double adjointPrimKinEnergy = theAdjointPrimary->GetKineticEnergy();
84 G4double adjointPrimP = theAdjointPrimary->GetTotalMomentum();
85
86 if(adjointPrimKinEnergy > GetHighEnergyLimit() * 0.999)
87 {
88 return;
89 }
90
91 // Sample secondary energy
92 G4double projectileKinEnergy =
93 SampleAdjSecEnergyFromCSMatrix(adjointPrimKinEnergy, isScatProjToProj);
94 CorrectPostStepWeight(fParticleChange, aTrack.GetWeight(),
95 adjointPrimKinEnergy, projectileKinEnergy,
96 isScatProjToProj);
97 // Caution!!! this weight correction should be always applied
98
99 // Kinematic:
100 // we consider a two body elastic scattering for the forward processes where
101 // the projectile knock on an e- at rest and gives it part of its energy
103 G4double projectileTotalEnergy = projectileM0 + projectileKinEnergy;
104 G4double projectileP2 =
105 projectileTotalEnergy * projectileTotalEnergy - projectileM0 * projectileM0;
106
107 // Companion
109 if(isScatProjToProj)
110 {
111 companionM0 = fAdjEquivDirectSecondPart->GetPDGMass();
112 }
113 G4double companionTotalEnergy =
114 companionM0 + projectileKinEnergy - adjointPrimKinEnergy;
115 G4double companionP2 =
116 companionTotalEnergy * companionTotalEnergy - companionM0 * companionM0;
117
118 // Projectile momentum
119 G4double P_parallel =
120 (adjointPrimP * adjointPrimP + projectileP2 - companionP2) /
121 (2. * adjointPrimP);
122 G4double P_perp = std::sqrt(projectileP2 - P_parallel * P_parallel);
123 G4ThreeVector dir_parallel = theAdjointPrimary->GetMomentumDirection();
124 G4double phi = G4UniformRand() * twopi;
125 G4ThreeVector projectileMomentum =
126 G4ThreeVector(P_perp * std::cos(phi), P_perp * std::sin(phi), P_parallel);
127 projectileMomentum.rotateUz(dir_parallel);
128
129 if(!isScatProjToProj)
130 { // kill the primary and add a secondary
131 fParticleChange->ProposeTrackStatus(fStopAndKill);
132 fParticleChange->AddSecondary(
133 new G4DynamicParticle(fAdjEquivDirectPrimPart, projectileMomentum));
134 }
135 else
136 {
137 fParticleChange->ProposeEnergy(projectileKinEnergy);
138 fParticleChange->ProposeMomentumDirection(projectileMomentum.unit());
139 }
140}
void RapidSampleSecondaries(const G4Track &aTrack, G4bool isScatProjToProj, G4ParticleChange *fParticleChange)
virtual void CorrectPostStepWeight(G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool isScatProjToProj)

References G4ParticleChange::AddSecondary(), G4VEmAdjointModel::CorrectPostStepWeight(), G4VEmAdjointModel::fAdjEquivDirectPrimPart, G4VEmAdjointModel::fAdjEquivDirectSecondPart, fStopAndKill, G4VEmAdjointModel::fUseMatrix, G4UniformRand, G4Track::GetDynamicParticle(), G4VEmAdjointModel::GetHighEnergyLimit(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetPDGMass(), G4DynamicParticle::GetTotalMomentum(), G4Track::GetWeight(), G4ParticleChange::ProposeEnergy(), G4ParticleChange::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), RapidSampleSecondaries(), CLHEP::Hep3Vector::rotateUz(), G4VEmAdjointModel::SampleAdjSecEnergyFromCSMatrix(), twopi, and CLHEP::Hep3Vector::unit().

◆ SelectCSMatrix()

void G4VEmAdjointModel::SelectCSMatrix ( G4bool  isScatProjToProj)
protectedinherited

Definition at line 527 of file G4VEmAdjointModel.cc.

528{
529 fCSMatrixUsed = 0;
533 { // Select Material
534 std::vector<G4double>* CS_Vs_Element = &fElementCSScatProjToProj;
536 if(!isScatProjToProj)
537 {
538 CS_Vs_Element = &fElementCSProdToProj;
540 }
541 G4double SumCS = 0.;
542 size_t ind = 0;
543 for(size_t i = 0; i < CS_Vs_Element->size(); ++i)
544 {
545 SumCS += (*CS_Vs_Element)[i];
546 if(G4UniformRand() <= SumCS / fLastCS)
547 {
548 ind = i;
549 break;
550 }
551 }
553 }
554}
size_t GetIndex() const
Definition: G4Element.hh:182
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:198
size_t GetIndex() const
Definition: G4Material.hh:256
std::vector< G4double > fElementCSScatProjToProj
std::vector< G4double > fElementCSProdToProj

References G4VEmAdjointModel::fCSMatrixUsed, G4VEmAdjointModel::fCurrentMaterial, G4VEmAdjointModel::fElementCSProdToProj, G4VEmAdjointModel::fElementCSScatProjToProj, G4VEmAdjointModel::fLastAdjointCSForProdToProj, G4VEmAdjointModel::fLastAdjointCSForScatProjToProj, G4VEmAdjointModel::fLastCS, G4VEmAdjointModel::fOneMatrixForAllElements, G4VEmAdjointModel::fUseMatrixPerElement, G4UniformRand, G4Material::GetElement(), G4Element::GetIndex(), and G4Material::GetIndex().

Referenced by G4VEmAdjointModel::SampleAdjSecEnergyFromCSMatrix().

◆ SetAdditionalWeightCorrectionFactorForPostStepOutsideModel()

void G4VEmAdjointModel::SetAdditionalWeightCorrectionFactorForPostStepOutsideModel ( G4double  factor)
inlineinherited

◆ SetAdjointEquivalentOfDirectPrimaryParticleDefinition()

void G4VEmAdjointModel::SetAdjointEquivalentOfDirectPrimaryParticleDefinition ( G4ParticleDefinition aPart)
inherited

◆ SetAdjointEquivalentOfDirectSecondaryParticleDefinition()

void G4VEmAdjointModel::SetAdjointEquivalentOfDirectSecondaryParticleDefinition ( G4ParticleDefinition aPart)
inlineinherited

Definition at line 165 of file G4VEmAdjointModel.hh.

167 {
169 }

References G4VEmAdjointModel::fAdjEquivDirectSecondPart.

◆ SetApplyCutInRange()

void G4VEmAdjointModel::SetApplyCutInRange ( G4bool  aBool)
inlineinherited

◆ SetCorrectWeightForPostStepInModel()

void G4VEmAdjointModel::SetCorrectWeightForPostStepInModel ( G4bool  aBool)
inlineinherited

◆ SetCSBiasingFactor()

virtual void G4VEmAdjointModel::SetCSBiasingFactor ( G4double  aVal)
inlinevirtualinherited

Definition at line 205 of file G4VEmAdjointModel.hh.

206 {
207 fCsBiasingFactor = aVal;
208 }

References G4VEmAdjointModel::fCsBiasingFactor.

◆ SetCSMatrices()

void G4VEmAdjointModel::SetCSMatrices ( std::vector< G4AdjointCSMatrix * > *  Vec1CSMatrix,
std::vector< G4AdjointCSMatrix * > *  Vec2CSMatrix 
)
inlineinherited

Definition at line 133 of file G4VEmAdjointModel.hh.

135 {
136 fCSMatrixProdToProjBackScat = Vec1CSMatrix;
137 fCSMatrixProjToProjBackScat = Vec2CSMatrix;
138 };
std::vector< G4AdjointCSMatrix * > * fCSMatrixProdToProjBackScat
std::vector< G4AdjointCSMatrix * > * fCSMatrixProjToProjBackScat

References G4VEmAdjointModel::fCSMatrixProdToProjBackScat, and G4VEmAdjointModel::fCSMatrixProjToProjBackScat.

◆ SetHighEnergyLimit()

void G4VEmAdjointModel::SetHighEnergyLimit ( G4double  aVal)
inherited

◆ SetLowEnergyLimit()

void G4VEmAdjointModel::SetLowEnergyLimit ( G4double  aVal)
inherited

◆ SetSecondPartOfSameType()

void G4VEmAdjointModel::SetSecondPartOfSameType ( G4bool  aBool)
inlineinherited

◆ SetUseMatrix()

void G4VEmAdjointModel::SetUseMatrix ( G4bool  aBool)
inlineinherited

◆ SetUseMatrixPerElement()

void G4VEmAdjointModel::SetUseMatrixPerElement ( G4bool  aBool)
inlineinherited

◆ SetUseOnlyOneMatrixForAllElements()

void G4VEmAdjointModel::SetUseOnlyOneMatrixForAllElements ( G4bool  aBool)
inlineinherited

Field Documentation

◆ fAdjEquivDirectPrimPart

G4ParticleDefinition* G4VEmAdjointModel::fAdjEquivDirectPrimPart = nullptr
protectedinherited

◆ fAdjEquivDirectSecondPart

G4ParticleDefinition* G4VEmAdjointModel::fAdjEquivDirectSecondPart = nullptr
protectedinherited

◆ fApplyCutInRange

G4bool G4VEmAdjointModel::fApplyCutInRange = true
protectedinherited

◆ fASelectedNucleus

G4int G4VEmAdjointModel::fASelectedNucleus = 0
protectedinherited

◆ fBraggDirectEMModel

G4VEmModel* G4AdjointhIonisationModel::fBraggDirectEMModel
private

◆ fCsBiasingFactor

G4double G4VEmAdjointModel::fCsBiasingFactor = 1.
protectedinherited

◆ fCSManager

G4AdjointCSManager* G4VEmAdjointModel::fCSManager
protectedinherited

◆ fCSMatrixProdToProjBackScat

std::vector<G4AdjointCSMatrix*>* G4VEmAdjointModel::fCSMatrixProdToProjBackScat = nullptr
protectedinherited

◆ fCSMatrixProjToProjBackScat

std::vector<G4AdjointCSMatrix*>* G4VEmAdjointModel::fCSMatrixProjToProjBackScat = nullptr
protectedinherited

◆ fCSMatrixUsed

size_t G4VEmAdjointModel::fCSMatrixUsed = 0
protectedinherited

◆ fCurrentCouple

G4MaterialCutsCouple* G4VEmAdjointModel::fCurrentCouple = nullptr
protectedinherited

◆ fCurrentMaterial

G4Material* G4VEmAdjointModel::fCurrentMaterial = nullptr
protectedinherited

◆ fDirectModel

G4VEmModel* G4VEmAdjointModel::fDirectModel = nullptr
protectedinherited

◆ fDirectPrimaryPart

G4ParticleDefinition* G4VEmAdjointModel::fDirectPrimaryPart = nullptr
protectedinherited

◆ fElementCSProdToProj

std::vector<G4double> G4VEmAdjointModel::fElementCSProdToProj
protectedinherited

◆ fElementCSScatProjToProj

std::vector<G4double> G4VEmAdjointModel::fElementCSScatProjToProj
protectedinherited

◆ fFormFact

G4double G4AdjointhIonisationModel::fFormFact = 0.
private

◆ fHighEnergyLimit

G4double G4VEmAdjointModel::fHighEnergyLimit = 0.
protectedinherited

◆ fInModelWeightCorr

G4bool G4VEmAdjointModel::fInModelWeightCorr
protectedinherited

◆ fKinEnergyProdForIntegration

G4double G4VEmAdjointModel::fKinEnergyProdForIntegration = 0.
protectedinherited

◆ fKinEnergyScatProjForIntegration

G4double G4VEmAdjointModel::fKinEnergyScatProjForIntegration = 0.
protectedinherited

◆ fLastAdjointCSForProdToProj

G4double G4VEmAdjointModel::fLastAdjointCSForProdToProj = 0.
protectedinherited

◆ fLastAdjointCSForScatProjToProj

G4double G4VEmAdjointModel::fLastAdjointCSForScatProjToProj = 0.
protectedinherited

◆ fLastCS

G4double G4VEmAdjointModel::fLastCS = 0.
protectedinherited

◆ fLowEnergyLimit

G4double G4VEmAdjointModel::fLowEnergyLimit = 0.
protectedinherited

◆ fMagMoment2

G4double G4AdjointhIonisationModel::fMagMoment2 = 0.
private

◆ fMass

G4double G4AdjointhIonisationModel::fMass = 0.
private

◆ fMassRatio

G4double G4AdjointhIonisationModel::fMassRatio = 0.
private

◆ fName

const G4String G4VEmAdjointModel::fName
protectedinherited

Definition at line 252 of file G4VEmAdjointModel.hh.

Referenced by G4VEmAdjointModel::GetName().

◆ fOneMatrixForAllElements

G4bool G4VEmAdjointModel::fOneMatrixForAllElements = false
protectedinherited

◆ fOneMinusRatio2

G4double G4AdjointhIonisationModel::fOneMinusRatio2 = 0.
private

◆ fOnePlusRatio2

G4double G4AdjointhIonisationModel::fOnePlusRatio2 = 0.
private

◆ fOutsideWeightFactor

G4double G4VEmAdjointModel::fOutsideWeightFactor = 1.
protectedinherited

◆ fPreStepEnergy

G4double G4VEmAdjointModel::fPreStepEnergy = 0.
protectedinherited

◆ fSecondPartSameType

G4bool G4VEmAdjointModel::fSecondPartSameType = false
protectedinherited

◆ fSelectedMaterial

G4Material* G4VEmAdjointModel::fSelectedMaterial = nullptr
protectedinherited

◆ fSpin

G4double G4AdjointhIonisationModel::fSpin = 0.
private

◆ fTcutPrim

G4double G4VEmAdjointModel::fTcutPrim = 0.
protectedinherited

Definition at line 279 of file G4VEmAdjointModel.hh.

◆ fTcutSecond

G4double G4VEmAdjointModel::fTcutSecond = 0.
protectedinherited

◆ fUseMatrix

G4bool G4VEmAdjointModel::fUseMatrix = false
protectedinherited

◆ fUseMatrixPerElement

G4bool G4VEmAdjointModel::fUseMatrixPerElement = false
protectedinherited

◆ fZSelectedNucleus

G4int G4VEmAdjointModel::fZSelectedNucleus = 0
protectedinherited

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