Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes
G4VEmAdjointModel Class Referenceabstract

#include <G4VEmAdjointModel.hh>

Inheritance diagram for G4VEmAdjointModel:
G4AdjointBremsstrahlungModel G4AdjointComptonModel G4AdjointIonIonisationModel G4AdjointPhotoElectricModel G4AdjointeIonisationModel G4AdjointhIonisationModel

Public Member Functions

virtual G4double AdjointCrossSection (const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool isScatProjToProj)
 
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.)
 
virtual G4double DiffCrossSectionPerAtomPrimToSecond (G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)
 
virtual G4double DiffCrossSectionPerVolumePrimToScatPrim (const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyScatProj)
 
virtual G4double DiffCrossSectionPerVolumePrimToSecond (const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
 
 G4VEmAdjointModel (const G4String &nam)
 
 G4VEmAdjointModel (G4VEmAdjointModel &)=delete
 
G4ParticleDefinitionGetAdjointEquivalentOfDirectPrimaryParticleDefinition ()
 
G4ParticleDefinitionGetAdjointEquivalentOfDirectSecondaryParticleDefinition ()
 
G4bool GetApplyCutInRange ()
 
G4double GetHighEnergyLimit ()
 
G4double GetLowEnergyLimit ()
 
G4String GetName ()
 
virtual G4double GetSecondAdjEnergyMaxForProdToProj (G4double primAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMaxForScatProjToProj (G4double primAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMinForProdToProj (G4double primAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMinForScatProjToProj (G4double primAdjEnergy, G4double tcut=0.)
 
G4bool GetSecondPartOfSameType ()
 
G4bool GetUseMatrix ()
 
G4bool GetUseMatrixPerElement ()
 
G4bool GetUseOnlyOneMatrixForAllElements ()
 
G4VEmAdjointModeloperator= (const G4VEmAdjointModel &right)=delete
 
virtual void SampleSecondaries (const G4Track &aTrack, G4bool isScatProjToProj, G4ParticleChange *fParticleChange)=0
 
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)
 
virtual ~G4VEmAdjointModel ()
 

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
 

Detailed Description

Definition at line 50 of file G4VEmAdjointModel.hh.

Constructor & Destructor Documentation

◆ G4VEmAdjointModel() [1/2]

G4VEmAdjointModel::G4VEmAdjointModel ( const G4String nam)
explicit

◆ ~G4VEmAdjointModel()

G4VEmAdjointModel::~G4VEmAdjointModel ( )
virtual

Definition at line 50 of file G4VEmAdjointModel.cc.

51{
54
57}
std::vector< G4AdjointCSMatrix * > * fCSMatrixProdToProjBackScat
std::vector< G4AdjointCSMatrix * > * fCSMatrixProjToProjBackScat

References fCSMatrixProdToProjBackScat, and fCSMatrixProjToProjBackScat.

◆ G4VEmAdjointModel() [2/2]

G4VEmAdjointModel::G4VEmAdjointModel ( G4VEmAdjointModel )
delete

Member Function Documentation

◆ AdjointCrossSection()

G4double G4VEmAdjointModel::AdjointCrossSection ( const G4MaterialCutsCouple aCouple,
G4double  primEnergy,
G4bool  isScatProjToProj 
)
virtual

Reimplemented in G4AdjointBremsstrahlungModel, G4AdjointComptonModel, G4AdjointhIonisationModel, and G4AdjointPhotoElectricModel.

Definition at line 60 of file G4VEmAdjointModel.cc.

63{
64 DefineCurrentMaterial(aCouple);
65 fPreStepEnergy = primEnergy;
66
67 std::vector<G4double>* CS_Vs_Element = &fElementCSProdToProj;
68 if(isScatProjToProj)
69 CS_Vs_Element = &fElementCSScatProjToProj;
70 fLastCS =
72 fTcutSecond, isScatProjToProj, *CS_Vs_Element);
73 if(isScatProjToProj)
75 else
77
78 return fLastCS;
79}
G4double ComputeAdjointCS(G4Material *aMaterial, G4VEmAdjointModel *aModel, G4double PrimEnergy, G4double Tcut, G4bool isScatProjToProj, std::vector< G4double > &AdjointCS_for_each_element)
G4double fLastAdjointCSForScatProjToProj
G4Material * fCurrentMaterial
std::vector< G4double > fElementCSScatProjToProj
void DefineCurrentMaterial(const G4MaterialCutsCouple *couple)
std::vector< G4double > fElementCSProdToProj
G4double fLastAdjointCSForProdToProj

References G4AdjointCSManager::ComputeAdjointCS(), DefineCurrentMaterial(), fCSManager, fCurrentMaterial, fElementCSProdToProj, fElementCSScatProjToProj, fLastAdjointCSForProdToProj, fLastAdjointCSForScatProjToProj, fLastCS, fPreStepEnergy, and fTcutSecond.

Referenced by G4AdjointBremsstrahlungModel::AdjointCrossSection(), G4AdjointComptonModel::AdjointCrossSection(), G4AdjointhIonisationModel::AdjointCrossSection(), G4AdjointCSManager::ComputeAdjointCS(), CorrectPostStepWeight(), G4VAdjointReverseReaction::GetMeanFreePath(), and G4AdjointForcedInteractionForGamma::PostStepDoIt().

◆ ComputeAdjointCrossSectionVectorPerAtomForScatProj()

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

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}
double G4double
Definition: G4Types.hh:83
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)
G4VEmAdjointModel(const G4String &nam)
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, DiffCrossSectionFunction2(), fASelectedNucleus, fKinEnergyScatProjForIntegration, fZSelectedNucleus, G4lrint(), G4VEmAdjointModel(), GetLowEnergyLimit(), GetSecondAdjEnergyMaxForScatProjToProj(), 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 
)

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, DiffCrossSectionFunction1(), fASelectedNucleus, fKinEnergyProdForIntegration, fZSelectedNucleus, G4lrint(), G4VEmAdjointModel(), GetSecondAdjEnergyMaxForProdToProj(), 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 
)

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 DiffCrossSectionFunction2(), fKinEnergyScatProjForIntegration, fSelectedMaterial, G4VEmAdjointModel(), GetLowEnergyLimit(), GetSecondAdjEnergyMaxForScatProjToProj(), 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 
)

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 DiffCrossSectionFunction1(), fKinEnergyProdForIntegration, fSelectedMaterial, G4VEmAdjointModel(), GetSecondAdjEnergyMaxForProdToProj(), GetSecondAdjEnergyMinForProdToProj(), and G4INCL::Math::min().

