Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/adjoint/include/G4VEmAdjointModel.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:  G4VEMAdjointModel
 28 //  Author:         L. Desorgher
 29 //  Organisation:   SpaceIT GmbH
 30 //
 31 //  Base class for Adjoint EM model. It is based on the use of direct
 32 //  G4VEmModel.
 33 ////////////////////////////////////////////////////////////////////////////////
 34 
 35 #ifndef G4VEmAdjointModel_h
 36 #define G4VEmAdjointModel_h 1
 37 
 38 #include "globals.hh"
 39 #include "G4ParticleDefinition.hh"
 40 #include "G4VEmModel.hh"
 41 
 42 class G4AdjointCSMatrix;
 43 class G4AdjointCSManager;
 44 class G4Material;
 45 class G4MaterialCutsCouple;
 46 class G4ParticleChange;
 47 class G4Region;
 48 class G4Track;
 49 
 50 class G4VEmAdjointModel
 51 {
 52  public:
 53   explicit G4VEmAdjointModel(const G4String& nam);
 54 
 55   virtual ~G4VEmAdjointModel();
 56 
 57   //------------------------------------------------------------------------
 58   // Virtual methods to be implemented for the sample secondaries concrete model
 59   //------------------------------------------------------------------------
 60 
 61   virtual void SampleSecondaries(const G4Track& aTrack, G4bool isScatProjToProj,
 62                                  G4ParticleChange* fParticleChange) = 0;
 63 
 64   //------------------------------------------------------------------------
 65   // Methods for adjoint processes
 66   //------------------------------------------------------------------------
 67 
 68   virtual G4double AdjointCrossSection(const G4MaterialCutsCouple* aCouple,
 69                                        G4double primEnergy,
 70                                        G4bool isScatProjToProj);
 71 
 72   // The implementation of the DiffCrossSection... here are correct for
 73   // energy loss process. For the photoelectric and Compton scattering
 74   // the method should be redefined
 75   virtual G4double DiffCrossSectionPerAtomPrimToSecond(
 76     G4double kinEnergyProj,  // kin energy of primary before interaction
 77     G4double kinEnergyProd,  // kinetic energy of the secondary particle
 78     G4double Z, G4double A = 0.);
 79 
 80   virtual G4double DiffCrossSectionPerAtomPrimToScatPrim(
 81     G4double kinEnergyProj,      // kin energy of primary before interaction
 82     G4double kinEnergyScatProj,  // kin energy of primary after interaction
 83     G4double Z, G4double A = 0.);
 84 
 85   virtual G4double DiffCrossSectionPerVolumePrimToSecond(
 86     const G4Material* aMaterial,
 87     G4double kinEnergyProj,  // kin energy of primary before interaction
 88     G4double kinEnergyProd   // kinetic energy of secondary particle
 89   );
 90 
 91   virtual G4double DiffCrossSectionPerVolumePrimToScatPrim(
 92     const G4Material* aMaterial,
 93     G4double kinEnergyProj,     // kin energy of primary before interaction
 94     G4double kinEnergyScatProj  // kinetic energy of primary after interaction
 95   );
 96 
 97   // Energy limits of adjoint secondary
 98   //------------------
 99 
100   virtual G4double GetSecondAdjEnergyMaxForScatProjToProj(
101     G4double primAdjEnergy);
102 
103   virtual G4double GetSecondAdjEnergyMinForScatProjToProj(
104     G4double primAdjEnergy, G4double tcut = 0.);
105 
106   virtual G4double GetSecondAdjEnergyMaxForProdToProj(G4double primAdjEnergy);
107 
108   virtual G4double GetSecondAdjEnergyMinForProdToProj(G4double primAdjEnergy);
109 
110   // Other Methods
111   //---------------
112 
113   void DefineCurrentMaterial(const G4MaterialCutsCouple* couple);
114 
115   std::vector<std::vector<double>*>
116   ComputeAdjointCrossSectionVectorPerAtomForSecond(G4double kinEnergyProd,
117                                                    G4double Z, G4double A = 0.,
118                                                    G4int nbin_pro_decade = 10);
119 
120   std::vector<std::vector<double>*>
121   ComputeAdjointCrossSectionVectorPerAtomForScatProj(
122     G4double kinEnergyProd, G4double Z, G4double A = 0.,
123     G4int nbin_pro_decade = 10);
124 
125   std::vector<std::vector<double>*>
126   ComputeAdjointCrossSectionVectorPerVolumeForSecond(
127     G4Material* aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade = 10);
128 
129   std::vector<std::vector<double>*>
130   ComputeAdjointCrossSectionVectorPerVolumeForScatProj(
131     G4Material* aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade = 10);
132 
133   inline void SetCSMatrices(std::vector<G4AdjointCSMatrix*>* Vec1CSMatrix,
134                             std::vector<G4AdjointCSMatrix*>* Vec2CSMatrix)
135   {
136     fCSMatrixProdToProjBackScat = Vec1CSMatrix;
137     fCSMatrixProjToProjBackScat = Vec2CSMatrix;
138   };
139 
140   inline G4ParticleDefinition*
141   GetAdjointEquivalentOfDirectPrimaryParticleDefinition() const
142   {
143     return fAdjEquivDirectPrimPart;
144   }
145 
146   inline G4ParticleDefinition*
147   GetAdjointEquivalentOfDirectSecondaryParticleDefinition() const
148   {
149     return fAdjEquivDirectSecondPart;
150   }
151 
152   inline G4double GetHighEnergyLimit() const { return fHighEnergyLimit; }
153 
154   inline G4double GetLowEnergyLimit() const { return fLowEnergyLimit; }
155 
156   void SetHighEnergyLimit(G4double aVal);
157 
158   void SetLowEnergyLimit(G4double aVal);
159 
160   inline void DefineDirectEMModel(G4VEmModel* aModel) { fDirectModel = aModel; }
161 
162   void SetAdjointEquivalentOfDirectPrimaryParticleDefinition(
163     G4ParticleDefinition* aPart);
164 
165   inline void SetAdjointEquivalentOfDirectSecondaryParticleDefinition(
166     G4ParticleDefinition* aPart)
167   {
168     fAdjEquivDirectSecondPart = aPart;
169   }
170 
171   inline void SetSecondPartOfSameType(G4bool aBool)
172   {
173     fSecondPartSameType = aBool;
174   }
175 
176   inline G4bool GetSecondPartOfSameType() const { return fSecondPartSameType; }
177 
178   inline void SetUseMatrix(G4bool aBool) { fUseMatrix = aBool; }
179 
180   inline void SetUseMatrixPerElement(G4bool aBool)
181   {
182     fUseMatrixPerElement = aBool;
183   }
184 
185   inline void SetUseOnlyOneMatrixForAllElements(G4bool aBool)
186   {
187     fOneMatrixForAllElements = aBool;
188   }
189 
190   inline void SetApplyCutInRange(G4bool aBool) { fApplyCutInRange = aBool; }
191 
192   inline G4bool GetUseMatrix() const { return fUseMatrix; }
193 
194   inline G4bool GetUseMatrixPerElement() const { return fUseMatrixPerElement; }
195 
196   inline G4bool GetUseOnlyOneMatrixForAllElements() const
197   {
198     return fOneMatrixForAllElements;
199   }
200 
201   inline G4bool GetApplyCutInRange() const { return fApplyCutInRange; }
202 
203   inline const G4String& GetName() const { return fName; }
204 
205   inline virtual void SetCSBiasingFactor(G4double aVal)
206   {
207     fCsBiasingFactor = aVal;
208   }
209 
210   inline void SetCorrectWeightForPostStepInModel(G4bool aBool)
211   {
212     fInModelWeightCorr = aBool;
213   }
214 
215   inline void SetAdditionalWeightCorrectionFactorForPostStepOutsideModel(
216     G4double factor)
217   {
218     fOutsideWeightFactor = factor;
219   }
220 
221   G4VEmAdjointModel(G4VEmAdjointModel&) = delete;
222   G4VEmAdjointModel& operator=(const G4VEmAdjointModel& right) = delete;
223 
224  protected:
225   G4double DiffCrossSectionFunction1(G4double kinEnergyProj);
226 
227   G4double DiffCrossSectionFunction2(G4double kinEnergyProj);
228 
229   // General methods to sample secondary energy
230   G4double SampleAdjSecEnergyFromCSMatrix(std::size_t MatrixIndex,
231                                           G4double prim_energy,
232                                           G4bool isScatProjToProj);
233 
234   G4double SampleAdjSecEnergyFromCSMatrix(G4double prim_energy,
235                                           G4bool isScatProjToProj);
236 
237   void SelectCSMatrix(G4bool isScatProjToProj);
238 
239   virtual G4double SampleAdjSecEnergyFromDiffCrossSectionPerAtom(
240     G4double prim_energy, G4bool isScatProjToProj);
241 
242   // Post  Step weight correction
243   virtual void CorrectPostStepWeight(G4ParticleChange* fParticleChange,
244                                      G4double old_weight,
245                                      G4double adjointPrimKinEnergy,
246                                      G4double projectileKinEnergy,
247                                      G4bool isScatProjToProj);
248 
249   G4AdjointCSManager* fCSManager;
250   G4VEmModel* fDirectModel = nullptr;
251 
252   const G4String fName;
253 
254   G4Material* fSelectedMaterial        = nullptr;
255   G4Material* fCurrentMaterial         = nullptr;
256   G4MaterialCutsCouple* fCurrentCouple = nullptr;
257 
258   // particle definition
259   G4ParticleDefinition* fAdjEquivDirectPrimPart   = nullptr;
260   G4ParticleDefinition* fAdjEquivDirectSecondPart = nullptr;
261   G4ParticleDefinition* fDirectPrimaryPart        = nullptr;
262 
263   // adjoint CS matrix for each element or material
264   std::vector<G4AdjointCSMatrix*>* fCSMatrixProdToProjBackScat = nullptr;
265   std::vector<G4AdjointCSMatrix*>* fCSMatrixProjToProjBackScat = nullptr;
266 
267   std::vector<G4double> fElementCSScatProjToProj;
268   std::vector<G4double> fElementCSProdToProj;
269 
270   G4double fKinEnergyProdForIntegration     = 0.;
271   G4double fKinEnergyScatProjForIntegration = 0.;
272 
273   G4double fLastCS                         = 0.;
274   G4double fLastAdjointCSForScatProjToProj = 0.;
275   G4double fLastAdjointCSForProdToProj     = 0.;
276 
277   G4double fPreStepEnergy = 0.;
278 
279   G4double fTcutPrim   = 0.;
280   G4double fTcutSecond = 0.;
281 
282   // Energy limits
283   G4double fHighEnergyLimit = 0.;
284   G4double fLowEnergyLimit  = 0.;
285 
286   // Cross Section biasing factor
287   G4double fCsBiasingFactor = 1.;
288 
289   // [1] This is needed for the forced interaction where part of the weight
290   // correction is given outside the model while the secondary are created in
291   // the model. The weight should be fixed before adding the secondary
292   G4double fOutsideWeightFactor = 1.;
293 
294   // Needed for CS integration at the initialisation phase
295   G4int fASelectedNucleus = 0;
296   G4int fZSelectedNucleus = 0;
297 
298   std::size_t fCSMatrixUsed = 0;  // Index of crosssection matrices used
299 
300   G4bool fSecondPartSameType = false;
301   G4bool fInModelWeightCorr =
302     false;  // correct_weight_for_post_step_in_model, see [1]
303 
304   G4bool fApplyCutInRange = true;
305 
306   // Type of Model with Matrix or not
307   G4bool fUseMatrix               = false;
308   G4bool fUseMatrixPerElement     = false;  // other possibility is per Material
309   G4bool fOneMatrixForAllElements = false;
310 };
311 
312 #endif
313