Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
G4BoldyshevTripletModel Class Reference

#include <G4BoldyshevTripletModel.hh>

Inheritance diagram for G4BoldyshevTripletModel:
G4VEmModel

Public Member Functions

 G4BoldyshevTripletModel (const G4ParticleDefinition *p=0, const G4String &nam="BoldyshevTriplet")
 
virtual ~G4BoldyshevTripletModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Member Functions

G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Detailed Description

Definition at line 42 of file G4BoldyshevTripletModel.hh.

Constructor & Destructor Documentation

G4BoldyshevTripletModel::G4BoldyshevTripletModel ( const G4ParticleDefinition p = 0,
const G4String nam = "BoldyshevTriplet" 
)

Definition at line 47 of file G4BoldyshevTripletModel.cc.

References python.hepunit::electron_mass_c2, G4cout, G4endl, python.hepunit::GeV, python.hepunit::MeV, and G4VEmModel::SetHighEnergyLimit().

49  :G4VEmModel(nam),fParticleChange(0),smallEnergy(4.*MeV),isInitialised(false),
50  crossSectionHandler(0),meanFreePathTable(0)
51 {
52  lowEnergyLimit = 4.0*electron_mass_c2;
53  highEnergyLimit = 100 * GeV;
54  SetHighEnergyLimit(highEnergyLimit);
55 
56  verboseLevel= 0;
57  // Verbosity scale:
58  // 0 = nothing
59  // 1 = warning for energy non-conservation
60  // 2 = details of energy budget
61  // 3 = calculation of cross sections, file openings, sampling of atoms
62  // 4 = entering in methods
63 
64  if(verboseLevel > 0) {
65  G4cout << "Triplet Gamma conversion is constructed " << G4endl
66  << "Energy range: "
67  << lowEnergyLimit / MeV << " MeV - "
68  << highEnergyLimit / GeV << " GeV"
69  << G4endl;
70  }
71 }
G4ParticleChangeForGamma * fParticleChange
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4GLOB_DLL std::ostream G4cout
float electron_mass_c2
Definition: hepunit.py:274
#define G4endl
Definition: G4ios.hh:61
G4BoldyshevTripletModel::~G4BoldyshevTripletModel ( )
virtual

Definition at line 75 of file G4BoldyshevTripletModel.cc.

76 {
77  if (crossSectionHandler) delete crossSectionHandler;
78 }

Member Function Documentation

G4double G4BoldyshevTripletModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 126 of file G4BoldyshevTripletModel.cc.

References G4VCrossSectionHandler::FindValue(), G4cout, and G4endl.

130 {
131  if (verboseLevel > 3) {
132  G4cout << "Calling ComputeCrossSectionPerAtom() of G4BoldyshevTripletModel"
133  << G4endl;
134  }
135  if (GammaEnergy < lowEnergyLimit || GammaEnergy > highEnergyLimit) return 0;
136 
137  G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
138  return cs;
139 }
int G4int
Definition: G4Types.hh:78
G4double FindValue(G4int Z, G4double e) const
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4BoldyshevTripletModel::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
protected
void G4BoldyshevTripletModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 83 of file G4BoldyshevTripletModel.cc.

References G4VCrossSectionHandler::Clear(), fParticleChange, G4cout, G4endl, G4VEmModel::GetParticleChangeForGamma(), python.hepunit::GeV, G4VEmModel::HighEnergyLimit(), G4VCrossSectionHandler::Initialise(), G4VCrossSectionHandler::LoadData(), G4VEmModel::LowEnergyLimit(), and python.hepunit::MeV.