Referenced by G4AdjointCSManager::BuildCrossSectionsModelAndMaterial().

◆ CorrectPostStepWeight()

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

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()
G4MaterialCutsCouple * fCurrentCouple
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool isScatProjToProj)
void SetSecondaryWeightByProcess(G4bool)
void SetParentWeightByProcess(G4bool)
void ProposeParentWeight(G4double finalWeight)

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

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

◆ DefineCurrentMaterial()

void G4VEmAdjointModel::DefineCurrentMaterial ( const G4MaterialCutsCouple couple)

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 G4AdjointElectron * AdjointElectron()
static G4AdjointGamma * AdjointGamma()
static G4AdjointPositron * AdjointPositron()
const G4Material * GetMaterial() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleDefinition * fAdjEquivDirectSecondPart

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

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

◆ DefineDirectEMModel()

void G4VEmAdjointModel::DefineDirectEMModel ( G4VEmModel aModel)
inline

Definition at line 160 of file G4VEmAdjointModel.hh.

160{ fDirectModel = aModel; }
G4VEmModel * fDirectModel

References fDirectModel.

◆ DiffCrossSectionFunction1()

G4double G4VEmAdjointModel::DiffCrossSectionFunction1 ( G4double  kinEnergyProj)
protected

Definition at line 155 of file G4VEmAdjointModel.cc.

