Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/include/G4eDPWACoulombScatteringModel.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 29 // GEANT4 Class header file
 30 //
 31 //
 32 // File name:     G4eDPWACoulombScatteringModel
 33 //
 34 // Author:        Mihaly Novak
 35 //
 36 // Creation date: 02.07.2020
 37 //
 38 // Modifications:
 39 //
 40 // Class Description:
 41 //
 42 // e-/e+ Coulomb scattering model based on numerical Differential Cross Sections
 43 // (DCS) obtained by Dirac Partial Wave Analysis (DPWA) and supplied by the
 44 // G4eDPWAElasticDCS class.
 45 // The model contains the possibility to incorporate the effects of angular
 46 // deflections of sub-threshold ionisation intercations when it's described by
 47 // the condensed history model. Note, this must be inactivated (by setting the
 48 // `isscpcor` input argument of the CTR to false) when ionisation is described
 49 // with a classical, event by event based simulation model instead of usign the
 50 // condensed history approach (otherwise, the corresponding angular defelctions
 51 // will be "double counted").
 52 //
 53 // -------------------------------------------------------------------
 54 
 55 
 56 
 57 #ifndef G4eDPWACoulombScatteringModel_h
 58 #define G4eDPWACoulombScatteringModel_h 1
 59 
 60 #include "G4VEmModel.hh"
 61 #include "globals.hh"
 62 
 63 class G4eDPWAElasticDCS;
 64 class G4ParticleChangeForGamma;
 65 class G4ParticleDefinition;
 66 class G4DataVector;
 67 
 68 class G4eDPWACoulombScatteringModel : public G4VEmModel {
 69 
 70 public:
 71 
 72   /**
 73    * Constructor.
 74    *
 75    * @param[in] ismixed  Indicates if the model is for mixed or for pure single
 76    *                     Coulomb scattering. Different type of tables are pre-
 77    *                     pared for sampling polar angle of Coulomb scattering
 78    *                     for mixed and for pure single scattering models: cosine
 79    *                     of the polar scattering angle can be sampled in a
 80    *                     restriced inteval (see mumin input parameter below).
 81    * @param[in] isscpcor Indicates if scattering power correction should be used.
 82    *                     Note, scattering power correction accounts the effects
 83    *                     angular deflections due to sub-threshold ionisations
 84    *                     when ionisation is described by using condensed history
 85    *                     model (should be active only in this case).
 86    * @param[in] mumin    When the model is used for mixed simulation, Coulomb
 87    *                     scatterings, resulting in a minimum t_c polar angular
 88    *                     deflection, modelled explicitly. Therefore, cross
 89    *                     sections are computed, and angular deflections are
 90    *                     sampled ina resricted [\theta_c,\pi] interval. The
 91    *                     minimum of this interval is determined by the mumin
 92    *                     parameter as:
 93    *                     \mu_{min} = \mu(\theta_c)=0.5[1-\cos(\theta_c)]
 94    */
 95   G4eDPWACoulombScatteringModel(G4bool ismixed=false, G4bool isscpcor=true,
 96                                 G4double mumin=0.0);
 97 
 98   ~G4eDPWACoulombScatteringModel() override;
 99 
100   //
101   // Interface methods:
102 
103   void     Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
104 
105   void     InitialiseLocal(const G4ParticleDefinition*, G4VEmModel*) override;
106 
107   G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition*, G4double ekin,
108                                       G4double Z, G4double A, G4double prodcut,
109                                       G4double emax) override;
110 
111   void     SampleSecondaries(std::vector<G4DynamicParticle*>*,
112                              const G4MaterialCutsCouple*,
113                              const G4DynamicParticle*,
114                              G4double tmin,
115                              G4double maxEnergy) override;
116 
117   G4double MinPrimaryEnergy(const G4Material*, const G4ParticleDefinition*, 
118                             G4double) override { return 10.0*CLHEP::eV; }
119 
120   void     SetTheDCS(G4eDPWAElasticDCS* theDCS) { fTheDCS = theDCS; }
121 
122   G4eDPWAElasticDCS* GetTheDCS() { return fTheDCS; }
123 
124 
125 private:
126 
127   // Indicates if the model is mixed: MSC for soft (theta<theta_c), Singe
128   // Scattering(SS) for hard scatterings(theta>theta_c). SS otherwise.
129   // Note, that while the model provides restricted (elastic and transport)
130   // cross sections, it's responsible to handle, i.e. provide final state,
131   // only for the Singe Scattering part in case of a mixed model.
132   G4bool                     fIsMixedModel;
133   // indicates if scattering power correction should be applied: correction due
134   // to deflection in case of sub-threshold, inelastic interactions -> only in
135   // case of condensed history simulation of inonisation!
136   G4bool                     fIsScpCorrection;
137   // mu(theta)=0.5[1-cos(theta)]: the model porvides final states \in [fMuMin,1]
138   G4double                   fMuMin;
139   // the object that provides cross sections and polar angle of scattering
140   G4eDPWAElasticDCS*         fTheDCS;
141   // particle change
142   G4ParticleChangeForGamma*  fParticleChange;
143 
144 };
145 
146 #endif
147