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

#include <G4PolarizedCompton.hh>

+ Inheritance diagram for G4PolarizedCompton:

Public Member Functions

 G4PolarizedCompton (const G4String &processName="pol-compt", G4ProcessType type=fElectromagnetic)
 
virtual ~G4PolarizedCompton ()
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void PrintInfo ()
 
void SetModel (const G4String &name)
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p)=0
 
virtual void PrintInfo ()=0
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void PrintInfoDefinition ()
 
void StartTracking (G4Track *)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
 
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double ComputeCrossSectionPerAtom (G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
 
G4double MeanFreePath (const G4Track &track)
 
G4double GetLambda (G4double &kinEnergy, const G4MaterialCutsCouple *couple)
 
void SetLambdaBinning (G4int nbins)
 
G4int LambdaBinning () const
 
void SetMinKinEnergy (G4double e)
 
G4double MinKinEnergy () const
 
void SetMaxKinEnergy (G4double e)
 
G4double MaxKinEnergy () const
 
void SetMinKinEnergyPrim (G4double e)
 
const G4PhysicsTableLambdaTable () const
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t &idxRegion) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=0)
 
void SetModel (G4VEmModel *, G4int index=1)
 
G4VEmModelModel (G4int index=1)
 
G4VEmModelEmModel (G4int index=1)
 
void SetEmModel (G4VEmModel *, G4int index=1)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false)
 
const G4ElementGetCurrentElement () const
 
void SetCrossSectionBiasingFactor (G4double f, G4bool flag=true)
 
G4double CrossSectionBiasingFactor () const
 
void ActivateForcedInteraction (G4double length=0.0, const G4String &r="", G4bool flag=true)
 
void ActivateSecondaryBiasing (const G4String &region, G4double factor, G4double energyLimit)
 
void SetPolarAngleLimit (G4double a)
 
G4double PolarAngleLimit () const
 
void SetLambdaFactor (G4double val)
 
void SetIntegral (G4bool val)
 
G4bool IsIntegral () const
 
void SetApplyCuts (G4bool val)
 
void SetBuildTableFlag (G4bool val)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Member Functions

virtual void InitialiseProcess (const G4ParticleDefinition *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildAsymmetryTable (const G4ParticleDefinition &part)
 
G4double ComputeAsymmetry (G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tAsymmetry)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
- Protected Member Functions inherited from G4VEmProcess
virtual void InitialiseProcess (const G4ParticleDefinition *)=0
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double &kinEnergy, size_t index)
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
G4double RecalculateLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple)
 
G4ParticleChangeForGammaGetParticleChange ()
 
void SetParticle (const G4ParticleDefinition *p)
 
void SetSecondaryParticle (const G4ParticleDefinition *p)
 
size_t CurrentMaterialCutsCoupleIndex () const
 
G4double GetGammaEnergyCut ()
 
G4double GetElectronEnergyCut ()
 
void SetStartFromNullFlag (G4bool val)
 
void SetSplineFlag (G4bool val)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VEmProcess
G4ParticleChangeForGamma fParticleChange
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 73 of file G4PolarizedCompton.hh.

Constructor & Destructor Documentation

◆ G4PolarizedCompton()

