source: trunk/source/processes/electromagnetic/utils/include/G4VEmModel.hh @ 1315

Last change on this file since 1315 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 18.7 KB
Line 
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// $Id: G4VEmModel.hh,v 1.75 2010/05/26 10:41:34 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class header file
32//
33//
34// File name:     G4VEmModel
35//
36// Author:        Vladimir Ivanchenko
37//
38// Creation date: 03.01.2002
39//
40// Modifications:
41//
42// 23-12-02 V.Ivanchenko change interface before move to cut per region
43// 24-01-03 Cut per region (V.Ivanchenko)
44// 13-02-03 Add name (V.Ivanchenko)
45// 25-02-03 Add sample theta and displacement (V.Ivanchenko)
46// 23-07-03 Replace G4Material by G4MaterialCutCouple in dE/dx and CrossSection
47//          calculation (V.Ivanchenko)
48// 01-03-04 L.Urban signature changed in SampleCosineTheta
49// 23-04-04 L.urban signature of SampleCosineTheta changed back
50// 17-11-04 Add method CrossSectionPerAtom (V.Ivanchenko)
51// 14-03-05 Reduce number of pure virtual methods and make inline part
52//          separate (V.Ivanchenko)
53// 24-03-05 Remove IsInCharge and add G4VParticleChange in the constructor (VI)
54// 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
55// 15-04-05 optimize internal interface for msc (V.Ivanchenko)
56// 08-05-05 A -> N (V.Ivanchenko)
57// 25-07-05 Move constructor and destructor to the body (V.Ivanchenko)
58// 02-02-06 ComputeCrossSectionPerAtom: default value A=0. (mma)
59// 06-02-06 add method ComputeMeanFreePath() (mma)
60// 07-03-06 Optimize msc methods (V.Ivanchenko)
61// 29-06-06 Add member currentElement and Get/Set methods (V.Ivanchenko)
62// 29-10-07 Added SampleScattering (V.Ivanchenko)
63// 15-07-08 Reorder class members and improve comments (VI)
64// 21-07-08 Added vector of G4ElementSelector and methods to use it (VI)
65// 12-09-08 Added methods GetParticleCharge, GetChargeSquareRatio,
66//          CorrectionsAlongStep, ActivateNuclearStopping (VI)
67// 16-02-09 Moved implementations of virtual methods to source (VI)
68// 07-04-09 Moved msc methods from G4VEmModel to G4VMscModel (VI)
69//
70// Class Description:
71//
72// Abstract interface to energy loss models
73
74// -------------------------------------------------------------------
75//
76
77#ifndef G4VEmModel_h
78#define G4VEmModel_h 1
79
80#include "globals.hh"
81#include "G4DynamicParticle.hh"
82#include "G4ParticleDefinition.hh"
83#include "G4MaterialCutsCouple.hh"
84#include "G4Material.hh"
85#include "G4Element.hh"
86#include "G4ElementVector.hh"
87#include "G4DataVector.hh"
88#include "G4VEmFluctuationModel.hh"
89#include "G4EmElementSelector.hh"
90#include "Randomize.hh"
91#include <vector>
92
93class G4PhysicsTable;
94class G4Region;
95class G4VParticleChange;
96class G4ParticleChangeForLoss;
97class G4ParticleChangeForGamma;
98class G4Track;
99
100class G4VEmModel
101{
102
103public:
104
105  G4VEmModel(const G4String& nam);
106
107  virtual ~G4VEmModel();
108
109  //------------------------------------------------------------------------
110  // Virtual methods to be implemented for any concrete model
111  //------------------------------------------------------------------------
112
113  virtual void Initialise(const G4ParticleDefinition*, 
114                          const G4DataVector&) = 0;
115
116  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
117                                 const G4MaterialCutsCouple*,
118                                 const G4DynamicParticle*,
119                                 G4double tmin = 0.0,
120                                 G4double tmax = DBL_MAX) = 0;
121
122  //------------------------------------------------------------------------
123  // Methods with standard implementation; may be overwritten if needed
124  //------------------------------------------------------------------------
125
126  // main method to compute dEdx
127  virtual G4double ComputeDEDXPerVolume(const G4Material*,
128                                        const G4ParticleDefinition*,
129                                        G4double kineticEnergy,
130                                        G4double cutEnergy = DBL_MAX);
131
132  // main method to compute cross section per Volume
133  virtual G4double CrossSectionPerVolume(const G4Material*,
134                                         const G4ParticleDefinition*,
135                                         G4double kineticEnergy,
136                                         G4double cutEnergy = 0.0,
137                                         G4double maxEnergy = DBL_MAX);
138
139  // main method to compute cross section per atom
140  virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
141                                              G4double kinEnergy, 
142                                              G4double Z, 
143                                              G4double A = 0., /* amu */
144                                              G4double cutEnergy = 0.0,
145                                              G4double maxEnergy = DBL_MAX);
146                                                                     
147  // min cut in kinetic energy allowed by the model
148  virtual G4double MinEnergyCut(const G4ParticleDefinition*,
149                                const G4MaterialCutsCouple*);
150
151  // Compute effective ion charge square
152  virtual G4double ChargeSquareRatio(const G4Track&);
153
154  // Compute effective ion charge square
155  virtual G4double GetChargeSquareRatio(const G4ParticleDefinition*,
156                                        const G4Material*,
157                                        G4double kineticEnergy);
158
159  // Compute ion charge
160  virtual G4double GetParticleCharge(const G4ParticleDefinition*,
161                                     const G4Material*,
162                                     G4double kineticEnergy);
163
164  // add correction to energy loss and compute non-ionizing energy loss
165  virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
166                                    const G4DynamicParticle*,
167                                    G4double& eloss,
168                                    G4double& niel,
169                                    G4double length);
170
171  // sample PIXE deexcitation
172  virtual void SampleDeexcitationAlongStep(const G4Material*,
173                                           const G4Track&,
174                                           G4double& eloss);
175
176  // initilisation at run time for a given material
177  virtual void SetupForMaterial(const G4ParticleDefinition*,
178                                const G4Material*,
179                                G4double kineticEnergy);
180
181  // add a region for the model
182  virtual void DefineForRegion(const G4Region*);
183
184protected:
185
186  // initialisation of the ParticleChange for the model
187  G4ParticleChangeForLoss* GetParticleChangeForLoss();
188
189  // initialisation of the ParticleChange for the model
190  G4ParticleChangeForGamma* GetParticleChangeForGamma();
191
192  // kinematically allowed max kinetic energy of a secondary
193  virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
194                                      G4double kineticEnergy); 
195
196public:
197
198  //------------------------------------------------------------------------
199  // Generic methods common to all models
200  //------------------------------------------------------------------------
201
202  // should be called at initialisation to build element selectors
203  void InitialiseElementSelectors(const G4ParticleDefinition*, 
204                                  const G4DataVector&);
205
206  // dEdx per unit length
207  inline G4double ComputeDEDX(const G4MaterialCutsCouple*,
208                              const G4ParticleDefinition*,
209                              G4double kineticEnergy,
210                              G4double cutEnergy = DBL_MAX);
211
212  // cross section per volume
213  inline G4double CrossSection(const G4MaterialCutsCouple*,
214                               const G4ParticleDefinition*,
215                               G4double kineticEnergy,
216                               G4double cutEnergy = 0.0,
217                               G4double maxEnergy = DBL_MAX);
218
219  // compute mean free path via cross section per volume
220  inline G4double ComputeMeanFreePath(const G4ParticleDefinition*,
221                                      G4double kineticEnergy,
222                                      const G4Material*,   
223                                      G4double cutEnergy = 0.0,
224                                      G4double maxEnergy = DBL_MAX);
225
226  // generic cross section per element
227  inline G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*,
228                                             const G4Element*,
229                                             G4double kinEnergy, 
230                                             G4double cutEnergy = 0.0,
231                                             G4double maxEnergy = DBL_MAX);
232
233  // select isotope in order to have precise mass of the nucleus
234  inline G4int SelectIsotopeNumber(const G4Element*);
235
236  // atom can be selected effitiantly if element selectors are initialised
237  inline const G4Element* SelectRandomAtom(const G4MaterialCutsCouple*,
238                                           const G4ParticleDefinition*,
239                                           G4double kineticEnergy,
240                                           G4double cutEnergy = 0.0,
241                                           G4double maxEnergy = DBL_MAX);
242
243  // to select atom cross section per volume is recomputed for each element
244  const G4Element* SelectRandomAtom(const G4Material*,
245                                    const G4ParticleDefinition*,
246                                    G4double kineticEnergy,
247                                    G4double cutEnergy = 0.0,
248                                    G4double maxEnergy = DBL_MAX);
249
250  //------------------------------------------------------------------------
251  // Get/Set methods
252  //------------------------------------------------------------------------
253
254  inline G4VEmFluctuationModel* GetModelOfFluctuations();
255
256  inline G4double HighEnergyLimit() const;
257
258  inline G4double LowEnergyLimit() const;
259
260  inline G4double PolarAngleLimit() const;
261
262  inline G4double SecondaryThreshold() const;
263
264  inline G4bool LPMFlag() const;
265
266  inline G4bool DeexcitationFlag() const;
267
268  inline void SetHighEnergyLimit(G4double);
269
270  inline void SetLowEnergyLimit(G4double);
271
272  inline void SetActivationHighEnergyLimit(G4double);
273
274  inline void SetActivationLowEnergyLimit(G4double);
275
276  inline G4bool IsActive(G4double kinEnergy);
277
278  inline void SetPolarAngleLimit(G4double);
279
280  inline void SetSecondaryThreshold(G4double);
281
282  inline void SetLPMFlag(G4bool val);
283
284  inline void SetDeexcitationFlag(G4bool val);
285
286  inline void ActivateNuclearStopping(G4bool);
287
288  inline G4double MaxSecondaryKinEnergy(const G4DynamicParticle* dynParticle);
289
290  inline const G4String& GetName() const;
291
292  inline void SetParticleChange(G4VParticleChange*, G4VEmFluctuationModel*);
293
294  inline void SetCurrentCouple(const G4MaterialCutsCouple*);
295
296  inline const G4Element* GetCurrentElement() const;
297
298protected:
299
300  inline const G4MaterialCutsCouple* CurrentCouple() const;
301
302  inline void SetCurrentElement(const G4Element*);
303
304private:
305
306  //  hide assignment operator
307  G4VEmModel & operator=(const  G4VEmModel &right);
308  G4VEmModel(const  G4VEmModel&);
309
310  // ======== Parameters of the class fixed at construction =========
311
312  G4VEmFluctuationModel* fluc;
313  const G4String   name;
314
315  // ======== Parameters of the class fixed at initialisation =======
316
317  G4double        lowLimit;
318  G4double        highLimit;
319  G4double        eMinActive;
320  G4double        eMaxActive;
321  G4double        polarAngleLimit;
322  G4double        secondaryThreshold;
323  G4bool          theLPMflag;
324
325  G4int           nSelectors;
326  std::vector<G4EmElementSelector*> elmSelectors;
327
328protected:
329
330  G4VParticleChange*  pParticleChange;
331  G4bool              nuclearStopping;
332
333  // ======== Cashed values - may be state dependent ================
334
335private:
336
337  const G4MaterialCutsCouple* currentCouple;
338  const G4Element*            currentElement;
339
340  G4int                  nsec;
341  G4bool                 flagDeexcitation;
342  std::vector<G4double>  xsec;
343
344};
345
346// ======== Run time inline methods ================
347
348inline void G4VEmModel::SetCurrentCouple(const G4MaterialCutsCouple* p)
349{
350  currentCouple = p;
351}
352
353//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
354
355inline const G4MaterialCutsCouple* G4VEmModel::CurrentCouple() const
356{
357  return currentCouple;
358}
359
360//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
361
362inline void G4VEmModel::SetCurrentElement(const G4Element* elm)
363{
364  currentElement = elm;
365}
366
367//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
368
369inline const G4Element* G4VEmModel::GetCurrentElement() const
370{
371  return currentElement;
372}
373
374//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
375
376inline 
377G4double G4VEmModel::MaxSecondaryKinEnergy(const G4DynamicParticle* dynPart)
378{
379  return MaxSecondaryEnergy(dynPart->GetDefinition(),
380                            dynPart->GetKineticEnergy());
381}
382
383//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
384
385inline G4double G4VEmModel::ComputeDEDX(const G4MaterialCutsCouple* c,
386                                        const G4ParticleDefinition* p,
387                                        G4double kinEnergy,
388                                        G4double cutEnergy)
389{
390  currentCouple = c;
391  return ComputeDEDXPerVolume(c->GetMaterial(),p,kinEnergy,cutEnergy);
392}
393
394//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
395
396inline G4double G4VEmModel::CrossSection(const G4MaterialCutsCouple* c,
397                                         const G4ParticleDefinition* p,
398                                         G4double kinEnergy,
399                                         G4double cutEnergy,
400                                         G4double maxEnergy)
401{
402  currentCouple = c;
403  return CrossSectionPerVolume(c->GetMaterial(),p,kinEnergy,cutEnergy,maxEnergy);
404}
405
406//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
407
408inline G4double G4VEmModel::ComputeMeanFreePath(const G4ParticleDefinition* p,
409                                                G4double ekin,
410                                                const G4Material* material,     
411                                                G4double emin,
412                                                G4double emax)
413{
414  G4double mfp = DBL_MAX;
415  G4double cross = CrossSectionPerVolume(material,p,ekin,emin,emax);
416  if (cross > DBL_MIN) { mfp = 1./cross; }
417  return mfp;
418}
419
420//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
421
422inline G4double G4VEmModel::ComputeCrossSectionPerAtom(
423                const G4ParticleDefinition* part,
424                const G4Element* elm,
425                G4double kinEnergy, 
426                G4double cutEnergy,
427                G4double maxEnergy)
428{
429  currentElement = elm;
430  return ComputeCrossSectionPerAtom(part,kinEnergy,elm->GetZ(),elm->GetN(),
431                                    cutEnergy,maxEnergy);
432}
433
434//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
435
436inline 
437const G4Element* G4VEmModel::SelectRandomAtom(const G4MaterialCutsCouple* couple,
438                                              const G4ParticleDefinition* p,
439                                              G4double kinEnergy,
440                                              G4double cutEnergy,
441                                              G4double maxEnergy)
442{
443  currentCouple = couple;
444  if(nSelectors > 0) {
445    currentElement = 
446      elmSelectors[couple->GetIndex()]->SelectRandomAtom(kinEnergy);
447  } else {
448    currentElement = SelectRandomAtom(couple->GetMaterial(),p,kinEnergy,
449                                      cutEnergy,maxEnergy);
450  }
451  return currentElement;
452}
453
454//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
455
456inline G4int G4VEmModel::SelectIsotopeNumber(const G4Element* elm)
457{
458  currentElement = elm;
459  G4int N = G4int(elm->GetN() + 0.5);
460  G4int ni = elm->GetNumberOfIsotopes();
461  if(ni > 0) {
462    G4int idx = 0;
463    if(ni > 1) {
464      G4double* ab = elm->GetRelativeAbundanceVector();
465      G4double x = G4UniformRand();
466      for(; idx<ni; ++idx) {
467        x -= ab[idx];
468        if (x <= 0.0) { break; }
469      }
470      if(idx >= ni) { idx = ni - 1; }
471    }
472    N = elm->GetIsotope(idx)->GetN();
473  }
474  return N;
475}
476
477// ======== Get/Set inline methods used at initialisation ================
478
479inline G4VEmFluctuationModel* G4VEmModel::GetModelOfFluctuations()
480{
481  return fluc;
482}
483
484//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
485
486inline G4double G4VEmModel::HighEnergyLimit() const
487{
488  return highLimit;
489}
490
491//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
492
493inline G4double G4VEmModel::LowEnergyLimit() const
494{
495  return lowLimit;
496}
497
498//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
499
500inline G4double G4VEmModel::PolarAngleLimit() const
501{
502  return polarAngleLimit;
503}
504
505//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
506
507inline G4double G4VEmModel::SecondaryThreshold() const
508{
509  return secondaryThreshold;
510}
511
512//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
513
514inline G4bool G4VEmModel::LPMFlag() const 
515{
516  return theLPMflag;
517}
518
519//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
520
521inline G4bool G4VEmModel::DeexcitationFlag() const 
522{
523  return flagDeexcitation;
524}
525
526//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
527
528inline void G4VEmModel::SetHighEnergyLimit(G4double val)
529{
530  highLimit = val;
531}
532
533//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
534
535inline void G4VEmModel::SetLowEnergyLimit(G4double val)
536{
537  lowLimit = val;
538}
539
540//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
541
542inline void G4VEmModel::SetActivationHighEnergyLimit(G4double val)
543{
544  eMaxActive = val;
545}
546
547//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
548
549inline void G4VEmModel::SetActivationLowEnergyLimit(G4double val)
550{
551  eMinActive = val;
552}
553
554//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
555
556inline G4bool G4VEmModel::IsActive(G4double kinEnergy)
557{
558  return (kinEnergy >= eMinActive && kinEnergy <= eMaxActive);
559}
560
561//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
562
563inline void G4VEmModel::SetPolarAngleLimit(G4double val)
564{
565  polarAngleLimit = val;
566}
567
568//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
569
570inline void G4VEmModel::SetSecondaryThreshold(G4double val) 
571{
572  secondaryThreshold = val;
573}
574
575//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
576
577inline void G4VEmModel::SetLPMFlag(G4bool val) 
578{
579  theLPMflag = val;
580}
581
582//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
583
584inline void G4VEmModel::SetDeexcitationFlag(G4bool val) 
585{
586  flagDeexcitation = val;
587}
588
589//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
590
591inline void G4VEmModel::ActivateNuclearStopping(G4bool val)
592{
593  nuclearStopping = val;
594}
595
596//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
597
598inline const G4String& G4VEmModel::GetName() const 
599{
600  return name;
601}
602
603//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
604
605inline void G4VEmModel::SetParticleChange(G4VParticleChange* p, 
606                                          G4VEmFluctuationModel* f = 0)
607{
608  if(p && pParticleChange != p) { pParticleChange = p; }
609  fluc = f;
610}
611
612//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
613
614#endif
615
Note: See TracBrowser for help on using the repository browser.