Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:21

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 // -----------------------------------------------------------------------------
0028 //
0029 // GEANT4 Class header file
0030 //
0031 // File name:     G4GoudsmitSaundersonTable
0032 //
0033 // Author:        Mihaly Novak / (Omrane Kadri)
0034 //
0035 // Creation date: 20.02.2009
0036 //
0037 // Class description:
0038 //   Class to handle multiple scattering angular distributions precomputed by
0039 //   using Kawrakow-Bielajew Goudsmit-Saunderson MSC model based on the screened
0040 //   Rutherford DCS for elastic scattering of electrons/positrons [1,2]. This
0041 //   class is used by G4GoudsmitSaundersonMscModel to sample the angular
0042 //   deflection of electrons/positrons after travelling a given path.
0043 //
0044 // Modifications:
0045 // 04.03.2009 V.Ivanchenko cleanup and format according to Geant4 EM style
0046 // 18.05.2015 M. Novak This class has been completely replaced (only the original
0047 //            class name was kept; class description was also inserted):
0048 //            A new version of Kawrakow-Bielajew Goudsmit-Saunderson MSC model
0049 //            based on the screened Rutherford DCS for elastic scattering of
0050 //            electrons/positrons has been introduced[1,2]. The corresponding MSC
0051 //            angular distributions over a 2D parameter grid have been recomputed
0052 //            and the CDFs are now stored in a variable transformed (smooth) form
0053 //            together with the corresponding rational interpolation parameters.
0054 //            The new version is several times faster, more robust and accurate
0055 //            compared to the earlier version (G4GoudsmitSaundersonMscModel class
0056 //            that use these data has been also completely replaced)
0057 // 28.04.2017 M. Novak: the GS angular distributions has been recomputed, the
0058 //            data size has been reduced from 16 MB down to 5 MB by using a new
0059 //            representation, the class has been modified significantly due to
0060 //            this new data representation.
0061 // 23.08.2017 M. Novak: Added funtionality to handle Mott-correction to the
0062 //            base GS angular distributions and some other factors (screening
0063 //            parameter, first and second moments) when Mott-correction is
0064 //            activated in the GS-MSC model.
0065 //
0066 // References:
0067 //   [1] A.F.Bielajew, NIMB, 111 (1996) 195-208
0068 //   [2] I.Kawrakow, A.F.Bielajew, NIMB 134(1998) 325-336
0069 //
0070 // -----------------------------------------------------------------------------
0071 
0072 
0073 #ifndef G4GoudsmitSaundersonTable_h
0074 #define G4GoudsmitSaundersonTable_h 1
0075 
0076 #include <vector>
0077 
0078 #include "G4Types.hh"
0079 
0080 class G4GSMottCorrection;
0081 class G4MaterialCutsCouple;
0082 
0083 class G4GoudsmitSaundersonTable {
0084 
0085 public:
0086   G4GoudsmitSaundersonTable(G4bool iselectron);
0087  ~G4GoudsmitSaundersonTable();
0088 
0089   void Initialise(G4double lownergylimit, G4double highenergylimit);
0090 
0091   // structure to store one GS transformed angular distribution (for a given s/lambda_el,s/lambda_elG1)
0092   struct GSMSCAngularDtr {
0093     G4int     fNumData;    // # of data points
0094     G4double *fUValues;    // array of transformed variables
0095     G4double *fParamA;     // array of interpolation parameters a
0096     G4double *fParamB;     // array of interpolation parameters b
0097   };
0098 
0099   void   LoadMSCData();
0100 
0101   G4bool   Sampling(G4double lambdaval, G4double qval,  G4double scra,
0102                     G4double &cost,     G4double &sint, G4double lekin,
0103                     G4double beta2,     G4int matindx, GSMSCAngularDtr **gsDtr,
0104                     G4int &mcekini, G4int &mcdelti, G4double &transfPar,
0105                     G4bool isfirst);
0106 
0107   G4double SampleCosTheta(G4double lambdaval, G4double qval,  G4double scra,
0108                           G4double lekin,     G4double beta2, G4int matindx,
0109                           GSMSCAngularDtr **gsDtr, G4int &mcekini, G4int &mcdelti,
0110                           G4double &transfPar, G4bool isfirst);
0111 
0112   G4double SampleGSSRCosTheta(const GSMSCAngularDtr* gsDrt, G4double transfpar);
0113 
0114   G4double SingleScattering(G4double lambdaval, G4double scra, G4double lekin,
0115                             G4double beta2, G4int matindx);
0116 
0117   GSMSCAngularDtr* GetGSAngularDtr(G4double scra, G4double &lambdaval,
0118                                    G4double &qval, G4double &transfpar);
0119 
0120   // material dependent MSC parameters (computed at initialisation) regarding
0121   // Moliere's screening parameter
0122   G4double GetMoliereBc(G4int matindx)  { return gMoliereBc[matindx];  }
0123 
0124   G4double GetMoliereXc2(G4int matindx) { return gMoliereXc2[matindx]; }
0125 
0126   void     GetMottCorrectionFactors(G4double logekin, G4double beta2,
0127                                     G4int matindx, G4double &mcToScr,
0128                                     G4double &mcToQ1, G4double &mcToG2PerG1);
0129 
0130   // set option to activate/inactivate Mott-correction
0131   void     SetOptionMottCorrection(G4bool val) { fIsMottCorrection = val; }
0132   // set option to activate/inactivate PWA-correction
0133   void     SetOptionPWACorrection(G4bool val) { fIsPWACorrection = val; }
0134 
0135   // this method returns with the scattering power correction (to avoid double counting of sub-threshold deflections)
0136   // interpolated from tables prepared at initialisation
0137   G4double ComputeScatteringPowerCorrection(const G4MaterialCutsCouple *matcut, G4double ekin);
0138 
0139   void     InitSCPCorrection();
0140 
0141 private:
0142   // initialisation of material dependent Moliere's MSC parameters
0143   void InitMoliereMSCParams();
0144 
0145 
0146  private:
0147    static G4bool             gIsInitialised;       // are the precomputed angular distributions already loaded in?
0148    static constexpr G4int    gLAMBNUM = 64;        // # L=s/lambda_el in [fLAMBMIN,fLAMBMAX]
0149    static constexpr G4int    gQNUM1   = 15;        // # Q=s/lambda_el G1 in [fQMIN1,fQMAX1] in the 1-st Q grid
0150    static constexpr G4int    gQNUM2   = 32;        // # Q=s/lambda_el G1 in [fQMIN2,fQMAX2] in the 2-nd Q grid
0151    static constexpr G4int    gNUMSCR1 = 201;       // # of screening parameters in the A(G1) function
0152    static constexpr G4int    gNUMSCR2 = 51;        // # of screening parameters in the A(G1) function
0153    static constexpr G4double gLAMBMIN = 1.0;       // minimum s/lambda_el
0154    static constexpr G4double gLAMBMAX = 100000.0;  // maximum s/lambda_el
0155    static constexpr G4double gQMIN1   = 0.001;     // minimum s/lambda_el G1 in the 1-st Q grid
0156    static constexpr G4double gQMAX1   = 0.99;      // maximum s/lambda_el G1 in the 1-st Q grid
0157    static constexpr G4double gQMIN2   = 0.99;      // minimum s/lambda_el G1 in the 2-nd Q grid
0158    static constexpr G4double gQMAX2   = 7.99;      // maximum s/lambda_el G1 in the 2-nd Q grid
0159    //
0160    G4bool   fIsElectron;          // GS-table for e- (for e+ otherwise)
0161    G4bool   fIsMottCorrection;    // flag to indicate if Mott-correction was requested to be used
0162    G4bool   fIsPWACorrection;     // flag to indicate is PWA corrections were requested to be used
0163    G4double fLogLambda0;          // ln(gLAMBMIN)
0164    G4double fLogDeltaLambda;      // ln(gLAMBMAX/gLAMBMIN)/(gLAMBNUM-1)
0165    G4double fInvLogDeltaLambda;   // 1/[ln(gLAMBMAX/gLAMBMIN)/(gLAMBNUM-1)]
0166    G4double fInvDeltaQ1;          // 1/[(gQMAX1-gQMIN1)/(gQNUM1-1)]
0167    G4double fDeltaQ2;             // [(gQMAX2-gQMIN2)/(gQNUM2-1)]
0168    G4double fInvDeltaQ2;          // 1/[(gQMAX2-gQMIN2)/(gQNUM2-1)]
0169    //
0170    G4double fLowEnergyLimit;
0171    G4double fHighEnergyLimit;
0172    //
0173    int      fNumSPCEbinPerDec;    // scattering power correction energy grid bins per decade
0174    struct SCPCorrection {
0175      bool   fIsUse;               //
0176      double fPrCut;               // sec. e- production cut energy
0177      double fLEmin;               // log min energy
0178      double fILDel;               // inverse log delta kinetic energy
0179      //std::vector<double> fVEkin;  // scattering power correction energies
0180      std::vector<double> fVSCPC;  // scattering power correction vector
0181    };
0182    std::vector<SCPCorrection*>  fSCPCPerMatCuts;
0183 
0184 
0185    // vector to store all GS transformed angular distributions (cumputed based on the Screened-Rutherford DCS)
0186    static std::vector<GSMSCAngularDtr*> gGSMSCAngularDistributions1;
0187    static std::vector<GSMSCAngularDtr*> gGSMSCAngularDistributions2;
0188 
0189    //@{
0190    /** Precomputed \f$ b_lambda_{c} $\f and \f$ \chi_c^{2} $\f material dependent
0191    *   Moliere parameters that can be used to compute the screening parameter,
0192    *   the elastic scattering cross section (or \f$ \lambda_{e} $\f) under the
0193    *   screened Rutherford cross section approximation. (These are used in
0194    *   G4GoudsmitSaundersonMscModel if fgIsUsePWATotalXsecData is FALSE.)
0195    */
0196    static std::vector<double> gMoliereBc;
0197    static std::vector<double> gMoliereXc2;
0198    //
0199    //
0200    G4GSMottCorrection   *fMottCorrection;
0201 };
0202 
0203 #endif