Geant4  10.02.p03
G4PolarizedCompton Class Reference

#include <G4PolarizedCompton.hh>

Inheritance diagram for G4PolarizedCompton:
Collaboration diagram for G4PolarizedCompton:

Public Member Functions

 G4PolarizedCompton (const G4String &processName="pol-compt", G4ProcessType type=fElectromagnetic)
 
virtual ~G4PolarizedCompton ()
 
virtual 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 void ProcessDescription (std::ostream &outFile) const
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void StartTracking (G4Track *)
 
virtual G4VParticleChange * PostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
virtual 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)
 
void SetMinKinEnergy (G4double e)
 
void SetMinKinEnergyPrim (G4double e)
 
void SetMaxKinEnergy (G4double e)
 
G4PhysicsTableLambdaTable () const
 
G4PhysicsTableLambdaTablePrim () const
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t &idxRegion) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=0)
 
G4VEmModelEmModel (G4int index=1) const
 
void SetEmModel (G4VEmModel *, G4int index=1)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
const G4ElementGetCurrentElement () const
 
const G4VEmModelGetCurrentModel () 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 SetIntegral (G4bool val)
 
void SetBuildTableFlag (G4bool val)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * AlongStepDoIt (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
 
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)
 
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 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
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

virtual void InitialiseProcess (const G4ParticleDefinition *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
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 G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double &kinEnergy, size_t index)
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
G4int LambdaBinning () const
 
G4double MinKinEnergy () const
 
G4double MaxKinEnergy () const
 
G4double PolarAngleLimit () const
 
G4bool IsIntegral () const
 
G4double RecalculateLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple)
 
G4ParticleChangeForGamma * GetParticleChange ()
 
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)
 
const G4ElementGetTargetElement () const
 
const G4IsotopeGetTargetIsotope () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Private Member Functions

void CleanTable ()
 
void BuildAsymmetryTable (const G4ParticleDefinition &part)
 
G4double ComputeAsymmetry (G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tAsymmetry)
 
G4double ComputeSaturationFactor (const G4Track &aTrack)
 
G4PolarizedComptonoperator= (const G4PolarizedCompton &right)
 
 G4PolarizedCompton (const G4PolarizedCompton &)
 

Private Attributes

G4bool buildAsymmetryTable
 
G4bool useAsymmetryTable
 
G4bool isInitialised
 
G4int mType
 
G4PolarizedComptonModelemModel
 
G4ThreeVector targetPolarization
 

Static Private Attributes

static G4PhysicsTabletheAsymmetryTable = nullptr
 

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
 
G4VParticleChange * pParticleChange
 
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 72 of file G4PolarizedCompton.hh.

Constructor & Destructor Documentation

◆ G4PolarizedCompton() [1/2]

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

Definition at line 73 of file G4PolarizedCompton.cc.

74  :
75  G4VEmProcess (processName, type),
76  buildAsymmetryTable(true),
77  useAsymmetryTable(true),
78  isInitialised(false),
79  mType(10),
80  targetPolarization(0.0,0.0,0.0)
81 {
83  SetBuildTableFlag(true);
87  SetSplineFlag(true);
88  emModel = nullptr;
89 }
static const double MeV
Definition: G4SIunits.hh:211
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
Definition: G4VEmProcess.cc:91
void SetBuildTableFlag(G4bool val)
void SetSplineFlag(G4bool val)
G4PolarizedComptonModel * emModel
void SetStartFromNullFlag(G4bool val)
void SetMinKinEnergyPrim(G4double e)
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
void SetSecondaryParticle(const G4ParticleDefinition *p)
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4ThreeVector targetPolarization
Here is the call graph for this function:

◆ ~G4PolarizedCompton()

G4PolarizedCompton::~G4PolarizedCompton ( )
virtual

Definition at line 93 of file G4PolarizedCompton.cc.

94 {
95  CleanTable();
96 }
Here is the call graph for this function:

◆ G4PolarizedCompton() [2/2]

G4PolarizedCompton::G4PolarizedCompton ( const G4PolarizedCompton )
private

Member Function Documentation

◆ BuildAsymmetryTable()

void G4PolarizedCompton::BuildAsymmetryTable ( const G4ParticleDefinition part)
private

Definition at line 285 of file G4PolarizedCompton.cc.