156{
157 G4double bias_factor =
159
161 {
165 bias_factor;
166 }
167 else
168 {
171 bias_factor;
172 }
173}
virtual G4double DiffCrossSectionPerVolumePrimToSecond(const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
virtual G4double DiffCrossSectionPerAtomPrimToSecond(G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)

References DiffCrossSectionPerAtomPrimToSecond(), DiffCrossSectionPerVolumePrimToSecond(), fASelectedNucleus, fCsBiasingFactor, fKinEnergyProdForIntegration, fSelectedMaterial, fUseMatrixPerElement, and fZSelectedNucleus.

Referenced by ComputeAdjointCrossSectionVectorPerAtomForSecond(), and ComputeAdjointCrossSectionVectorPerVolumeForSecond().

◆ DiffCrossSectionFunction2()

G4double G4VEmAdjointModel::DiffCrossSectionFunction2 ( G4double  kinEnergyProj)
protected

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 DiffCrossSectionPerAtomPrimToScatPrim(), DiffCrossSectionPerVolumePrimToScatPrim(), fASelectedNucleus, fCsBiasingFactor, fKinEnergyScatProjForIntegration, fSelectedMaterial, fUseMatrixPerElement, and fZSelectedNucleus.

Referenced by ComputeAdjointCrossSectionVectorPerAtomForScatProj(), and ComputeAdjointCrossSectionVectorPerVolumeForScatProj().

◆ DiffCrossSectionPerAtomPrimToScatPrim()

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

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, DiffCrossSectionPerAtomPrimToSecond(), and Z.

Referenced by DiffCrossSectionFunction2(), and SampleAdjSecEnergyFromDiffCrossSectionPerAtom().

◆ DiffCrossSectionPerAtomPrimToSecond()

G4double G4VEmAdjointModel::DiffCrossSectionPerAtomPrimToSecond ( G4double  kinEnergyProj,
G4double  kinEnergyProd,
G4double  Z,
G4double  A = 0. 
)
virtual

Reimplemented in G4AdjointComptonModel, G4AdjointeIonisationModel, G4AdjointhIonisationModel, and G4AdjointIonIonisationModel.

Definition at line 82 of file G4VEmAdjointModel.cc.

84{
85 G4double dSigmadEprod = 0.;
86 G4double Emax_proj = GetSecondAdjEnergyMaxForProdToProj(kinEnergyProd);
87 G4double Emin_proj = GetSecondAdjEnergyMinForProdToProj(kinEnergyProd);
88
89 // the produced particle should have a kinetic energy less than the projectile
90 if(kinEnergyProj > Emin_proj && kinEnergyProj <= Emax_proj)
91 {
92 G4double E1 = kinEnergyProd;
93 G4double E2 = kinEnergyProd * 1.000001;
95 fDirectPrimaryPart, kinEnergyProj, Z, A, E1, 1.e20);
97 fDirectPrimaryPart, kinEnergyProj, Z, A, E2, 1.e20);
98
99 dSigmadEprod = (sigma1 - sigma2) / (E2 - E1);
100 }
101 return dSigmadEprod;
102}
G4ParticleDefinition * fDirectPrimaryPart
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(), fDirectModel, fDirectPrimaryPart, GetSecondAdjEnergyMaxForProdToProj(), GetSecondAdjEnergyMinForProdToProj(), and Z.

Referenced by DiffCrossSectionFunction1(), DiffCrossSectionPerAtomPrimToScatPrim(), and SampleAdjSecEnergyFromDiffCrossSectionPerAtom().

◆ DiffCrossSectionPerVolumePrimToScatPrim()

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

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 DiffCrossSectionPerVolumePrimToSecond().

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

◆ DiffCrossSectionPerVolumePrimToSecond()

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

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(), fDirectModel, fDirectPrimaryPart, GetSecondAdjEnergyMaxForProdToProj(), and GetSecondAdjEnergyMinForProdToProj().

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

◆ GetAdjointEquivalentOfDirectPrimaryParticleDefinition()

G4ParticleDefinition * G4VEmAdjointModel::GetAdjointEquivalentOfDirectPrimaryParticleDefinition ( )
inline

Definition at line 141 of file G4VEmAdjointModel.hh.

142 {
144 }
G4ParticleDefinition * fAdjEquivDirectPrimPart

References fAdjEquivDirectPrimPart.

Referenced by G4AdjointCSManager::ComputeTotalAdjointCS().

◆ GetAdjointEquivalentOfDirectSecondaryParticleDefinition()

G4ParticleDefinition * G4VEmAdjointModel::GetAdjointEquivalentOfDirectSecondaryParticleDefinition ( )
inline

Definition at line 147 of file G4VEmAdjointModel.hh.

148 {
150 }

References fAdjEquivDirectSecondPart.

Referenced by G4AdjointCSManager::ComputeTotalAdjointCS().

◆ GetApplyCutInRange()

G4bool G4VEmAdjointModel::GetApplyCutInRange ( )
inline

◆ GetHighEnergyLimit()

G4double G4VEmAdjointModel::GetHighEnergyLimit ( )
inline

◆ GetLowEnergyLimit()

G4double G4VEmAdjointModel::GetLowEnergyLimit ( )
inline

◆ GetName()

G4String G4VEmAdjointModel::GetName ( )
inline

Definition at line 203 of file G4VEmAdjointModel.hh.

203{ return fName; }

References fName.

◆ GetSecondAdjEnergyMaxForProdToProj()