G4PolarizedCompton::G4PolarizedCompton ( const G4String processName = "pol-compt",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 70 of file G4PolarizedCompton.cc.

71 :
72 G4VEmProcess (processName, type),
73 buildAsymmetryTable(true),
74 useAsymmetryTable(true),
75 isInitialised(false),
76 selectedModel(0),
77 mType(10),
78 theAsymmetryTable(NULL)
79{
81 SetMinKinEnergy(0.1*keV);
82 SetMaxKinEnergy(100.0*GeV);
84 emModel = 0;
85}
@ fComptonScattering
void SetMinKinEnergy(G4double e)
void SetLambdaBinning(G4int nbins)
void SetMaxKinEnergy(G4double e)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403

◆ ~G4PolarizedCompton()

G4PolarizedCompton::~G4PolarizedCompton ( )
virtual

Definition at line 89 of file G4PolarizedCompton.cc.

90{
91 if (theAsymmetryTable) {
92 theAsymmetryTable->clearAndDestroy();
93 delete theAsymmetryTable;
94 }
95}
void clearAndDestroy()

Member Function Documentation

◆ BuildAsymmetryTable()

void G4PolarizedCompton::BuildAsymmetryTable ( const G4ParticleDefinition part)
protected

Definition at line 303 of file G4PolarizedCompton.cc.

304{
305 // Access to materials
306 const G4ProductionCutsTable* theCoupleTable=
308 size_t numOfCouples = theCoupleTable->GetTableSize();
309 for(size_t i=0; i<numOfCouples; ++i) {
310 if (!theAsymmetryTable) break;
311 if (theAsymmetryTable->GetFlag(i)) {
312
313 // create physics vector and fill it
314 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
315 // use same parameters as for lambda
316 G4PhysicsVector* aVector = LambdaPhysicsVector(couple);
317 // modelManager->FillLambdaVector(aVector, couple, startFromNull);
318
319 for (G4int j = 0 ; j < LambdaBinning() ; ++j ) {
320 G4double lowEdgeEnergy = aVector->GetLowEdgeEnergy(j);
321 G4double tasm=0.;
322 G4double asym = ComputeAsymmetry(lowEdgeEnergy, couple, part, 0., tasm);
323 aVector->PutValue(j,asym);
324 }
325
326 G4PhysicsTableHelper::SetPhysicsVector(theAsymmetryTable, i, aVector);
327 }
328 }
329
330}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
static void SetPhysicsVector(G4PhysicsTable *physTable, size_t idx, G4PhysicsVector *vec)
G4bool GetFlag(size_t i) const
virtual G4double GetLowEdgeEnergy(size_t binNumber) const
void PutValue(size_t index, G4double theValue)
G4double ComputeAsymmetry(G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tAsymmetry)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *)
G4int LambdaBinning() const

Referenced by BuildPhysicsTable().

◆ BuildPhysicsTable()

void G4PolarizedCompton::BuildPhysicsTable ( const G4ParticleDefinition part)
protectedvirtual

Reimplemented from G4VProcess.

Definition at line 292 of file G4PolarizedCompton.cc.

293{
294 // *** build (unpolarized) cross section tables (Lambda)
296 if(buildAsymmetryTable)
298}
void BuildAsymmetryTable(const G4ParticleDefinition &part)
void BuildPhysicsTable(const G4ParticleDefinition &)

◆ ComputeAsymmetry()

G4double G4PolarizedCompton::ComputeAsymmetry ( G4double  energy,
const G4MaterialCutsCouple couple,
const G4ParticleDefinition particle,
G4double  cut,
G4double tAsymmetry 
)
protected

Definition at line 335 of file G4PolarizedCompton.cc.

340{
341 G4double lAsymmetry = 0.0;
342 tAsymmetry=0;
343
344 //
345 // calculate polarized cross section
346 //
347 G4ThreeVector thePolarization=G4ThreeVector(0.,0.,1.);
348 emModel->SetTargetPolarization(thePolarization);
349 emModel->SetBeamPolarization(thePolarization);
350 G4double sigma2=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
351
352 //
353 // calculate unpolarized cross section
354 //
355 thePolarization=G4ThreeVector();
356 emModel->SetTargetPolarization(thePolarization);
357 emModel->SetBeamPolarization(thePolarization);
358 G4double sigma0=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
359
360 // determine assymmetries
361 if (sigma0>0.) {
362 lAsymmetry=sigma2/sigma0-1.;
363 }
364 return lAsymmetry;
365}
CLHEP::Hep3Vector G4ThreeVector
void SetTargetPolarization(const G4ThreeVector &pTarget)
void SetBeamPolarization(const G4ThreeVector &pBeam)
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:418

Referenced by BuildAsymmetryTable().

◆ GetMeanFreePath()

G4double G4PolarizedCompton::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
protectedvirtual

Implements G4VDiscreteProcess.

Definition at line 140 of file G4PolarizedCompton.cc.