286 {
287  // cleanup old, initialise new table
288  CleanTable();
291 
292  // Access to materials
293  const G4ProductionCutsTable* theCoupleTable=
295  size_t numOfCouples = theCoupleTable->GetTableSize();
296  if(!theAsymmetryTable) { return; }
297  G4int nbins = LambdaBinning();
298  G4double emin = MinKinEnergy();
300  G4PhysicsLogVector* aVector = 0;
301  G4PhysicsLogVector* bVector = 0;
302 
303  for(size_t i=0; i<numOfCouples; ++i) {
304  if (theAsymmetryTable->GetFlag(i)) {
305 
306  // create physics vector and fill it
307  const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
308  // use same parameters as for lambda
309  if(!aVector) {
310  aVector = new G4PhysicsLogVector(emin, emax, nbins);
311  aVector->SetSpline(true);
312  bVector = aVector;
313  } else {
314  bVector = new G4PhysicsLogVector(*aVector);
315  }
316 
317  for (G4int j = 0; j <= nbins; ++j ) {
318  G4double energy = bVector->Energy(j);
319  G4double tasm=0.;
320  G4double asym = ComputeAsymmetry(energy, couple, part, 0., tasm);
321  bVector->PutValue(j,asym);
322  }
324  }
325  }
326 }
static void SetPhysicsVector(G4PhysicsTable *physTable, size_t idx, G4PhysicsVector *vec)
G4double MinKinEnergy() const
G4double MaxKinEnergy() const
G4int LambdaBinning() const
int G4int
Definition: G4Types.hh:78
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
void SetSpline(G4bool)
double energy
Definition: plottest35.C:25
void PutValue(size_t index, G4double theValue)
static const G4double emax
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
G4double ComputeAsymmetry(G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tAsymmetry)
G4bool GetFlag(size_t i) const
G4double Energy(size_t index) const
double G4double
Definition: G4Types.hh:76
static G4PhysicsTable * theAsymmetryTable
Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildPhysicsTable()

void G4PolarizedCompton::BuildPhysicsTable ( const G4ParticleDefinition part)
protectedvirtual

Reimplemented from G4VEmProcess.

Definition at line 270 of file G4PolarizedCompton.cc.

271 {
272  // *** build (unpolarized) cross section tables (Lambda)
274  if(buildAsymmetryTable && emModel) {
275  G4bool isMaster = true;
276  const G4PolarizedCompton* masterProcess =
277  static_cast<const G4PolarizedCompton*>(GetMasterProcess());
278  if(masterProcess && masterProcess != this) { isMaster = false; }
279  if(isMaster) { BuildAsymmetryTable(part); }
280  }
281 }
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
void BuildAsymmetryTable(const G4ParticleDefinition &part)
G4PolarizedComptonModel * emModel
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
bool G4bool
Definition: G4Types.hh:79
Here is the call graph for this function:

◆ CleanTable()

void G4PolarizedCompton::CleanTable ( )
private

Definition at line 100 of file G4PolarizedCompton.cc.

101 {
102  if( theAsymmetryTable) {
104  delete theAsymmetryTable;
105  theAsymmetryTable = nullptr;
106  }
107 }
static G4PhysicsTable * theAsymmetryTable
void clearAndDestroy()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeAsymmetry()

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

Definition at line 330 of file G4PolarizedCompton.cc.

335 {
336  G4double lAsymmetry = 0.0;
337  tAsymmetry=0;
338 
339  //
340  // calculate polarized cross section
341  //
342  G4ThreeVector thePolarization=G4ThreeVector(0.,0.,1.);
343  emModel->SetTargetPolarization(thePolarization);
344  emModel->SetBeamPolarization(thePolarization);
345  G4double sigma2=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
346 
347  //
348  // calculate unpolarized cross section
349  //
350  thePolarization=G4ThreeVector();
351  emModel->SetTargetPolarization(thePolarization);
352  emModel->SetBeamPolarization(thePolarization);
353  G4double sigma0=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
354 
355  // determine assymmetries
356  if (sigma0 > 0.) {
357  lAsymmetry = sigma2/sigma0-1.;
358  }
359  return lAsymmetry;
360 }
CLHEP::Hep3Vector G4ThreeVector
G4PolarizedComptonModel * emModel
double energy
Definition: plottest35.C:25
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:501
void SetBeamPolarization(const G4ThreeVector &pBeam)
double G4double
Definition: G4Types.hh:76
void SetTargetPolarization(const G4ThreeVector &pTarget)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeSaturationFactor()

G4double G4PolarizedCompton::ComputeSaturationFactor ( const G4Track &  aTrack)
private