G4double G4VEmAdjointModel::GetSecondAdjEnergyMaxForProdToProj ( G4double  primAdjEnergy)
virtual

◆ GetSecondAdjEnergyMaxForScatProjToProj()

G4double G4VEmAdjointModel::GetSecondAdjEnergyMaxForScatProjToProj ( G4double  primAdjEnergy)
virtual

◆ GetSecondAdjEnergyMinForProdToProj()

G4double G4VEmAdjointModel::GetSecondAdjEnergyMinForProdToProj ( G4double  primAdjEnergy)
virtual

◆ GetSecondAdjEnergyMinForScatProjToProj()

G4double G4VEmAdjointModel::GetSecondAdjEnergyMinForScatProjToProj ( G4double  primAdjEnergy,
G4double  tcut = 0. 
)
virtual

◆ GetSecondPartOfSameType()

G4bool G4VEmAdjointModel::GetSecondPartOfSameType ( )
inline

◆ GetUseMatrix()

G4bool G4VEmAdjointModel::GetUseMatrix ( )
inline

Definition at line 192 of file G4VEmAdjointModel.hh.

192{ return fUseMatrix; }

References fUseMatrix.

Referenced by G4AdjointCSManager::ComputeAdjointCS().

◆ GetUseMatrixPerElement()

G4bool G4VEmAdjointModel::GetUseMatrixPerElement ( )
inline

Definition at line 194 of file G4VEmAdjointModel.hh.

194{ return fUseMatrixPerElement; }

References fUseMatrixPerElement.

Referenced by G4AdjointCSManager::ComputeAdjointCS().

◆ GetUseOnlyOneMatrixForAllElements()

G4bool G4VEmAdjointModel::GetUseOnlyOneMatrixForAllElements ( )
inline

Definition at line 196 of file G4VEmAdjointModel.hh.

197 {
199 }

References fOneMatrixForAllElements.

Referenced by G4AdjointCSManager::ComputeAdjointCS().

◆ operator=()

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

◆ SampleAdjSecEnergyFromCSMatrix() [1/2]

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

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 fCSMatrixUsed, SampleAdjSecEnergyFromCSMatrix(), and SelectCSMatrix().

◆ SampleAdjSecEnergyFromCSMatrix() [2/2]

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

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 Emax
static const G4double dE
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
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, fApplyCutInRange, G4AdjointInterpolator::FindPositionForLogVector(), fSecondPartSameType, fTcutSecond, G4cout, G4endl, G4UniformRand, G4AdjointCSMatrix::GetData(), G4AdjointInterpolator::GetInstance(), G4AdjointCSMatrix::GetLogPrimEnergyVector(), GetSecondAdjEnergyMaxForProdToProj(), GetSecondAdjEnergyMaxForScatProjToProj(), GetSecondAdjEnergyMinForProdToProj(), GetSecondAdjEnergyMinForScatProjToProj(), G4AdjointInterpolator::Interpolate(), G4AdjointInterpolator::InterpolateForLogVector(), G4AdjointCSMatrix::IsScatProjToProj(), G4AdjointInterpolator::LinearInterpolation(), G4INCL::Math::max(), and G4INCL::Math::min().

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

◆ SampleAdjSecEnergyFromDiffCrossSectionPerAtom()

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

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 DiffCrossSectionPerAtomPrimToScatPrim(), DiffCrossSectionPerAtomPrimToSecond(), Emax, Emin, fTcutSecond, G4UniformRand, GetSecondAdjEnergyMaxForProdToProj(), GetSecondAdjEnergyMaxForScatProjToProj(), and GetSecondAdjEnergyMinForProdToProj().

◆ SampleSecondaries()

virtual void G4VEmAdjointModel::SampleSecondaries ( const G4Track aTrack,
G4bool  isScatProjToProj,
G4ParticleChange fParticleChange 
)
pure virtual

◆ SelectCSMatrix()

void G4VEmAdjointModel::SelectCSMatrix ( G4bool  isScatProjToProj)
protected

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

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

Referenced by SampleAdjSecEnergyFromCSMatrix().

◆ SetAdditionalWeightCorrectionFactorForPostStepOutsideModel()

void G4VEmAdjointModel::SetAdditionalWeightCorrectionFactorForPostStepOutsideModel ( G4double  factor)
inline

Definition at line 215 of file G4VEmAdjointModel.hh.

217 {
218 fOutsideWeightFactor = factor;
219 }

References fOutsideWeightFactor.