144{
145 // *** get unploarised mean free path from lambda table ***
146 G4double mfp = G4VEmProcess::GetMeanFreePath(aTrack, previousStepSize, condition);
147
148
149 if (theAsymmetryTable && useAsymmetryTable) {
150 // *** get asymmetry, if target is polarized ***
151 const G4DynamicParticle* aDynamicGamma = aTrack.GetDynamicParticle();
152 const G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
153 const G4StokesVector GammaPolarization = aTrack.GetPolarization();
154 const G4ParticleMomentum GammaDirection0 = aDynamicGamma->GetMomentumDirection();
155
156 G4Material* aMaterial = aTrack.GetMaterial();
157 G4VPhysicalVolume* aPVolume = aTrack.GetVolume();
158 G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
159
160 // G4Material* bMaterial = aLVolume->GetMaterial();
162
163 const G4bool VolumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
164 G4StokesVector ElectronPolarization = polarizationManger->GetVolumePolarization(aLVolume);
165
166 if (!VolumeIsPolarized || mfp == DBL_MAX) return mfp;
167
168 if (verboseLevel>=2) {
169
170 G4cout << " Mom " << GammaDirection0 << G4endl;
171 G4cout << " Polarization " << GammaPolarization << G4endl;
172 G4cout << " MaterialPol. " << ElectronPolarization << G4endl;
173 G4cout << " Phys. Volume " << aPVolume->GetName() << G4endl;
174 G4cout << " Log. Volume " << aLVolume->GetName() << G4endl;
175 G4cout << " Material " << aMaterial << G4endl;
176 }
177
179 G4PhysicsVector * aVector=(*theAsymmetryTable)(midx);
180
181 G4double asymmetry=0;
182 if (aVector) {
183 G4bool isOutRange;
184 asymmetry = aVector->GetValue(GammaEnergy, isOutRange);
185 } else {
186 G4cout << " MaterialIndex " << midx << " is out of range \n";
187 asymmetry=0;
188 }
189
190 // we have to determine angle between particle motion
191 // and target polarisation here
192 // circ pol * Vec(ElectronPol)*Vec(PhotonMomentum)
193 // both vectors in global reference frame
194
195 G4double pol=ElectronPolarization*GammaDirection0;
196
197 G4double polProduct = GammaPolarization.p3() * pol;
198 mfp *= 1. / ( 1. + polProduct * asymmetry );
199
200 if (verboseLevel>=2) {
201 G4cout << " MeanFreePath: " << mfp / mm << " mm " << G4endl;
202 G4cout << " Asymmetry: " << asymmetry << G4endl;
203 G4cout << " PolProduct: " << polProduct << G4endl;
204 }
205 }
206
207 return mfp;
208}
G4double condition(const G4ErrorSymMatrix &m)
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4String GetName() const
G4double GetValue(G4double theEnergy, G4bool &isOutRange)
bool IsPolarized(G4LogicalVolume *lVol) const
static G4PolarizationManager * GetInstance()
const G4ThreeVector & GetVolumePolarization(G4LogicalVolume *lVol) const
G4double p3() const
G4VPhysicalVolume * GetVolume() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
const G4ThreeVector & GetPolarization() const
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
size_t CurrentMaterialCutsCoupleIndex() const
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4int verboseLevel
Definition: G4VProcess.hh:368
#define DBL_MAX
Definition: templates.hh:83

◆ InitialiseProcess()

void G4PolarizedCompton::InitialiseProcess ( const G4ParticleDefinition )
protectedvirtual

Implements G4VEmProcess.

Definition at line 99 of file G4PolarizedCompton.cc.

100{
101 if(!isInitialised) {
102 isInitialised = true;
103 SetBuildTableFlag(true);
105 G4double emin = MinKinEnergy();
106 G4double emax = MaxKinEnergy();
107 emModel = new G4PolarizedComptonModel();
108 if(0 == mType) selectedModel = new G4KleinNishinaCompton();
109 else if(10 == mType) selectedModel = emModel;
110 selectedModel->SetLowEnergyLimit(emin);
111 selectedModel->SetHighEnergyLimit(emax);
112 AddEmModel(1, selectedModel);
113 }
114}
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetBuildTableFlag(G4bool val)
void SetSecondaryParticle(const G4ParticleDefinition *p)
G4double MaxKinEnergy() const
G4double MinKinEnergy() const

◆ IsApplicable()

G4bool G4PolarizedCompton::IsApplicable ( const G4ParticleDefinition p)
inlinevirtual

