Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/adjoint/include/G4AdjointCSManager.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 //  Class:    G4AdjointCSManager
 28 //  Author:         L. Desorgher
 29 //  Organisation:   SpaceIT GmbH
 30 //
 31 // Class is responsible for the management of all adjoint cross section
 32 // matrices, and for the computation of the total forward and adjoint cross
 33 // sections. Total adjoint and forward cross sections are needed to correct the
 34 // weight of a particle after a tracking step or after the occurrence of a
 35 // reverse reaction. It is also used to sample an adjoint secondary from a
 36 // given adjoint cross section matrix.
 37 //
 38 ////////////////////////////////////////////////////////////////////////////////
 39 
 40 #ifndef G4AdjointCSManager_h
 41 #define G4AdjointCSManager_h 1
 42 
 43 #include "globals.hh"
 44 #include "G4AdjointCSMatrix.hh"
 45 #include "G4ThreadLocalSingleton.hh"
 46 
 47 #include <vector>
 48 
 49 class G4Element;
 50 class G4Material;
 51 class G4MaterialCutsCouple;
 52 class G4ParticleDefinition;
 53 class G4PhysicsTable;
 54 class G4VEmProcess;
 55 class G4VEmAdjointModel;
 56 class G4VEnergyLossProcess;
 57 
 58 class G4AdjointCSManager
 59 {
 60   friend class G4ThreadLocalSingleton<G4AdjointCSManager>;
 61 
 62  public:
 63   ~G4AdjointCSManager();
 64   static G4AdjointCSManager* GetAdjointCSManager();
 65 
 66   G4int GetNbProcesses();
 67 
 68   // Registration of the different models and processes
 69 
 70   std::size_t RegisterEmAdjointModel(G4VEmAdjointModel*);
 71 
 72   void RegisterEmProcess(G4VEmProcess* aProcess,
 73                          G4ParticleDefinition* aPartDef);
 74 
 75   void RegisterEnergyLossProcess(G4VEnergyLossProcess* aProcess,
 76                                  G4ParticleDefinition* aPartDef);
 77 
 78   void RegisterAdjointParticle(G4ParticleDefinition* aPartDef);
 79 
 80   // Building of the CS Matrices and Total Forward and Adjoint LambdaTables
 81   void BuildCrossSectionMatrices();
 82 
 83   void BuildTotalSigmaTables();
 84 
 85   // Get TotalCrossSections form Total Lambda Tables, Needed for Weight
 86   // correction and scaling of the
 87   G4double GetTotalAdjointCS(G4ParticleDefinition* aPartDef, G4double Ekin,
 88                              const G4MaterialCutsCouple* aCouple);
 89 
 90   G4double GetTotalForwardCS(G4ParticleDefinition* aPartDef, G4double Ekin,
 91                              const G4MaterialCutsCouple* aCouple);
 92 
 93   G4double GetAdjointSigma(G4double Ekin_nuc, std::size_t index_model,
 94                            G4bool is_scat_proj_to_proj,
 95                            const G4MaterialCutsCouple* aCouple);
 96 
 97   void GetEminForTotalCS(G4ParticleDefinition* aPartDef,
 98                          const G4MaterialCutsCouple* aCouple,
 99                          G4double& emin_adj, G4double& emin_fwd);
100 
101   void GetMaxFwdTotalCS(G4ParticleDefinition* aPartDef,
102                         const G4MaterialCutsCouple* aCouple,
103                         G4double& e_sigma_max, G4double& sigma_max);
104 
105   void GetMaxAdjTotalCS(G4ParticleDefinition* aPartDef,
106                         const G4MaterialCutsCouple* aCouple,
107                         G4double& e_sigma_max, G4double& sigma_max);
108 
109   // CrossSection Correction 1 or FwdCS/AdjCS following the G4boolean value of
110   // forward_CS_is_used and forward_CS_mode
111   G4double GetCrossSectionCorrection(G4ParticleDefinition* aPartDef,
112                                      G4double PreStepEkin,
113                                      const G4MaterialCutsCouple* aCouple,
114                                      G4bool& fwd_is_used);
115 
116   // Cross section mode
117   inline void SetFwdCrossSectionMode(G4bool aBool) { fForwardCSMode = aBool; }
118 
119   // Weight correction
120   G4double GetContinuousWeightCorrection(G4ParticleDefinition* aPartDef,
121                                          G4double PreStepEkin,
122                                          G4double AfterStepEkin,
123                                          const G4MaterialCutsCouple* aCouple,
124                                          G4double step_length);
125 
126   G4double GetPostStepWeightCorrection();
127 
128   // called by the adjoint model to get the CS, if not otherwise specified
129   G4double ComputeAdjointCS(G4Material* aMaterial, G4VEmAdjointModel* aModel,
130                             G4double PrimEnergy, G4double Tcut,
131                             G4bool isScatProjToProj,
132                             std::vector<G4double>& AdjointCS_for_each_element);
133 
134   // called by the adjoint model to sample secondary energy from the CS matrix
135   G4Element* SampleElementFromCSMatrices(G4Material* aMaterial,
136                                          G4VEmAdjointModel* aModel,
137                                          G4double PrimEnergy, G4double Tcut,
138                                          G4bool isScatProjToProj);
139 
140   // Total Adjoint CS is computed at initialisation phase
141   G4double ComputeTotalAdjointCS(const G4MaterialCutsCouple* aMatCutCouple,
142                                  G4ParticleDefinition* aPart,
143                                  G4double PrimEnergy);
144 
145   G4ParticleDefinition* GetAdjointParticleEquivalent(
146     G4ParticleDefinition* theFwdPartDef);
147 
148   G4ParticleDefinition* GetForwardParticleEquivalent(
149     G4ParticleDefinition* theAdjPartDef);
150 
151   // inline
152   inline void SetIon(G4ParticleDefinition* adjIon, G4ParticleDefinition* fwdIon)
153   {
154     fAdjIon = adjIon;
155     fFwdIon = fwdIon;
156   }
157 
158  private:
159   G4AdjointCSManager();
160 
161   void DefineCurrentMaterial(const G4MaterialCutsCouple* couple);
162 
163   void DefineCurrentParticle(const G4ParticleDefinition* aPartDef);
164 
165   G4double ComputeAdjointCS(G4double aPrimEnergy,
166                             G4AdjointCSMatrix* anAdjointCSMatrix,
167                             G4double Tcut);
168 
169   std::vector<G4AdjointCSMatrix*> BuildCrossSectionsModelAndElement(
170     G4VEmAdjointModel* aModel, G4int Z, G4int A, G4int nbin_pro_decade);
171 
172   std::vector<G4AdjointCSMatrix*> BuildCrossSectionsModelAndMaterial(
173     G4VEmAdjointModel* aModel, G4Material* aMaterial, G4int nbin_pro_decade);
174 
175   static constexpr G4double fTmin = 0.1 * CLHEP::keV;
176   static constexpr G4double fTmax = 100. * CLHEP::TeV;
177   // fNbins chosen to avoid error
178   // in the CS value close to CS jump. (For example at Tcut)
179   static constexpr G4int fNbins = 320;
180 
181   static G4ThreadLocal G4AdjointCSManager* fInstance;
182 
183   // only one ion can be considered by simulation
184   G4ParticleDefinition* fAdjIon = nullptr;
185   G4ParticleDefinition* fFwdIon = nullptr;
186 
187   G4MaterialCutsCouple* fCurrentCouple = nullptr;
188   G4Material* fCurrentMaterial         = nullptr;
189 
190   // x dim is for G4VAdjointEM*, y dim is for elements
191   std::vector<std::vector<G4AdjointCSMatrix*>>
192     fAdjointCSMatricesForScatProjToProj;
193 
194   std::vector<std::vector<G4AdjointCSMatrix*>> fAdjointCSMatricesForProdToProj;
195 
196   std::vector<G4VEmAdjointModel*> fAdjointModels;
197 
198   std::vector<std::size_t> fIndexOfAdjointEMModelInAction;
199   std::vector<G4bool> fIsScatProjToProj;
200   std::vector<std::vector<G4double>> fLastAdjointCSVsModelsAndElements;
201 
202   // total adjoint and total forward cross section table in function of material
203   // and in function of adjoint particle type
204   std::vector<G4PhysicsTable*> fTotalFwdSigmaTable;
205   std::vector<G4PhysicsTable*> fTotalAdjSigmaTable;
206 
207   // Sigma table for each G4VAdjointEMModel
208   std::vector<G4PhysicsTable*> fSigmaTableForAdjointModelScatProjToProj;
209   std::vector<G4PhysicsTable*> fSigmaTableForAdjointModelProdToProj;
210 
211   std::vector<std::vector<G4double>> fEminForFwdSigmaTables;
212   std::vector<std::vector<G4double>> fEminForAdjSigmaTables;
213   std::vector<std::vector<G4double>> fEkinofFwdSigmaMax;
214   std::vector<std::vector<G4double>> fEkinofAdjSigmaMax;
215 
216   // list of forward G4VEmProcess and of G4VEnergyLossProcess for the different
217   // adjoint particle
218   std::vector<std::vector<G4VEmProcess*>*> fForwardProcesses;
219   std::vector<std::vector<G4VEnergyLossProcess*>*> fForwardLossProcesses;
220 
221   // list of adjoint particles considered
222   std::vector<G4ParticleDefinition*> fAdjointParticlesInAction;
223 
224   G4double fMassRatio              = 1.;  // ion
225   G4double fLastCSCorrectionFactor = 1.;
226 
227   G4ParticleDefinition* fCurrentParticleDef = nullptr;
228   std::size_t fCurrentParticleIndex = 0;
229   std::size_t fCurrentMatIndex      = 0;
230 
231   G4bool fCSMatricesBuilt = false;
232   G4bool fSigmaTableBuilt = false;
233   G4bool fForwardCSUsed   = true;
234   G4bool fForwardCSMode   = true;
235   // Two CS mode are possible:
236   // 1) fForwardCSMode = false, the Adjoint CS are used as it is implying
237   //    an AlongStep Weight Correction.
238   // 2) fForwardCSMode = true, the Adjoint CS are scaled to have the total
239   //    adjoint CS equal to the fwd one implying a PostStep Weight Correction.
240   // For energies where the total Fwd CS or the total adjoint CS are zero,
241   // the scaling is not possible and fForwardCSUsed is set to false
242 };
243 #endif
244