Referenced by G4AdjointForcedInteractionForGamma::PostStepDoIt().

◆ SetAdjointEquivalentOfDirectPrimaryParticleDefinition()

void G4VEmAdjointModel::SetAdjointEquivalentOfDirectPrimaryParticleDefinition ( G4ParticleDefinition aPart)

Definition at line 729 of file G4VEmAdjointModel.cc.

731{
735 else if(fAdjEquivDirectPrimPart->GetParticleName() == "adj_gamma")
737}
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
const G4String & GetParticleName() const

References G4Electron::Electron(), fAdjEquivDirectPrimPart, fDirectPrimaryPart, G4Gamma::Gamma(), and G4ParticleDefinition::GetParticleName().

◆ SetAdjointEquivalentOfDirectSecondaryParticleDefinition()

void G4VEmAdjointModel::SetAdjointEquivalentOfDirectSecondaryParticleDefinition ( G4ParticleDefinition aPart)
inline

Definition at line 165 of file G4VEmAdjointModel.hh.

167 {
169 }

References fAdjEquivDirectSecondPart.

◆ SetApplyCutInRange()

void G4VEmAdjointModel::SetApplyCutInRange ( G4bool  aBool)
inline

◆ SetCorrectWeightForPostStepInModel()

void G4VEmAdjointModel::SetCorrectWeightForPostStepInModel ( G4bool  aBool)
inline

Definition at line 210 of file G4VEmAdjointModel.hh.

211 {
212 fInModelWeightCorr = aBool;
213 }

References fInModelWeightCorr.

Referenced by G4AdjointForcedInteractionForGamma::PostStepDoIt().

◆ SetCSBiasingFactor()

virtual void G4VEmAdjointModel::SetCSBiasingFactor ( G4double  aVal)
inlinevirtual

Definition at line 205 of file G4VEmAdjointModel.hh.

206 {
207 fCsBiasingFactor = aVal;
208 }

References fCsBiasingFactor.

◆ SetCSMatrices()

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

Definition at line 133 of file G4VEmAdjointModel.hh.

135 {
136 fCSMatrixProdToProjBackScat = Vec1CSMatrix;
137 fCSMatrixProjToProjBackScat = Vec2CSMatrix;
138 };

References fCSMatrixProdToProjBackScat, and fCSMatrixProjToProjBackScat.

◆ SetHighEnergyLimit()

void G4VEmAdjointModel::SetHighEnergyLimit ( G4double  aVal)

Definition at line 713 of file G4VEmAdjointModel.cc.

714{
715 fHighEnergyLimit = aVal;
716 if(fDirectModel)
718}
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767

References fDirectModel, fHighEnergyLimit, and G4VEmModel::SetHighEnergyLimit().

◆ SetLowEnergyLimit()

void G4VEmAdjointModel::SetLowEnergyLimit ( G4double  aVal)

Definition at line 721 of file G4VEmAdjointModel.cc.

722{
723 fLowEnergyLimit = aVal;
724 if(fDirectModel)
726}
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774

References fDirectModel, fLowEnergyLimit, and G4VEmModel::SetLowEnergyLimit().

◆ SetSecondPartOfSameType()

void G4VEmAdjointModel::SetSecondPartOfSameType ( G4bool  aBool)
inline

◆ SetUseMatrix()

void G4VEmAdjointModel::SetUseMatrix ( G4bool  aBool)
inline

◆ SetUseMatrixPerElement()

void G4VEmAdjointModel::SetUseMatrixPerElement ( G4bool  aBool)
inline

◆ SetUseOnlyOneMatrixForAllElements()

void G4VEmAdjointModel::SetUseOnlyOneMatrixForAllElements ( G4bool  aBool)
inline

Definition at line 185 of file G4VEmAdjointModel.hh.

186 {
188 }

References fOneMatrixForAllElements.

Referenced by G4AdjointComptonModel::G4AdjointComptonModel().

Field Documentation

◆ fAdjEquivDirectPrimPart

G4ParticleDefinition* G4VEmAdjointModel::fAdjEquivDirectPrimPart = nullptr
protected

◆ fAdjEquivDirectSecondPart

G4ParticleDefinition* G4VEmAdjointModel::fAdjEquivDirectSecondPart = nullptr
protected

◆ fApplyCutInRange

G4bool G4VEmAdjointModel::fApplyCutInRange = true
protected

◆ fASelectedNucleus

G4int G4VEmAdjointModel::fASelectedNucleus = 0
protected