Implements G4VEmProcess.

Definition at line 138 of file G4PolarizedCompton.hh.

139{
140 return (&p == G4Gamma::Gamma());
141}
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86

◆ PostStepGetPhysicalInteractionLength()

G4double G4PolarizedCompton::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
protectedvirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 210 of file G4PolarizedCompton.cc.

214{
215 // *** get unploarised mean free path from lambda table ***
217
218
219 if (theAsymmetryTable && useAsymmetryTable) {
220 // *** get asymmetry, if target is polarized ***
221 const G4DynamicParticle* aDynamicGamma = aTrack.GetDynamicParticle();
222 const G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
223 const G4StokesVector GammaPolarization = aTrack.GetPolarization();
224 const G4ParticleMomentum GammaDirection0 = aDynamicGamma->GetMomentumDirection();
225
226 G4Material* aMaterial = aTrack.GetMaterial();
227 G4VPhysicalVolume* aPVolume = aTrack.GetVolume();
228 G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
229
230 // G4Material* bMaterial = aLVolume->GetMaterial();
232
233 const G4bool VolumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
234 G4StokesVector ElectronPolarization = polarizationManger->GetVolumePolarization(aLVolume);
235
236 if (!VolumeIsPolarized || mfp == DBL_MAX) return mfp;
237
238 if (verboseLevel>=2) {
239
240 G4cout << " Mom " << GammaDirection0 << G4endl;
241 G4cout << " Polarization " << GammaPolarization << G4endl;
242 G4cout << " MaterialPol. " << ElectronPolarization << G4endl;
243 G4cout << " Phys. Volume " << aPVolume->GetName() << G4endl;
244 G4cout << " Log. Volume " << aLVolume->GetName() << G4endl;
245 G4cout << " Material " << aMaterial << G4endl;
246 }
247
249 G4PhysicsVector * aVector=(*theAsymmetryTable)(midx);
250
251 G4double asymmetry=0;
252 if (aVector) {
253 G4bool isOutRange;
254 asymmetry = aVector->GetValue(GammaEnergy, isOutRange);
255 } else {
256 G4cout << " MaterialIndex " << midx << " is out of range \n";
257 asymmetry=0;
258 }
259
260 // we have to determine angle between particle motion
261 // and target polarisation here
262 // circ pol * Vec(ElectronPol)*Vec(PhotonMomentum)
263 // both vectors in global reference frame
264
265 G4double pol=ElectronPolarization*GammaDirection0;
266
267 G4double polProduct = GammaPolarization.p3() * pol;
268 mfp *= 1. / ( 1. + polProduct * asymmetry );
269
270 if (verboseLevel>=2) {
271 G4cout << " MeanFreePath: " << mfp / mm << " mm " << G4endl;
272 G4cout << " Asymmetry: " << asymmetry << G4endl;
273 G4cout << " PolProduct: " << polProduct << G4endl;
274 }
275 }
276
277 return mfp;
278}
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)

◆ PreparePhysicsTable()

void G4PolarizedCompton::PreparePhysicsTable ( const G4ParticleDefinition part)
protectedvirtual

Reimplemented from G4VProcess.

Definition at line 282 of file G4PolarizedCompton.cc.

283{
285 if(buildAsymmetryTable)
286 theAsymmetryTable = G4PhysicsTableHelper::PreparePhysicsTable(theAsymmetryTable);
287}
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
void PreparePhysicsTable(const G4ParticleDefinition &)

◆ PrintInfo()

void G4PolarizedCompton::PrintInfo ( )
virtual

Implements G4VEmProcess.

Definition at line 118 of file G4PolarizedCompton.cc.

119{
120 G4cout << " Total cross sections has a good parametrisation"
121 << " from 10 KeV to (100/Z) GeV"
122 << "\n Sampling according " << selectedModel->GetName() << " model"
123 << G4endl;
124}
const G4String & GetName() const
Definition: G4VEmModel.hh:655

◆ SetModel()

void G4PolarizedCompton::SetModel ( const G4String name)

Definition at line 128 of file G4PolarizedCompton.cc.

129{
130 if(ss == "Klein-Nishina") mType = 0;
131 if(ss == "Polarized-Compton") mType = 10;
132}

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