85 {
86  if (verboseLevel > 3)
87  G4cout << "Calling G4BoldyshevTripletModel::Initialise()" << G4endl;
88 
89  if (crossSectionHandler)
90  {
91  crossSectionHandler->Clear();
92  delete crossSectionHandler;
93  }
94 
95  // Read data tables for all materials
96 
97  crossSectionHandler = new G4CrossSectionHandler();
98  crossSectionHandler->Initialise(0,lowEnergyLimit,100.*GeV,400);
99  G4String crossSectionFile = "tripdata/pp-trip-cs-"; // here only pair in electron field cs should be used
100  crossSectionHandler->LoadData(crossSectionFile);
101 
102  //
103 
104  if (verboseLevel > 0) {
105  G4cout << "Loaded cross section files for Livermore GammaConversion" << G4endl;
106  G4cout << "To obtain the total cross section this should be used only " << G4endl
107  << "in connection with G4NuclearGammaConversion " << G4endl;
108  }
109 
110  if (verboseLevel > 0) {
111  G4cout << "Livermore Electron Gamma Conversion model is initialized " << G4endl
112  << "Energy range: "
113  << LowEnergyLimit() / MeV << " MeV - "
114  << HighEnergyLimit() / GeV << " GeV"
115  << G4endl;
116  }
117 
118  if(isInitialised) return;
120  isInitialised = true;
121 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4ParticleChangeForGamma * fParticleChange
G4GLOB_DLL std::ostream G4cout
void Initialise(G4VDataSetAlgorithm *interpolation=0, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int numberOfBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)
void LoadData(const G4String &dataFile)
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
void G4BoldyshevTripletModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 143 of file G4BoldyshevTripletModel.cc.

References test::a, test::b, plottest35::c1, G4Electron::Electron(), python.hepunit::electron_mass_c2, fParticleChange, fStopAndKill, G4cout, G4endl, G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4INCL::Math::max(), n, python.hepunit::pi, G4Positron::Positron(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), CLHEP::Hep3Vector::rotateUz(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), and python.hepunit::twopi.

148 {
149 
150 // The energies of the secondary particles are sampled using
151 // a modified Wheeler-Lamb model (see PhysRevD 7 (1973), 26)
152 
153  if (verboseLevel > 3)
154  G4cout << "Calling SampleSecondaries() of G4BoldyshevTripletModel" << G4endl;
155 
156  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
157  G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
158 
159  G4double epsilon ;
161 
162  G4double positronTotEnergy, electronTotEnergy, thetaEle, thetaPos;
163  G4double ener_re=0., theta_re, phi_re, phi;
164 
165  // Calculo de theta - elecron de recoil
166 
167  G4double energyThreshold = sqrt(2.)*electron_mass_c2; // -> momentumThreshold_N = 1
168  energyThreshold = 1.1*electron_mass_c2;
169  // G4cout << energyThreshold << G4endl;
170 
171  G4double momentumThreshold_c = sqrt(energyThreshold * energyThreshold - electron_mass_c2*electron_mass_c2); // momentun in MeV/c unit
172  G4double momentumThreshold_N = momentumThreshold_c/electron_mass_c2; // momentun in mc unit
173 
174  // Calculation of recoil electron production
175 
176  G4double SigmaTot = (28./9.) * std::log ( 2.* photonEnergy / electron_mass_c2 ) - 218. / 27. ;
177  G4double X_0 = 2. * ( sqrt(momentumThreshold_N*momentumThreshold_N + 1) -1 );
178  G4double SigmaQ = (82./27. - (14./9.) * log (X_0) + 4./15.*X_0 - 0.0348 * X_0 * X_0);
179  G4double recoilProb = G4UniformRand();
180  //G4cout << "SIGMA TOT " << SigmaTot << " " << "SigmaQ " << SigmaQ << " " << SigmaQ/SigmaTot << " " << recoilProb << G4endl;
181 
182  if (recoilProb >= SigmaQ/SigmaTot) // create electron recoil
183  {
184 
185  G4double cosThetaMax = ( ( energyThreshold - electron_mass_c2 ) / (momentumThreshold_c) + electron_mass_c2*
186  ( energyThreshold + electron_mass_c2 ) / (photonEnergy*momentumThreshold_c) );
187 
188  if (cosThetaMax > 1) G4cout << "ERRORE " << G4endl;
189 
190  G4double r1;
191  G4double r2;
192  G4double are, bre, loga, f1_re, greject, cost;
193 
194  do {
195  r1 = G4UniformRand();
196  r2 = G4UniformRand();
197  // cost = (pow(4./enern,0.5*r1)) ;
198  cost = pow(cosThetaMax,r1);
199  theta_re = acos(cost);
200  are = 1./(14.*cost*cost);
201  bre = (1.-5.*cost*cost)/(2.*cost);
202  loga = log((1.+ cost)/(1.- cost));
203  f1_re = 1. - bre*loga;
204 
205  if ( theta_re >= 4.47*CLHEP::pi/180.)
206  {
207  greject = are*f1_re;
208  } else {
209  greject = 1. ;
210  }
211  } while(greject < r2);
212 
213  // Calculo de phi - elecron de recoil
214 
215  G4double r3, r4, rt;
216 
217  do {
218 
219  r3 = G4UniformRand();
220  r4 = G4UniformRand();
221  phi_re = twopi*r3 ;
222  G4double sint2 = 1. - cost*cost ;
223  G4double fp = 1. - sint2*loga/(2.*cost) ;
224  rt = (1.-cos(2.*phi_re)*fp/f1_re)/(2.*pi) ;
225 
226  } while(rt < r4);
227 
228  // Calculo de la energia - elecron de recoil - relacion momento maximo <-> angulo
229 
230  G4double S = electron_mass_c2*(2.* photonEnergy + electron_mass_c2);
233  *(S - electron_mass_c2*electron_mass_c2)*sin(theta_re)*sin(theta_re);
234  ener_re = electron_mass_c2 * (S + electron_mass_c2*electron_mass_c2)/sqrt(D2);
235 
236  // New Recoil energy calculation
237 
238  G4double momentum_recoil = 2* (electron_mass_c2) * (std::cos(theta_re)/(std::sin(phi_re)*std::sin(phi_re)));
239  G4double ener_recoil = sqrt( momentum_recoil*momentum_recoil + electron_mass_c2*electron_mass_c2);
240  ener_re = ener_recoil;
241 
242  // G4cout << "electron de retroceso " << ener_re << " " << theta_re << " " << phi_re << G4endl;
243 
244  // Recoil electron creation
245  G4double dxEle_re=sin(theta_re)*std::cos(phi_re),dyEle_re=sin(theta_re)*std::sin(phi_re), dzEle_re=cos(theta_re);
246 
247  G4double electronRKineEnergy = std::max(0.,ener_re - electron_mass_c2) ;
248 
249  G4ThreeVector electronRDirection (dxEle_re, dyEle_re, dzEle_re);
250  electronRDirection.rotateUz(photonDirection);
251 
253  electronRDirection,
254  electronRKineEnergy);
255  fvect->push_back(particle3);
256 
257  }
258  else
259  {
260  // deposito la energia ener_re - electron_mass_c2
261  // G4cout << "electron de retroceso " << ener_re << G4endl;
263  }
264 
265  // Depaola (2004) suggested distribution for e+e- energy
266 
267  // G4double t = 0.5*asinh(momentumThreshold_N);
268  G4double t = 0.5*log(momentumThreshold_N + sqrt(momentumThreshold_N*momentumThreshold_N+1));
269 
270  G4cout << 0.5*asinh(momentumThreshold_N) << " " << t << G4endl;
271 
272  G4double J1 = 0.5*(t*cosh(t)/sinh(t) - log(2.*sinh(t)));
273  G4double J2 = (-2./3.)*log(2.*sinh(t)) + t*cosh(t)/sinh(t) + (sinh(t)-t*pow(cosh(t),3))/(3.*pow(sinh(t),3));
274  G4double b = 2.*(J1-J2)/J1;
275 
276  G4double n = 1 - b/6.;
277  G4double re=0.;
278  re = G4UniformRand();
279  G4double a = 0.;
280  G4double b1 = 16. - 3.*b - 36.*b*re*n + 36.*b*pow(re,2.)*pow(n,2.) +
281  6.*pow(b,2.)*re*n;
282  a = pow((b1/b),0.5);
283  G4double c1 = (-6. + 12.*re*n + b + 2*a)*pow(b,2.);
284  epsilon = (pow(c1,1./3.))/(2.*b) + (b-4.)/(2.*pow(c1,1./3.))+0.5;
285 
286  G4double photonEnergy1 = photonEnergy - ener_re ; // resto al foton la energia del electron de retro.
287  positronTotEnergy = epsilon*photonEnergy1;
288  electronTotEnergy = photonEnergy1 - positronTotEnergy; // temporarly
289 
290  G4double momento_e = sqrt(electronTotEnergy*electronTotEnergy -
292  G4double momento_p = sqrt(positronTotEnergy*positronTotEnergy -
294 
295  thetaEle = acos((sqrt(p0*p0/(momento_e*momento_e) +1.)- p0/momento_e)) ;
296  thetaPos = acos((sqrt(p0*p0/(momento_p*momento_p) +1.)- p0/momento_p)) ;
297  phi = twopi * G4UniformRand();
298 
299  G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
300  G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
301 
302 
303  // Kinematics of the created pair:
304  // the electron and positron are assumed to have a symetric angular
305  // distribution with respect to the Z axis along the parent photon
306 
307  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
308 
309  // SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
310 
311  G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
312  electronDirection.rotateUz(photonDirection);
313 
315  electronDirection,
316  electronKineEnergy);
317 
318  // The e+ is always created (even with kinetic energy = 0) for further annihilation
319  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
320 
321  // SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
322 
323  G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
324  positronDirection.rotateUz(photonDirection);
325 
326  // Create G4DynamicParticle object for the particle2
328  positronDirection, positronKineEnergy);
329  // Fill output vector
330 
331 
332  fvect->push_back(particle1);
333  fvect->push_back(particle2);
334 
335 
336 
337 
338  // kill incident photon
341 
342 }
G4double GetKineticEnergy() const
G4ParticleChangeForGamma * fParticleChange
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
float electron_mass_c2
Definition: hepunit.py:274
const G4int n
static G4Positron * Positron()
Definition: G4Positron.cc:94
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
tuple c1
Definition: plottest35.py:14

Field Documentation

G4ParticleChangeForGamma* G4BoldyshevTripletModel::fParticleChange
protected

Definition at line 70 of file G4BoldyshevTripletModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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