Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:59:20

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 ////////////////////////////////////////////////////////////////////////////////
0027 //  Class:  G4VEMAdjointModel
0028 //  Author:         L. Desorgher
0029 //  Organisation:   SpaceIT GmbH
0030 //
0031 //  Base class for Adjoint EM model. It is based on the use of direct
0032 //  G4VEmModel.
0033 ////////////////////////////////////////////////////////////////////////////////
0034 
0035 #ifndef G4VEmAdjointModel_h
0036 #define G4VEmAdjointModel_h 1
0037 
0038 #include "globals.hh"
0039 #include "G4ParticleDefinition.hh"
0040 #include "G4VEmModel.hh"
0041 
0042 class G4AdjointCSMatrix;
0043 class G4AdjointCSManager;
0044 class G4Material;
0045 class G4MaterialCutsCouple;
0046 class G4ParticleChange;
0047 class G4Region;
0048 class G4Track;
0049 
0050 class G4VEmAdjointModel
0051 {
0052  public:
0053   explicit G4VEmAdjointModel(const G4String& nam);
0054 
0055   virtual ~G4VEmAdjointModel();
0056 
0057   //------------------------------------------------------------------------
0058   // Virtual methods to be implemented for the sample secondaries concrete model
0059   //------------------------------------------------------------------------
0060 
0061   virtual void SampleSecondaries(const G4Track& aTrack, G4bool isScatProjToProj,
0062                                  G4ParticleChange* fParticleChange) = 0;
0063 
0064   //------------------------------------------------------------------------
0065   // Methods for adjoint processes
0066   //------------------------------------------------------------------------
0067 
0068   virtual G4double AdjointCrossSection(const G4MaterialCutsCouple* aCouple,
0069                                        G4double primEnergy,
0070                                        G4bool isScatProjToProj);
0071 
0072   // The implementation of the DiffCrossSection... here are correct for
0073   // energy loss process. For the photoelectric and Compton scattering
0074   // the method should be redefined
0075   virtual G4double DiffCrossSectionPerAtomPrimToSecond(
0076     G4double kinEnergyProj,  // kin energy of primary before interaction
0077     G4double kinEnergyProd,  // kinetic energy of the secondary particle
0078     G4double Z, G4double A = 0.);
0079 
0080   virtual G4double DiffCrossSectionPerAtomPrimToScatPrim(
0081     G4double kinEnergyProj,      // kin energy of primary before interaction
0082     G4double kinEnergyScatProj,  // kin energy of primary after interaction
0083     G4double Z, G4double A = 0.);
0084 
0085   virtual G4double DiffCrossSectionPerVolumePrimToSecond(
0086     const G4Material* aMaterial,
0087     G4double kinEnergyProj,  // kin energy of primary before interaction
0088     G4double kinEnergyProd   // kinetic energy of secondary particle
0089   );
0090 
0091   virtual G4double DiffCrossSectionPerVolumePrimToScatPrim(
0092     const G4Material* aMaterial,
0093     G4double kinEnergyProj,     // kin energy of primary before interaction
0094     G4double kinEnergyScatProj  // kinetic energy of primary after interaction
0095   );
0096 
0097   // Energy limits of adjoint secondary
0098   //------------------
0099 
0100   virtual G4double GetSecondAdjEnergyMaxForScatProjToProj(
0101     G4double primAdjEnergy);
0102 
0103   virtual G4double GetSecondAdjEnergyMinForScatProjToProj(
0104     G4double primAdjEnergy, G4double tcut = 0.);
0105 
0106   virtual G4double GetSecondAdjEnergyMaxForProdToProj(G4double primAdjEnergy);
0107 
0108   virtual G4double GetSecondAdjEnergyMinForProdToProj(G4double primAdjEnergy);
0109 
0110   // Other Methods
0111   //---------------
0112 
0113   void DefineCurrentMaterial(const G4MaterialCutsCouple* couple);
0114 
0115   std::vector<std::vector<double>*>
0116   ComputeAdjointCrossSectionVectorPerAtomForSecond(G4double kinEnergyProd,
0117                                                    G4double Z, G4double A = 0.,
0118                                                    G4int nbin_pro_decade = 10);
0119 
0120   std::vector<std::vector<double>*>
0121   ComputeAdjointCrossSectionVectorPerAtomForScatProj(
0122     G4double kinEnergyProd, G4double Z, G4double A = 0.,
0123     G4int nbin_pro_decade = 10);
0124 
0125   std::vector<std::vector<double>*>
0126   ComputeAdjointCrossSectionVectorPerVolumeForSecond(
0127     G4Material* aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade = 10);
0128 
0129   std::vector<std::vector<double>*>
0130   ComputeAdjointCrossSectionVectorPerVolumeForScatProj(
0131     G4Material* aMaterial, G4double kinEnergyProd, G4int nbin_pro_decade = 10);
0132 
0133   inline void SetCSMatrices(std::vector<G4AdjointCSMatrix*>* Vec1CSMatrix,
0134                             std::vector<G4AdjointCSMatrix*>* Vec2CSMatrix)
0135   {
0136     fCSMatrixProdToProjBackScat = Vec1CSMatrix;
0137     fCSMatrixProjToProjBackScat = Vec2CSMatrix;
0138   };
0139 
0140   inline G4ParticleDefinition*
0141   GetAdjointEquivalentOfDirectPrimaryParticleDefinition()
0142   {
0143     return fAdjEquivDirectPrimPart;
0144   }
0145 
0146   inline G4ParticleDefinition*
0147   GetAdjointEquivalentOfDirectSecondaryParticleDefinition()
0148   {
0149     return fAdjEquivDirectSecondPart;
0150   }
0151 
0152   inline G4double GetHighEnergyLimit() { return fHighEnergyLimit; }
0153 
0154   inline G4double GetLowEnergyLimit() { return fLowEnergyLimit; }
0155 
0156   void SetHighEnergyLimit(G4double aVal);
0157 
0158   void SetLowEnergyLimit(G4double aVal);
0159 
0160   inline void DefineDirectEMModel(G4VEmModel* aModel) { fDirectModel = aModel; }
0161 
0162   void SetAdjointEquivalentOfDirectPrimaryParticleDefinition(
0163     G4ParticleDefinition* aPart);
0164 
0165   inline void SetAdjointEquivalentOfDirectSecondaryParticleDefinition(
0166     G4ParticleDefinition* aPart)
0167   {
0168     fAdjEquivDirectSecondPart = aPart;
0169   }
0170 
0171   inline void SetSecondPartOfSameType(G4bool aBool)
0172   {
0173     fSecondPartSameType = aBool;
0174   }
0175 
0176   inline G4bool GetSecondPartOfSameType() { return fSecondPartSameType; }
0177 
0178   inline void SetUseMatrix(G4bool aBool) { fUseMatrix = aBool; }
0179 
0180   inline void SetUseMatrixPerElement(G4bool aBool)
0181   {
0182     fUseMatrixPerElement = aBool;
0183   }
0184 
0185   inline void SetUseOnlyOneMatrixForAllElements(G4bool aBool)
0186   {
0187     fOneMatrixForAllElements = aBool;
0188   }
0189 
0190   inline void SetApplyCutInRange(G4bool aBool) { fApplyCutInRange = aBool; }
0191 
0192   inline G4bool GetUseMatrix() { return fUseMatrix; }
0193 
0194   inline G4bool GetUseMatrixPerElement() { return fUseMatrixPerElement; }
0195 
0196   inline G4bool GetUseOnlyOneMatrixForAllElements()
0197   {
0198     return fOneMatrixForAllElements;
0199   }
0200 
0201   inline G4bool GetApplyCutInRange() { return fApplyCutInRange; }
0202 
0203   inline G4String GetName() { return fName; }
0204 
0205   inline virtual void SetCSBiasingFactor(G4double aVal)
0206   {
0207     fCsBiasingFactor = aVal;
0208   }
0209 
0210   inline void SetCorrectWeightForPostStepInModel(G4bool aBool)
0211   {
0212     fInModelWeightCorr = aBool;
0213   }
0214 
0215   inline void SetAdditionalWeightCorrectionFactorForPostStepOutsideModel(
0216     G4double factor)
0217   {
0218     fOutsideWeightFactor = factor;
0219   }
0220 
0221   G4VEmAdjointModel(G4VEmAdjointModel&) = delete;
0222   G4VEmAdjointModel& operator=(const G4VEmAdjointModel& right) = delete;
0223 
0224  protected:
0225   G4double DiffCrossSectionFunction1(G4double kinEnergyProj);
0226 
0227   G4double DiffCrossSectionFunction2(G4double kinEnergyProj);
0228 
0229   // General methods to sample secondary energy
0230   G4double SampleAdjSecEnergyFromCSMatrix(std::size_t MatrixIndex,
0231                                           G4double prim_energy,
0232                                           G4bool isScatProjToProj);
0233 
0234   G4double SampleAdjSecEnergyFromCSMatrix(G4double prim_energy,
0235                                           G4bool isScatProjToProj);
0236 
0237   void SelectCSMatrix(G4bool isScatProjToProj);
0238 
0239   virtual G4double SampleAdjSecEnergyFromDiffCrossSectionPerAtom(
0240     G4double prim_energy, G4bool isScatProjToProj);
0241 
0242   // Post  Step weight correction
0243   virtual void CorrectPostStepWeight(G4ParticleChange* fParticleChange,
0244                                      G4double old_weight,
0245                                      G4double adjointPrimKinEnergy,
0246                                      G4double projectileKinEnergy,
0247                                      G4bool isScatProjToProj);
0248 
0249   G4AdjointCSManager* fCSManager;
0250   G4VEmModel* fDirectModel = nullptr;
0251 
0252   const G4String fName;
0253 
0254   G4Material* fSelectedMaterial        = nullptr;
0255   G4Material* fCurrentMaterial         = nullptr;
0256   G4MaterialCutsCouple* fCurrentCouple = nullptr;
0257 
0258   // particle definition
0259   G4ParticleDefinition* fAdjEquivDirectPrimPart   = nullptr;
0260   G4ParticleDefinition* fAdjEquivDirectSecondPart = nullptr;
0261   G4ParticleDefinition* fDirectPrimaryPart        = nullptr;
0262 
0263   // adjoint CS matrix for each element or material
0264   std::vector<G4AdjointCSMatrix*>* fCSMatrixProdToProjBackScat = nullptr;
0265   std::vector<G4AdjointCSMatrix*>* fCSMatrixProjToProjBackScat = nullptr;
0266 
0267   std::vector<G4double> fElementCSScatProjToProj;
0268   std::vector<G4double> fElementCSProdToProj;
0269 
0270   G4double fKinEnergyProdForIntegration     = 0.;
0271   G4double fKinEnergyScatProjForIntegration = 0.;
0272 
0273   G4double fLastCS                         = 0.;
0274   G4double fLastAdjointCSForScatProjToProj = 0.;
0275   G4double fLastAdjointCSForProdToProj     = 0.;
0276 
0277   G4double fPreStepEnergy = 0.;
0278 
0279   G4double fTcutPrim   = 0.;
0280   G4double fTcutSecond = 0.;
0281 
0282   // Energy limits
0283   G4double fHighEnergyLimit = 0.;
0284   G4double fLowEnergyLimit  = 0.;
0285 
0286   // Cross Section biasing factor
0287   G4double fCsBiasingFactor = 1.;
0288 
0289   // [1] This is needed for the forced interaction where part of the weight
0290   // correction is given outside the model while the secondary are created in
0291   // the model. The weight should be fixed before adding the secondary
0292   G4double fOutsideWeightFactor = 1.;
0293 
0294   // Needed for CS integration at the initialisation phase
0295   G4int fASelectedNucleus = 0;
0296   G4int fZSelectedNucleus = 0;
0297 
0298   std::size_t fCSMatrixUsed = 0;  // Index of crosssection matrices used
0299 
0300   G4bool fSecondPartSameType = false;
0301   G4bool fInModelWeightCorr =
0302     false;  // correct_weight_for_post_step_in_model, see [1]
0303 
0304   G4bool fApplyCutInRange = true;
0305 
0306   // Type of Model with Matrix or not
0307   G4bool fUseMatrix               = false;
0308   G4bool fUseMatrixPerElement     = false;  // other possibility is per Material
0309   G4bool fOneMatrixForAllElements = false;
0310 };
0311 
0312 #endif