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

#include <G4AdjointhIonisationModel.hh>

+ Inheritance diagram for G4AdjointhIonisationModel:

Public Member Functions

 G4AdjointhIonisationModel (G4ParticleDefinition *projectileDefinition)
 
virtual ~G4AdjointhIonisationModel ()
 
virtual void SampleSecondaries (const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)
 
void RapidSampleSecondaries (const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)
 
virtual G4double DiffCrossSectionPerAtomPrimToSecond (G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)
 
virtual G4double AdjointCrossSection (const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
 
virtual G4double GetSecondAdjEnergyMaxForScatProjToProjCase (G4double PrimAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMinForScatProjToProjCase (G4double PrimAdjEnergy, G4double Tcut=0)
 
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase (G4double PrimAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMinForProdToProjCase (G4double PrimAdjEnergy)
 
- Public Member Functions inherited from G4VEmAdjointModel
 G4VEmAdjointModel (const G4String &nam)
 
virtual ~G4VEmAdjointModel ()
 
virtual void SampleSecondaries (const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)=0
 
virtual G4double AdjointCrossSection (const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
 
virtual G4double GetAdjointCrossSection (const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)
 
virtual G4double DiffCrossSectionPerAtomPrimToSecond (G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)
 
virtual G4double DiffCrossSectionPerAtomPrimToScatPrim (G4double kinEnergyProj, G4double kinEnergyScatProj, G4double Z, G4double A=0.)
 
virtual G4double DiffCrossSectionPerVolumePrimToSecond (const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyProd)
 
virtual G4double DiffCrossSectionPerVolumePrimToScatPrim (const G4Material *aMaterial, G4double kinEnergyProj, G4double kinEnergyScatProj)
 
virtual G4double GetSecondAdjEnergyMaxForScatProjToProjCase (G4double PrimAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMinForScatProjToProjCase (G4double PrimAdjEnergy, G4double Tcut=0)
 
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase (G4double PrimAdjEnergy)
 
virtual G4double GetSecondAdjEnergyMinForProdToProjCase (G4double PrimAdjEnergy)
 
void DefineCurrentMaterial (const G4MaterialCutsCouple *couple)
 
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerAtomForSecond (G4double kinEnergyProd, G4double Z, G4double A=0., G4int nbin_pro_decade=10)
 
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerAtomForScatProj (G4double kinEnergyProd, G4double Z, G4double A=0., G4int nbin_pro_decade=10)
 
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerVolumeForSecond (G4Material *aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade=10)
 
std::vector< std::vector< double > * > ComputeAdjointCrossSectionVectorPerVolumeForScatProj (G4Material *aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade=10)
 
void SetCSMatrices (std::vector< G4AdjointCSMatrix * > *Vec1CSMatrix, std::vector< G4AdjointCSMatrix * > *Vec2CSMatrix)
 
G4ParticleDefinitionGetAdjointEquivalentOfDirectPrimaryParticleDefinition ()
 
G4ParticleDefinitionGetAdjointEquivalentOfDirectSecondaryParticleDefinition ()
 
G4double GetHighEnergyLimit ()
 
G4double GetLowEnergyLimit ()
 
void SetHighEnergyLimit (G4double aVal)
 
void SetLowEnergyLimit (G4double aVal)
 
void DefineDirectEMModel (G4VEmModel *aModel)
 
void SetAdjointEquivalentOfDirectPrimaryParticleDefinition (G4ParticleDefinition *aPart)
 
void SetAdjointEquivalentOfDirectSecondaryParticleDefinition (G4ParticleDefinition *aPart)
 
void SetSecondPartOfSameType (G4bool aBool)
 
G4bool GetSecondPartOfSameType ()
 
void SetUseMatrix (G4bool aBool)
 
void SetUseMatrixPerElement (G4bool aBool)
 
void SetUseOnlyOneMatrixForAllElements (G4bool aBool)
 
void SetApplyCutInRange (G4bool aBool)
 
G4bool GetUseMatrix ()
 
G4bool GetUseMatrixPerElement ()
 
G4bool GetUseOnlyOneMatrixForAllElements ()
 
G4bool GetApplyCutInRange ()
 
G4String GetName ()
 
virtual void SetCSBiasingFactor (G4double aVal)
 
void SetCorrectWeightForPostStepInModel (G4bool aBool)
 
void SetAdditionalWeightCorrectionFactorForPostStepOutsideModel (G4double factor)
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmAdjointModel
G4double DiffCrossSectionFunction1 (G4double kinEnergyProj)
 
G4double DiffCrossSectionFunction2 (G4double kinEnergyProj)
 
G4double DiffCrossSectionPerVolumeFunctionForIntegrationOverEkinProj (G4double EkinProd)
 
G4double SampleAdjSecEnergyFromCSMatrix (size_t MatrixIndex, G4double prim_energy, G4bool IsScatProjToProjCase)
 
G4double SampleAdjSecEnergyFromCSMatrix (G4double prim_energy, G4bool IsScatProjToProjCase)
 
void SelectCSMatrix (G4bool IsScatProjToProjCase)
 
virtual G4double SampleAdjSecEnergyFromDiffCrossSectionPerAtom (G4double prim_energy, G4bool IsScatProjToProjCase)
 
virtual void CorrectPostStepWeight (G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool IsScatProjToProjCase)
 
- Protected Attributes inherited from G4VEmAdjointModel
G4VEmModeltheDirectEMModel
 
G4VParticleChangepParticleChange
 
const G4String name
 
G4int ASelectedNucleus
 
G4int ZSelectedNucleus
 
G4MaterialSelectedMaterial
 
G4double kinEnergyProdForIntegration
 
G4double kinEnergyScatProjForIntegration
 
G4double kinEnergyProjForIntegration
 
std::vector< G4AdjointCSMatrix * > * pOnCSMatrixForProdToProjBackwardScattering
 
std::vector< G4AdjointCSMatrix * > * pOnCSMatrixForScatProjToProjBackwardScattering
 
std::vector< G4doubleCS_Vs_ElementForScatProjToProjCase
 
std::vector< G4doubleCS_Vs_ElementForProdToProjCase
 
G4double lastCS
 
G4double lastAdjointCSForScatProjToProjCase
 
G4double lastAdjointCSForProdToProjCase
 
G4ParticleDefinitiontheAdjEquivOfDirectPrimPartDef
 
G4ParticleDefinitiontheAdjEquivOfDirectSecondPartDef
 
G4ParticleDefinitiontheDirectPrimaryPartDef
 
G4bool second_part_of_same_type
 
G4double preStepEnergy
 
G4MaterialcurrentMaterial
 
G4MaterialCutsCouplecurrentCouple
 
size_t currentMaterialIndex
 
size_t currentCoupleIndex
 
G4double currentTcutForDirectPrim
 
G4double currentTcutForDirectSecond
 
G4bool ApplyCutInRange
 
G4double mass_ratio_product
 
G4double mass_ratio_projectile
 
G4double HighEnergyLimit
 
G4double LowEnergyLimit
 
G4double CS_biasing_factor
 
G4bool UseMatrix
 
G4bool UseMatrixPerElement
 
G4bool UseOnlyOneMatrixForAllElements
 
size_t indexOfUsedCrossSectionMatrix
 
size_t model_index
 
G4bool correct_weight_for_post_step_in_model
 
G4double additional_weight_correction_factor_for_post_step_outside_model
 

Detailed Description

Definition at line 71 of file G4AdjointhIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4AdjointhIonisationModel()

G4AdjointhIonisationModel::G4AdjointhIonisationModel ( G4ParticleDefinition projectileDefinition)

Definition at line 45 of file G4AdjointhIonisationModel.cc.

45 :
46 G4VEmAdjointModel("Adjoint_hIonisation")
47{
48
49
50
51 UseMatrix =true;
53 ApplyCutInRange = true;
57
58 //The direct EM Modfel is taken has BetheBloch it is only used for the computation
59 // of the differential cross section.
60 //The Bragg model could be used as an alternative as it offers the same differential cross section
61
62 theDirectEMModel = new G4BetheBlochModel(projectileDefinition);
63 theBraggDirectEMModel = new G4BraggModel(projectileDefinition);
65
66 theDirectPrimaryPartDef = projectileDefinition;
68 if (projectileDefinition == G4Proton::Proton()) {
70
71 }
72
73 DefineProjectileProperty();
74}
static G4AdjointElectron * AdjointElectron()
static G4AdjointProton * AdjointProton()
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4VEmModel * theDirectEMModel
G4ParticleDefinition * theDirectPrimaryPartDef
G4bool UseOnlyOneMatrixForAllElements
G4ParticleDefinition * theAdjEquivOfDirectSecondPartDef
G4ParticleDefinition * theAdjEquivOfDirectPrimPartDef

◆ ~G4AdjointhIonisationModel()

G4AdjointhIonisationModel::~G4AdjointhIonisationModel ( )
virtual

Definition at line 77 of file G4AdjointhIonisationModel.cc.

78{;}

Member Function Documentation

◆ AdjointCrossSection()

G4double G4AdjointhIonisationModel::AdjointCrossSection ( const G4MaterialCutsCouple aCouple,
G4double  primEnergy,
G4bool  IsScatProjToProjCase 
)
virtual

Reimplemented from G4VEmAdjointModel.

Definition at line 437 of file G4AdjointhIonisationModel.cc.

440{
441 if (UseMatrix) return G4VEmAdjointModel::AdjointCrossSection(aCouple,primEnergy,IsScatProjToProjCase);
442 DefineCurrentMaterial(aCouple);
443
444
445 G4double Cross=currentMaterial->GetElectronDensity()*twopi_mc2_rcl2*mass;
446
447 if (!IsScatProjToProjCase ){
448 G4double Emax_proj = GetSecondAdjEnergyMaxForProdToProjCase(primEnergy);
449 G4double Emin_proj = GetSecondAdjEnergyMinForProdToProjCase(primEnergy);
450 if (Emax_proj>Emin_proj && primEnergy > currentTcutForDirectSecond) {
451 Cross*=(1./Emin_proj -1./Emax_proj)/primEnergy;
452 }
453 else Cross=0.;
454
455
456
457
458
459
460 }
461 else {
464 G4double diff1=Emin_proj-primEnergy;
465 G4double diff2=Emax_proj-primEnergy;
466 G4double t1=(1./diff1+1./Emin_proj-1./diff2-1./Emax_proj)/primEnergy;
467 //G4double t2=2.*std::log(diff2*Emin_proj/Emax_proj/diff1)/primEnergy/primEnergy;
468 G4double t2=2.*std::log(Emax_proj/Emin_proj)/primEnergy/primEnergy;
469 Cross*=(t1+t2);
470
471
472 }
473 lastCS =Cross;
474 return Cross;
475}
double G4double
Definition: G4Types.hh:83
virtual G4double GetSecondAdjEnergyMaxForScatProjToProjCase(G4double PrimAdjEnergy)
virtual G4double GetSecondAdjEnergyMinForProdToProjCase(G4double PrimAdjEnergy)
virtual G4double GetSecondAdjEnergyMinForScatProjToProjCase(G4double PrimAdjEnergy, G4double Tcut=0)
virtual G4double GetSecondAdjEnergyMaxForProdToProjCase(G4double PrimAdjEnergy)
G4double GetElectronDensity() const
Definition: G4Material.hh:215
void DefineCurrentMaterial(const G4MaterialCutsCouple *couple)
G4Material * currentMaterial
G4double currentTcutForDirectSecond
virtual G4double AdjointCrossSection(const G4MaterialCutsCouple *aCouple, G4double primEnergy, G4bool IsScatProjToProjCase)

◆ DiffCrossSectionPerAtomPrimToSecond()

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

Reimplemented from G4VEmAdjointModel.

Definition at line 305 of file G4AdjointhIonisationModel.cc.

310{//Probably that here the Bragg Model should be also used for kinEnergyProj/nuc<2MeV
311
312
313
314 G4double dSigmadEprod=0;
315 G4double Emax_proj = GetSecondAdjEnergyMaxForProdToProjCase(kinEnergyProd);
316 G4double Emin_proj = GetSecondAdjEnergyMinForProdToProjCase(kinEnergyProd);
317
318
319 if (kinEnergyProj>Emin_proj && kinEnergyProj<=Emax_proj){ //the produced particle should have a kinetic energy smaller than the projectile
320 G4double Tmax=kinEnergyProj;
321
322 G4double E1=kinEnergyProd;
323 G4double E2=kinEnergyProd*1.000001;
324 G4double dE=(E2-E1);
325 G4double sigma1,sigma2;
326 if (kinEnergyProj >2.*MeV){
327 sigma1=theDirectEMModel->ComputeCrossSectionPerAtom(theDirectPrimaryPartDef,kinEnergyProj,Z,A ,E1,1.e20);
328 sigma2=theDirectEMModel->ComputeCrossSectionPerAtom(theDirectPrimaryPartDef,kinEnergyProj,Z,A ,E2,1.e20);
329 }
330 else {
331 sigma1=theBraggDirectEMModel->ComputeCrossSectionPerAtom(theDirectPrimaryPartDef,kinEnergyProj,Z,A ,E1,1.e20);
332 sigma2=theBraggDirectEMModel->ComputeCrossSectionPerAtom(theDirectPrimaryPartDef,kinEnergyProj,Z,A ,E2,1.e20);
333 }
334
335
336 dSigmadEprod=(sigma1-sigma2)/dE;
337 if (dSigmadEprod>1.) {
338 G4cout<<"sigma1 "<<kinEnergyProj/MeV<<'\t'<<kinEnergyProd/MeV<<'\t'<<sigma1<<G4endl;
339 G4cout<<"sigma2 "<<kinEnergyProj/MeV<<'\t'<<kinEnergyProd/MeV<<'\t'<<sigma2<<G4endl;
340 G4cout<<"dsigma "<<kinEnergyProj/MeV<<'\t'<<kinEnergyProd/MeV<<'\t'<<dSigmadEprod<<G4endl;
341
342 }
343
344
345
346 //correction of differential cross section at high energy to correct for the suppression of particle at secondary at high
347 //energy used in the Bethe Bloch Model. This correction consist to multiply by g the probability function used
348 //to test the rejection of a secondary
349 //-------------------------
350
351 //Source code taken from G4BetheBlochModel::SampleSecondaries
352
353 G4double deltaKinEnergy = kinEnergyProd;
354
355 //Part of the taken code
356 //----------------------
357
358
359
360 // projectile formfactor - suppresion of high energy
361 // delta-electron production at high energy
362 G4double x = formfact*deltaKinEnergy;
363 if(x > 1.e-6) {
364
365
366 G4double totEnergy = kinEnergyProj + mass;
367 G4double etot2 = totEnergy*totEnergy;
368 G4double beta2 = kinEnergyProj*(kinEnergyProj + 2.0*mass)/etot2;
369 G4double f;
370 G4double f1 = 0.0;
371 f = 1.0 - beta2*deltaKinEnergy/Tmax;
372 if( 0.5 == spin ) {
373 f1 = 0.5*deltaKinEnergy*deltaKinEnergy/etot2;
374 f += f1;
375 }
376 G4double x1 = 1.0 + x;
377 G4double gg = 1.0/(x1*x1);
378 if( 0.5 == spin ) {
379 G4double x2 = 0.5*electron_mass_c2*deltaKinEnergy/(mass*mass);
380 gg *= (1.0 + magMoment2*(x2 - f1/f)/(1.0 + x2));
381 }
382 if(gg > 1.0) {
383 G4cout << "### G4BetheBlochModel in Adjoint Sim WARNING: g= " << g
384 << G4endl;
385 gg=1.;
386 }
387 //G4cout<<"gg"<<gg<<G4endl;
388 dSigmadEprod*=gg;
389 }
390
391 }
392
393 return dSigmadEprod;
394}
#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:359

Referenced by RapidSampleSecondaries().

◆ GetSecondAdjEnergyMaxForProdToProjCase()

G4double G4AdjointhIonisationModel::GetSecondAdjEnergyMaxForProdToProjCase ( G4double  PrimAdjEnergy)
virtual

◆ GetSecondAdjEnergyMaxForScatProjToProjCase()

G4double G4AdjointhIonisationModel::GetSecondAdjEnergyMaxForScatProjToProjCase ( G4double  PrimAdjEnergy)
virtual

Reimplemented from G4VEmAdjointModel.

Definition at line 478 of file G4AdjointhIonisationModel.cc.

479{
480 G4double Tmax=PrimAdjEnergy*one_plus_ratio_2/(one_minus_ratio_2-2.*ratio*PrimAdjEnergy/mass);
481 return Tmax;
482}

Referenced by AdjointCrossSection(), and RapidSampleSecondaries().

◆ GetSecondAdjEnergyMinForProdToProjCase()

G4double G4AdjointhIonisationModel::GetSecondAdjEnergyMinForProdToProjCase ( G4double  PrimAdjEnergy)
virtual

Reimplemented from G4VEmAdjointModel.

Definition at line 495 of file G4AdjointhIonisationModel.cc.

496{ G4double Tmin= (2*PrimAdjEnergy-4*mass + std::sqrt(4.*PrimAdjEnergy*PrimAdjEnergy +16.*mass*mass + 8.*PrimAdjEnergy*mass*(1/ratio +ratio)))/4.;
497 return Tmin;
498}

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

◆ GetSecondAdjEnergyMinForScatProjToProjCase()

G4double G4AdjointhIonisationModel::GetSecondAdjEnergyMinForScatProjToProjCase ( G4double  PrimAdjEnergy,
G4double  Tcut = 0 
)
virtual

Reimplemented from G4VEmAdjointModel.

Definition at line 485 of file G4AdjointhIonisationModel.cc.

486{ return PrimAdjEnergy+Tcut;
487}

Referenced by AdjointCrossSection(), and RapidSampleSecondaries().

◆ RapidSampleSecondaries()

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

Definition at line 159 of file G4AdjointhIonisationModel.cc.

162{
163
164 const G4DynamicParticle* theAdjointPrimary =aTrack.GetDynamicParticle();
166
167
168 G4double adjointPrimKinEnergy = theAdjointPrimary->GetKineticEnergy();
169 G4double adjointPrimP =theAdjointPrimary->GetTotalMomentum();
170
171 if (adjointPrimKinEnergy>HighEnergyLimit*0.999){
172 return;
173 }
174
175 G4double projectileKinEnergy =0.;
176 G4double eEnergy=0.;
177 G4double newCS=currentMaterial->GetElectronDensity()*twopi_mc2_rcl2*mass;
178 if (!IsScatProjToProjCase){//1/E^2 distribution
179
180 eEnergy=adjointPrimKinEnergy;
181 G4double Emax = GetSecondAdjEnergyMaxForProdToProjCase(adjointPrimKinEnergy);
182 G4double Emin= GetSecondAdjEnergyMinForProdToProjCase(adjointPrimKinEnergy);
183 if (Emin>=Emax) return;
184 G4double a=1./Emax;
185 G4double b=1./Emin;
186 newCS=newCS*(b-a)/eEnergy;
187
188 projectileKinEnergy =1./(b- (b-a)*G4UniformRand());
189
190
191 }
192 else { G4double Emax = GetSecondAdjEnergyMaxForScatProjToProjCase(adjointPrimKinEnergy);
194 if (Emin>=Emax) return;
195 G4double diff1=Emin-adjointPrimKinEnergy;
196 G4double diff2=Emax-adjointPrimKinEnergy;
197
198 G4double t1=adjointPrimKinEnergy*(1./diff1-1./diff2);
199 G4double t2=adjointPrimKinEnergy*(1./Emin-1./Emax);
200 /*G4double f31=diff1/Emin;
201 G4double f32=diff2/Emax/f31;*/
202 G4double t3=2.*std::log(Emax/Emin);
203 G4double sum_t=t1+t2+t3;
204 newCS=newCS*sum_t/adjointPrimKinEnergy/adjointPrimKinEnergy;
205 G4double t=G4UniformRand()*sum_t;
206 if (t <=t1 ){
207 G4double q= G4UniformRand()*t1/adjointPrimKinEnergy ;
208 projectileKinEnergy =adjointPrimKinEnergy +1./(1./diff1-q);
209
210 }
211 else if (t <=t2 ) {
212 G4double q= G4UniformRand()*t2/adjointPrimKinEnergy;
213 projectileKinEnergy =1./(1./Emin-q);
214 }
215 else {
216 projectileKinEnergy=Emin*std::pow(Emax/Emin,G4UniformRand());
217
218 }
219 eEnergy=projectileKinEnergy-adjointPrimKinEnergy;
220
221
222 }
223
224
225
226 G4double diffCS_perAtom_Used=twopi_mc2_rcl2*mass*adjointPrimKinEnergy/projectileKinEnergy/projectileKinEnergy/eEnergy/eEnergy;
227
228
229
230 //Weight correction
231 //-----------------------
232 //First w_corr is set to the ratio between adjoint total CS and fwd total CS
234
235 //G4cout<<w_corr<<G4endl;
236 w_corr*=newCS/lastCS;
237 //G4cout<<w_corr<<G4endl;
238 //Then another correction is needed due to the fact that a biaised differential CS has been used rather than the one consistent with the direct model
239 //Here we consider the true diffCS as the one obtained by the numerical differentiation over Tcut of the direct CS
240
241 G4double diffCS = DiffCrossSectionPerAtomPrimToSecond(projectileKinEnergy, eEnergy,1,1);
242 w_corr*=diffCS/diffCS_perAtom_Used;
243 //G4cout<<w_corr<<G4endl;
244
245 G4double new_weight = aTrack.GetWeight()*w_corr;
246 fParticleChange->SetParentWeightByProcess(false);
247 fParticleChange->SetSecondaryWeightByProcess(false);
248 fParticleChange->ProposeParentWeight(new_weight);
249
250
251
252
253 //Kinematic:
254 //we consider a two body elastic scattering for the forward processes where the projectile knock on an e- at rest and gives
255 // him part of its energy
256 //----------------------------------------------------------------------------------------
257
259 G4double projectileTotalEnergy = projectileM0+projectileKinEnergy;
260 G4double projectileP2 = projectileTotalEnergy*projectileTotalEnergy - projectileM0*projectileM0;
261
262
263
264 //Companion
265 //-----------
267 if (IsScatProjToProjCase) {
269 }
270 G4double companionTotalEnergy =companionM0+ projectileKinEnergy-adjointPrimKinEnergy;
271 G4double companionP2 = companionTotalEnergy*companionTotalEnergy - companionM0*companionM0;
272
273
274 //Projectile momentum
275 //--------------------
276 G4double P_parallel = (adjointPrimP*adjointPrimP + projectileP2 - companionP2)/(2.*adjointPrimP);
277 G4double P_perp = std::sqrt( projectileP2 - P_parallel*P_parallel);
278 G4ThreeVector dir_parallel=theAdjointPrimary->GetMomentumDirection();
279 G4double phi =G4UniformRand()*2.*3.1415926;
280 G4ThreeVector projectileMomentum = G4ThreeVector(P_perp*std::cos(phi),P_perp*std::sin(phi),P_parallel);
281 projectileMomentum.rotateUz(dir_parallel);
282
283
284
285 if (!IsScatProjToProjCase ){ //kill the primary and add a secondary
286 fParticleChange->ProposeTrackStatus(fStopAndKill);
287 fParticleChange->AddSecondary(new G4DynamicParticle(theAdjEquivOfDirectPrimPartDef,projectileMomentum));
288 //G4cout<<"projectileMomentum "<<projectileMomentum<<G4endl;
289 }
290 else {
291 fParticleChange->ProposeEnergy(projectileKinEnergy);
292 fParticleChange->ProposeMomentumDirection(projectileMomentum.unit());
293 }
294
295
296
297
298
299
300
301}
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
G4double GetPostStepWeightCorrection()
static G4AdjointCSManager * GetAdjointCSManager()
virtual G4double DiffCrossSectionPerAtomPrimToSecond(G4double kinEnergyProj, G4double kinEnergyProd, G4double Z, G4double A=0.)
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)
void SetSecondaryWeightByProcess(G4bool)
void SetParentWeightByProcess(G4bool)
void ProposeParentWeight(G4double finalWeight)

Referenced by SampleSecondaries().

◆ SampleSecondaries()

void G4AdjointhIonisationModel::SampleSecondaries ( const G4Track aTrack,
G4bool  IsScatProjToProjCase,
G4ParticleChange fParticleChange 
)
virtual

Implements G4VEmAdjointModel.

Definition at line 83 of file G4AdjointhIonisationModel.cc.

86{
87 if (!UseMatrix) return RapidSampleSecondaries(aTrack,IsScatProjToProjCase,fParticleChange);
88
89 const G4DynamicParticle* theAdjointPrimary =aTrack.GetDynamicParticle();
90
91 //Elastic inverse scattering
92 //---------------------------------------------------------
93 G4double adjointPrimKinEnergy = theAdjointPrimary->GetKineticEnergy();
94 G4double adjointPrimP =theAdjointPrimary->GetTotalMomentum();
95
96 if (adjointPrimKinEnergy>HighEnergyLimit*0.999){
97 return;
98 }
99
100 //Sample secondary energy
101 //-----------------------
102 G4double projectileKinEnergy = SampleAdjSecEnergyFromCSMatrix(adjointPrimKinEnergy, IsScatProjToProjCase);
103 CorrectPostStepWeight(fParticleChange,
104 aTrack.GetWeight(),
105 adjointPrimKinEnergy,
106 projectileKinEnergy,
107 IsScatProjToProjCase); //Caution !!!this weight correction should be always applied
108
109
110 //Kinematic:
111 //we consider a two body elastic scattering for the forward processes where the projectile knock on an e- at rest and gives
112 // him part of its energy
113 //----------------------------------------------------------------------------------------
114
116 G4double projectileTotalEnergy = projectileM0+projectileKinEnergy;
117 G4double projectileP2 = projectileTotalEnergy*projectileTotalEnergy - projectileM0*projectileM0;
118
119
120
121 //Companion
122 //-----------
124 if (IsScatProjToProjCase) {
126 }
127 G4double companionTotalEnergy =companionM0+ projectileKinEnergy-adjointPrimKinEnergy;
128 G4double companionP2 = companionTotalEnergy*companionTotalEnergy - companionM0*companionM0;
129
130
131 //Projectile momentum
132 //--------------------
133 G4double P_parallel = (adjointPrimP*adjointPrimP + projectileP2 - companionP2)/(2.*adjointPrimP);
134 G4double P_perp = std::sqrt( projectileP2 - P_parallel*P_parallel);
135 G4ThreeVector dir_parallel=theAdjointPrimary->GetMomentumDirection();
136 G4double phi =G4UniformRand()*2.*3.1415926;
137 G4ThreeVector projectileMomentum = G4ThreeVector(P_perp*std::cos(phi),P_perp*std::sin(phi),P_parallel);
138 projectileMomentum.rotateUz(dir_parallel);
139
140
141
142 if (!IsScatProjToProjCase ){ //kill the primary and add a secondary
143 fParticleChange->ProposeTrackStatus(fStopAndKill);
144 fParticleChange->AddSecondary(new G4DynamicParticle(theAdjEquivOfDirectPrimPartDef,projectileMomentum));
145 //G4cout<<"projectileMomentum "<<projectileMomentum<<G4endl;
146 }
147 else {
148 fParticleChange->ProposeEnergy(projectileKinEnergy);
149 fParticleChange->ProposeMomentumDirection(projectileMomentum.unit());
150 }
151
152
153
154
155}
void RapidSampleSecondaries(const G4Track &aTrack, G4bool IsScatProjToProjCase, G4ParticleChange *fParticleChange)
virtual void CorrectPostStepWeight(G4ParticleChange *fParticleChange, G4double old_weight, G4double adjointPrimKinEnergy, G4double projectileKinEnergy, G4bool IsScatProjToProjCase)
G4double SampleAdjSecEnergyFromCSMatrix(size_t MatrixIndex, G4double prim_energy, G4bool IsScatProjToProjCase)

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