Geant4-11
G4PolarizedIonisation.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// -------------------------------------------------------------------
27//
28// Geant4 Class file
29//
30// File name: G4PolarizedIonisation
31//
32// Author: A.Schaelicke on base of Vladimir Ivanchenko code
33
35
36#include "G4Electron.hh"
37#include "G4EmParameters.hh"
42#include "G4Positron.hh"
44#include "G4StokesVector.hh"
45#include "G4SystemOfUnits.hh"
46#include "G4UnitsTable.hh"
48
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52 , fAsymmetryTable(nullptr)
53 , fTransverseAsymmetryTable(nullptr)
54 , fIsElectron(true)
55 , fIsInitialised(false)
56{
57 verboseLevel = 0;
60 fFlucModel = nullptr;
61 fEmModel = nullptr;
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
68void G4PolarizedIonisation::ProcessDescription(std::ostream& out) const
69{
70 out << "Polarized version of G4eIonisation.\n";
71
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77{
79 {
81 delete fAsymmetryTable;
82 fAsymmetryTable = nullptr;
83 }
85 {
89 }
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94 const G4Material*,
95 G4double cut)
96{
97 G4double x = cut;
98 if(fIsElectron)
99 {
100 x += cut;
101 }
102 return x;
103}
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107{
108 return (&p == G4Electron::Electron() || &p == G4Positron::Positron());
109}
110//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
112 const G4ParticleDefinition* part, const G4ParticleDefinition*)
113{
114 if(!fIsInitialised)
115 {
116 if(part == G4Positron::Positron())
117 {
118 fIsElectron = false;
119 }
120
121 if(!FluctModel())
122 {
124 }
126
133
134 fIsInitialised = true;
135 }
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
140 G4double step,
141 G4ForceCondition* cond)
142{
143 // *** get unploarised mean free path from lambda table ***
144 G4double mfp = G4VEnergyLossProcess::GetMeanFreePath(track, step, cond);
146 {
147 mfp *= ComputeSaturationFactor(track);
148 }
149 if(verboseLevel >= 2)
150 {
151 G4cout << "G4PolarizedIonisation::MeanFreePath: " << mfp / mm << " mm "
152 << G4endl;
153 }
154 return mfp;
155}
156
157//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
159 const G4Track& track, G4double step, G4ForceCondition* cond)
160{
161 // save previous values
164
165 // *** get unpolarised mean free path from lambda table ***
166 // this changes theNumberOfInteractionLengthLeft and currentInteractionLength
168 track, step, cond);
169 G4double x0 = x;
170 G4double satFact = 1.;
171
172 // *** add corrections on polarisation ***
174 {
175 satFact = ComputeSaturationFactor(track);
176 G4double curLength = currentInteractionLength * satFact;
177 G4double prvLength = iLength * satFact;
178 if(nLength > 0.0)
179 {
181 std::max(nLength - step / prvLength, 0.0);
182 }
183 x = theNumberOfInteractionLengthLeft * curLength;
184 }
185 if(verboseLevel >= 2)
186 {
187 G4cout << "G4PolarizedIonisation::PostStepGPIL: " << std::setprecision(8)
188 << x / mm << " mm;" << G4endl
189 << " unpolarized value: " << std::setprecision(8)
190 << x0 / mm << " mm." << G4endl;
191 }
192 return x;
193}
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
197{
198 G4Material* aMaterial = track.GetMaterial();
199 G4VPhysicalVolume* aPVolume = track.GetVolume();
200 G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
201
202 G4PolarizationManager* polarizationManager =
204
205 const G4bool volumeIsPolarized = polarizationManager->IsPolarized(aLVolume);
206 G4StokesVector volPolarization =
207 polarizationManager->GetVolumePolarization(aLVolume);
208
209 G4double factor = 1.0;
210
211 if(volumeIsPolarized && !volPolarization.IsZero())
212 {
213 // *** get asymmetry, if target is polarized ***
214 const G4DynamicParticle* aDynamicPart = track.GetDynamicParticle();
215 const G4double energy = aDynamicPart->GetKineticEnergy();
216 const G4StokesVector polarization = G4StokesVector(track.GetPolarization());
217 const G4ParticleMomentum direction0 = aDynamicPart->GetMomentumDirection();
218
219 if(verboseLevel >= 2)
220 {
221 G4cout << "G4PolarizedIonisation::ComputeSaturationFactor: " << G4endl;
222 G4cout << " Energy(MeV) " << energy / MeV << G4endl;
223 G4cout << " Direction " << direction0 << G4endl;
224 G4cout << " Polarization " << polarization << G4endl;
225 G4cout << " MaterialPol. " << volPolarization << G4endl;
226 G4cout << " Phys. Volume " << aPVolume->GetName() << G4endl;
227 G4cout << " Log. Volume " << aLVolume->GetName() << G4endl;
228 G4cout << " Material " << aMaterial << G4endl;
229 }
230
231 size_t midx = CurrentMaterialCutsCoupleIndex();
232 const G4PhysicsVector* aVector = nullptr;
233 const G4PhysicsVector* bVector = nullptr;
234 if(midx < fAsymmetryTable->size())
235 {
236 aVector = (*fAsymmetryTable)(midx);
237 }
238 if(midx < fTransverseAsymmetryTable->size())
239 {
240 bVector = (*fTransverseAsymmetryTable)(midx);
241 }
242 if(aVector && bVector)
243 {
244 G4double lAsymmetry = aVector->Value(energy);
245 G4double tAsymmetry = bVector->Value(energy);
246 G4double polZZ = polarization.z() * (volPolarization * direction0);
247 G4double polXX =
248 polarization.x() *
249 (volPolarization * G4PolarizationHelper::GetParticleFrameX(direction0));
250 G4double polYY =
251 polarization.y() *
252 (volPolarization * G4PolarizationHelper::GetParticleFrameY(direction0));
253
254 factor /= (1. + polZZ * lAsymmetry + (polXX + polYY) * tAsymmetry);
255
256 if(verboseLevel >= 2)
257 {
258 G4cout << " Asymmetry: " << lAsymmetry << ", " << tAsymmetry
259 << G4endl;
260 G4cout << " PolProduct: " << polXX << ", " << polYY << ", " << polZZ
261 << G4endl;
262 G4cout << " Factor: " << factor << G4endl;
263 }
264 }
265 else
266 {
268 ed << "Problem with asymmetry tables: material index " << midx
269 << " is out of range or tables are not filled";
270 G4Exception("G4PolarizedIonisation::ComputeSaturationFactor", "em0048",
271 JustWarning, ed, "");
272 }
273 }
274 return factor;
275}
276
277//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
279{
280 // *** build DEDX and (unpolarized) cross section tables
282 G4bool master = true;
283 const G4PolarizedIonisation* masterProcess =
284 static_cast<const G4PolarizedIonisation*>(GetMasterProcess());
285 if(masterProcess && masterProcess != this)
286 {
287 master = false;
288 }
289 if(master)
290 {
292 }
293}
294
295//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
297 const G4ParticleDefinition& part)
298{
299 // cleanup old, initialise new table
300 CleanTables();
304
305 const G4ProductionCutsTable* theCoupleTable =
307 size_t numOfCouples = theCoupleTable->GetTableSize();
308
309 for(size_t j = 0; j < numOfCouples; ++j)
310 {
311 // get cut value
312 const G4MaterialCutsCouple* couple =
313 theCoupleTable->GetMaterialCutsCouple(j);
314
315 G4double cut = (*theCoupleTable->GetEnergyCutsVector(1))[j];
316
317 // create physics vectors then fill it (same parameters as lambda vector)
318 G4PhysicsVector* ptrVectorA = LambdaPhysicsVector(couple, cut);
319 G4PhysicsVector* ptrVectorB = LambdaPhysicsVector(couple, cut);
320 size_t bins = ptrVectorA->GetVectorLength();
321
322 for(size_t i = 0; i < bins; ++i)
323 {
324 G4double lowEdgeEnergy = ptrVectorA->Energy(i);
325 G4double tasm = 0.;
326 G4double asym = ComputeAsymmetry(lowEdgeEnergy, couple, part, cut, tasm);
327 ptrVectorA->PutValue(i, asym);
328 ptrVectorB->PutValue(i, tasm);
329 }
330 fAsymmetryTable->insertAt(j, ptrVectorA);
331 fTransverseAsymmetryTable->insertAt(j, ptrVectorB);
332 }
333}
334
335//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
337 G4double energy, const G4MaterialCutsCouple* couple,
338 const G4ParticleDefinition& aParticle, G4double cut, G4double& tAsymmetry)
339{
340 G4double lAsymmetry = 0.0;
341 tAsymmetry = 0.0;
342 if(fIsElectron)
343 {
344 lAsymmetry = tAsymmetry = -1.0;
345 }
346
347 // calculate polarized cross section
348 G4ThreeVector targetPolarization = G4ThreeVector(0., 0., 1.);
349 fEmModel->SetTargetPolarization(targetPolarization);
350 fEmModel->SetBeamPolarization(targetPolarization);
351 G4double sigma2 =
352 fEmModel->CrossSection(couple, &aParticle, energy, cut, energy);
353
354 // calculate transversely polarized cross section
355 targetPolarization = G4ThreeVector(1., 0., 0.);
356 fEmModel->SetTargetPolarization(targetPolarization);
357 fEmModel->SetBeamPolarization(targetPolarization);
358 G4double sigma3 =
359 fEmModel->CrossSection(couple, &aParticle, energy, cut, energy);
360
361 // calculate unpolarized cross section
362 targetPolarization = G4ThreeVector();
363 fEmModel->SetTargetPolarization(targetPolarization);
364 fEmModel->SetBeamPolarization(targetPolarization);
365 G4double sigma0 =
366 fEmModel->CrossSection(couple, &aParticle, energy, cut, energy);
367 // determine asymmetries
368 if(sigma0 > 0.)
369 {
370 lAsymmetry = sigma2 / sigma0 - 1.;
371 tAsymmetry = sigma3 / sigma0 - 1.;
372 }
373 if(std::fabs(lAsymmetry) > 1.)
374 {
376 ed << "G4PolarizedIonisation::ComputeAsymmetry : E(MeV)= " << energy
377 << " lAsymmetry= " << lAsymmetry << " (" << std::fabs(lAsymmetry) - 1.
378 << ")";
379 G4Exception("G4PolarizedIonisation::ComputeAsymmetry", "pol002",
380 JustWarning, ed);
381 }
382 if(std::fabs(tAsymmetry) > 1.)
383 {
385 ed << "G4PolarizedIonisation::ComputeAsymmetry : E(MeV)= " << energy
386 << " tAsymmetry= " << tAsymmetry << " (" << std::fabs(tAsymmetry) - 1.
387 << ")";
388 G4Exception("G4PolarizedIonisation::ComputeAsymmetry", "pol003",
389 JustWarning, ed);
390 }
391 return lAsymmetry;
392}
@ fIonisation
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
static constexpr double mm
Definition: G4SIunits.hh:95
static constexpr double MeV
Definition: G4SIunits.hh:200
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4EmParameters * Instance()
G4double MinKinEnergy() const
G4double MaxKinEnergy() const
const G4String & GetName() const
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
void clearAndDestroy()
void insertAt(std::size_t, G4PhysicsVector *)
void PutValue(const std::size_t index, const G4double value)
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
static G4ThreeVector GetParticleFrameY(const G4ThreeVector &)
static G4ThreeVector GetParticleFrameX(const G4ThreeVector &)
bool IsPolarized(G4LogicalVolume *lVol) const
const G4StokesVector GetVolumePolarization(G4LogicalVolume *lVol) const
static G4PolarizationManager * GetInstance()
void SetTargetPolarization(const G4ThreeVector &pTarget)
void SetBeamPolarization(const G4ThreeVector &pBeam)
G4double ComputeAsymmetry(G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tasm)
G4VEmFluctuationModel * fFlucModel
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4double ComputeSaturationFactor(const G4Track &aTrack)
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *, G4double cut) override
void BuildAsymmetryTables(const G4ParticleDefinition &part)
G4PolarizedIonisation(const G4String &name="pol-eIoni")
G4PhysicsTable * fTransverseAsymmetryTable
G4PhysicsTable * fAsymmetryTable
G4PolarizedIonisationModel * fEmModel
virtual void ProcessDescription(std::ostream &) const override
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
virtual ~G4PolarizedIonisation() override
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *) override
virtual G4bool IsApplicable(const G4ParticleDefinition &p) override
static G4Positron * Positron()
Definition: G4Positron.cc:93
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4bool IsZero() const
G4VPhysicalVolume * GetVolume() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
const G4ThreeVector & GetPolarization() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:774
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:539
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetFluctModel(G4VEmFluctuationModel *)
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *, G4double cut)
void ProcessDescription(std::ostream &outFile) const override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void SetEmModel(G4VEmModel *, G4int index=0)
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4VEmFluctuationModel * FluctModel() const
void BuildPhysicsTable(const G4ParticleDefinition &) override
size_t CurrentMaterialCutsCoupleIndex() const
void SetSecondaryParticle(const G4ParticleDefinition *p)
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4double currentInteractionLength
Definition: G4VProcess.hh:335
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:518
G4int verboseLevel
Definition: G4VProcess.hh:356
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4double energy(const ThreeVector &p, const G4double m)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const char * name(G4int ptype)
#define DBL_MAX
Definition: templates.hh:62