Definition at line 204 of file G4PolarizedCompton.cc.

205 {
206  G4double factor = 1.0;
207 
208  // *** get asymmetry, if target is polarized ***
209  const G4DynamicParticle* aDynamicGamma = aTrack.GetDynamicParticle();
210  const G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
211  const G4StokesVector GammaPolarization = aTrack.GetPolarization();
212  const G4ParticleMomentum GammaDirection0 = aDynamicGamma->GetMomentumDirection();
213 
214  G4Material* aMaterial = aTrack.GetMaterial();
215  G4VPhysicalVolume* aPVolume = aTrack.GetVolume();
216  G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
217 
218  // G4Material* bMaterial = aLVolume->GetMaterial();
220 
221  const G4bool VolumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
222  G4StokesVector ElectronPolarization = polarizationManger->GetVolumePolarization(aLVolume);
223 
224  if (VolumeIsPolarized) {
225 
226  if (verboseLevel>=2) {
227  G4cout << "G4PolarizedCompton::ComputeSaturationFactor: " << G4endl;
228  G4cout << " Mom " << GammaDirection0 << G4endl;
229  G4cout << " Polarization " << GammaPolarization << G4endl;
230  G4cout << " MaterialPol. " << ElectronPolarization << G4endl;
231  G4cout << " Phys. Volume " << aPVolume->GetName() << G4endl;
232  G4cout << " Log. Volume " << aLVolume->GetName() << G4endl;
233  G4cout << " Material " << aMaterial << G4endl;
234  }
235 
236  size_t midx = CurrentMaterialCutsCoupleIndex();
237  const G4PhysicsVector* aVector = nullptr;
238  if(midx < theAsymmetryTable->size()) {
239  aVector = (*theAsymmetryTable)(midx);
240  }
241  if (aVector) {
242  G4double asymmetry = aVector->Value(GammaEnergy);
243 
244  // we have to determine angle between particle motion
245  // and target polarisation here
246  // circ pol * Vec(ElectronPol)*Vec(PhotonMomentum)
247  // both vectors in global reference frame
248 
249  G4double pol = ElectronPolarization*GammaDirection0;
250  G4double polProduct = GammaPolarization.p3() * pol;
251  factor /= (1. + polProduct * asymmetry);
252  if (verboseLevel>=2) {
253  G4cout << " Asymmetry: " << asymmetry << G4endl;
254  G4cout << " PolProduct: " << polProduct << G4endl;
255  G4cout << " Factor: " << factor << G4endl;
256  }
257  } else {
259  ed << "Problem with asymmetry table: material index " << midx
260  << " is out of range or the table is not filled";
261  G4Exception("G4PolarizedComptonModel::ComputeSaturationFactor","em0048",
262  JustWarning, ed, "");
263  }
264  }
265  return factor;
266 }
G4int verboseLevel
Definition: G4VProcess.hh:368
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
const G4ThreeVector & GetVolumePolarization(G4LogicalVolume *lVol) const
static G4PolarizationManager * GetInstance()
size_t CurrentMaterialCutsCoupleIndex() const
G4double GetKineticEnergy() const
bool IsPolarized(G4LogicalVolume *lVol) const
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
bool G4bool
Definition: G4Types.hh:79
G4double p3() const
G4double Value(G4double theEnergy, size_t &lastidx) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4String & GetName() const
static const G4double factor
const G4ThreeVector & GetMomentumDirection() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4LogicalVolume * GetLogicalVolume() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetMeanFreePath()

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

Reimplemented from G4VEmProcess.

Definition at line 155 of file G4PolarizedCompton.cc.

158 {
159  // *** get unploarised mean free path from lambda table ***
160  G4double mfp = G4VEmProcess::GetMeanFreePath(aTrack, previousStepSize, condition);
161 
162  if (theAsymmetryTable && useAsymmetryTable && mfp < DBL_MAX) {
163  mfp *= ComputeSaturationFactor(aTrack);
164  }
165  if (verboseLevel>=2) {
166  G4cout << "G4PolarizedCompton::MeanFreePath: " << mfp / mm << " mm " << G4endl;
167  }
168  return mfp;
169 }
G4double condition(const G4ErrorSymMatrix &m)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double ComputeSaturationFactor(const G4Track &aTrack)
G4GLOB_DLL std::ostream G4cout
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static G4PhysicsTable * theAsymmetryTable
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:

◆ InitialiseProcess()

void G4PolarizedCompton::InitialiseProcess ( const G4ParticleDefinition )
protectedvirtual

Implements G4VEmProcess.

Definition at line 118 of file G4PolarizedCompton.cc.

119 {
120  if(!isInitialised) {
121  isInitialised = true;
122  if(0 == mType) {
123  if(!EmModel(1)) { SetEmModel(new G4KleinNishinaCompton(), 1); }
124  } else {
126  SetEmModel(emModel, 1);
127  }
129  EmModel(1)->SetLowEnergyLimit(param->MinKinEnergy());
131  AddEmModel(1, EmModel(1));
132  }
133 }
G4PolarizedComptonModel * emModel
G4VEmModel * EmModel(G4int index=1) const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
void SetEmModel(G4VEmModel *, G4int index=1)
G4double MinKinEnergy() const
G4double MaxKinEnergy() const
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
static G4EmParameters * Instance()
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
Here is the call graph for this function:

◆ IsApplicable()

G4bool G4PolarizedCompton::IsApplicable ( const G4ParticleDefinition p)
virtual

Implements G4VEmProcess.

Definition at line 111 of file G4PolarizedCompton.cc.

112 {
113  return (&p == G4Gamma::Gamma());
114 }
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
Here is the call graph for this function:

◆ operator=()

G4PolarizedCompton& G4PolarizedCompton::operator= ( const G4PolarizedCompton right)
private

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VEmProcess.

Definition at line 173 of file G4PolarizedCompton.cc.

177 {
178  // save previous value
180 
181  // *** compute uppolarized step limit ***
183  previousStepSize,
184  condition);
185 
186  // *** add corrections on polarisation ***
189  if(nLength > 0.0) {
191  std::max(nLength - previousStepSize/curLength, 0.0);
192  }
193  x = theNumberOfInteractionLengthLeft * curLength;
194  }
195  if (verboseLevel>=2) {
196  G4cout << "G4PolarizedCompton::PostStepGetPhysicalInteractionLength: "
197  << x/mm << " mm " << G4endl;
198  }
199  return x;
200 }
G4double condition(const G4ErrorSymMatrix &m)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4double ComputeSaturationFactor(const G4Track &aTrack)
G4GLOB_DLL std::ostream G4cout
G4double currentInteractionLength
Definition: G4VProcess.hh:297
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static G4PhysicsTable * theAsymmetryTable
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:

◆ PrintInfo()

void G4PolarizedCompton::PrintInfo ( )
virtual

Implements G4VEmProcess.

Definition at line 137 of file G4PolarizedCompton.cc.

138 {
139  G4cout << " Total cross sections has a good parametrisation"
140  << " from 10 KeV to (100/Z) GeV"
141  << "\n Sampling according " << EmModel(1)->GetName() << " model"
142  << G4endl;
143 }
const G4String & GetName() const
Definition: G4VEmModel.hh:795
G4VEmModel * EmModel(G4int index=1) const
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the call graph for this function:

◆ SetModel()

void G4PolarizedCompton::SetModel ( const G4String name)

Definition at line 147 of file G4PolarizedCompton.cc.

148 {
149  if(ss == "Klein-Nishina") { mType = 0; }
150  if(ss == "Polarized-Compton") { mType = 10; }
151 }

Member Data Documentation

◆ buildAsymmetryTable

G4bool G4PolarizedCompton::buildAsymmetryTable
private

Definition at line 123 of file G4PolarizedCompton.hh.

◆ emModel

G4PolarizedComptonModel* G4PolarizedCompton::emModel
private

Definition at line 130 of file G4PolarizedCompton.hh.

◆ isInitialised

G4bool G4PolarizedCompton::isInitialised
private

Definition at line 126 of file G4PolarizedCompton.hh.

◆ mType

G4int G4PolarizedCompton::mType
private

Definition at line 127 of file G4PolarizedCompton.hh.

◆ targetPolarization

G4ThreeVector G4PolarizedCompton::targetPolarization
private

Definition at line 132 of file G4PolarizedCompton.hh.

◆ theAsymmetryTable

G4PhysicsTable * G4PolarizedCompton::theAsymmetryTable = nullptr
staticprivate

Definition at line 131 of file G4PolarizedCompton.hh.

◆ useAsymmetryTable

G4bool G4PolarizedCompton::useAsymmetryTable
private

Definition at line 124 of file G4PolarizedCompton.hh.


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