◆ fCsBiasingFactor

G4double G4VEmAdjointModel::fCsBiasingFactor = 1.
protected

◆ fCSManager

G4AdjointCSManager* G4VEmAdjointModel::fCSManager
protected

◆ fCSMatrixProdToProjBackScat

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

Definition at line 264 of file G4VEmAdjointModel.hh.

Referenced by SetCSMatrices(), and ~G4VEmAdjointModel().

◆ fCSMatrixProjToProjBackScat

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

Definition at line 265 of file G4VEmAdjointModel.hh.

Referenced by SetCSMatrices(), and ~G4VEmAdjointModel().

◆ fCSMatrixUsed

size_t G4VEmAdjointModel::fCSMatrixUsed = 0
protected

Definition at line 298 of file G4VEmAdjointModel.hh.

Referenced by SampleAdjSecEnergyFromCSMatrix(), and SelectCSMatrix().

◆ fCurrentCouple

G4MaterialCutsCouple* G4VEmAdjointModel::fCurrentCouple = nullptr
protected

◆ fCurrentMaterial

G4Material* G4VEmAdjointModel::fCurrentMaterial = nullptr
protected

◆ fDirectModel

G4VEmModel* G4VEmAdjointModel::fDirectModel = nullptr
protected

◆ fDirectPrimaryPart

G4ParticleDefinition* G4VEmAdjointModel::fDirectPrimaryPart = nullptr
protected

◆ fElementCSProdToProj

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

Definition at line 268 of file G4VEmAdjointModel.hh.

Referenced by AdjointCrossSection(), and SelectCSMatrix().

◆ fElementCSScatProjToProj

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

Definition at line 267 of file G4VEmAdjointModel.hh.

Referenced by AdjointCrossSection(), and SelectCSMatrix().

◆ fHighEnergyLimit

G4double G4VEmAdjointModel::fHighEnergyLimit = 0.
protected

◆ fInModelWeightCorr

G4bool G4VEmAdjointModel::fInModelWeightCorr
protected

◆ fKinEnergyProdForIntegration

G4double G4VEmAdjointModel::fKinEnergyProdForIntegration = 0.
protected

◆ fKinEnergyScatProjForIntegration

G4double G4VEmAdjointModel::fKinEnergyScatProjForIntegration = 0.
protected

◆ fLastAdjointCSForProdToProj

G4double G4VEmAdjointModel::fLastAdjointCSForProdToProj = 0.
protected

◆ fLastAdjointCSForScatProjToProj

G4double G4VEmAdjointModel::fLastAdjointCSForScatProjToProj = 0.
protected

◆ fLastCS

G4double G4VEmAdjointModel::fLastCS = 0.
protected

◆ fLowEnergyLimit

G4double G4VEmAdjointModel::fLowEnergyLimit = 0.
protected

Definition at line 284 of file G4VEmAdjointModel.hh.

Referenced by GetLowEnergyLimit(), and SetLowEnergyLimit().

◆ fName

const G4String G4VEmAdjointModel::fName
protected

Definition at line 252 of file G4VEmAdjointModel.hh.

Referenced by GetName().

◆ fOneMatrixForAllElements

G4bool G4VEmAdjointModel::fOneMatrixForAllElements = false
protected

◆ fOutsideWeightFactor

G4double G4VEmAdjointModel::fOutsideWeightFactor = 1.
protected

◆ fPreStepEnergy

G4double G4VEmAdjointModel::fPreStepEnergy = 0.
protected

Definition at line 277 of file G4VEmAdjointModel.hh.

Referenced by AdjointCrossSection(), and CorrectPostStepWeight().

◆ fSecondPartSameType

G4bool G4VEmAdjointModel::fSecondPartSameType = false
protected

◆ fSelectedMaterial

G4Material* G4VEmAdjointModel::fSelectedMaterial = nullptr
protected

◆ fTcutPrim

G4double G4VEmAdjointModel::fTcutPrim = 0.
protected

Definition at line 279 of file G4VEmAdjointModel.hh.

◆ fTcutSecond

G4double G4VEmAdjointModel::fTcutSecond = 0.
protected

◆ fUseMatrix

G4bool G4VEmAdjointModel::fUseMatrix = false
protected

◆ fUseMatrixPerElement

G4bool G4VEmAdjointModel::fUseMatrixPerElement = false
protected

◆ fZSelectedNucleus

G4int G4VEmAdjointModel::fZSelectedNucleus = 0
protected

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