Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/utils/include/G4EmExtraParameters.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 // GEANT4 Class header file
 29 //
 30 // File name:     G4EmExtraParameters
 31 //
 32 // Author:        Vladimir Ivanchenko
 33 //                  
 34 // Creation date: 06.05.2019
 35 //
 36 // Class Description:
 37 //
 38 // An internal utility class, responsable for keeping parameters
 39 // for EM processes and models.
 40 //
 41 // It is initialized by the master thread but can be updated 
 42 // at any moment via G4EmParameters interface. It is not assumed
 43 // to be used for a direct initialisation
 44 //
 45 // -------------------------------------------------------------------
 46 //
 47 
 48 #ifndef G4EmExtraParameters_h
 49 #define G4EmExtraParameters_h 1
 50 
 51 #include "globals.hh"
 52 #include "G4ios.hh"
 53 #include "G4ThreeVector.hh"
 54 #include <vector>
 55 
 56 class G4EmExtraParametersMessenger;
 57 class G4VEnergyLossProcess;
 58 class G4VEmProcess;
 59 class G4ParticleDefinition;
 60 class G4VAtomDeexcitation;
 61 
 62 class G4EmExtraParameters
 63 {
 64 public:
 65 
 66   explicit G4EmExtraParameters();
 67 
 68   ~G4EmExtraParameters();
 69 
 70   void Initialise();
 71 
 72   G4bool GetDirectionalSplitting();
 73   void SetDirectionalSplitting(G4bool v);
 74 
 75   G4bool QuantumEntanglement();
 76   void SetQuantumEntanglement(G4bool v);
 77 
 78   void SetDirectionalSplittingRadius(G4double r);
 79   G4double GetDirectionalSplittingRadius();
 80 
 81   void SetDirectionalSplittingTarget(const G4ThreeVector& v);
 82   G4ThreeVector GetDirectionalSplittingTarget() const;
 83 
 84   void SetStepFunction(G4double v1, G4double v2);
 85   G4double GetStepFunctionP1() const;
 86   G4double GetStepFunctionP2() const;
 87  
 88   void SetStepFunctionMuHad(G4double v1, G4double v2);
 89   G4double GetStepFunctionMuHadP1() const;
 90   G4double GetStepFunctionMuHadP2() const;
 91 
 92   void SetStepFunctionLightIons(G4double v1, G4double v2);
 93   G4double GetStepFunctionLightIonsP1() const;
 94   G4double GetStepFunctionLightIonsP2() const;
 95 
 96   void SetStepFunctionIons(G4double v1, G4double v2);
 97   G4double GetStepFunctionIonsP1() const;
 98   G4double GetStepFunctionIonsP2() const;
 99 
100   void FillStepFunction(const G4ParticleDefinition*, G4VEnergyLossProcess*) const;
101 
102   // parameters per region or per process 
103   void AddPAIModel(const G4String& particle,
104                    const G4String& region,
105                    const G4String& type);
106   const std::vector<G4String>& ParticlesPAI() const;
107   const std::vector<G4String>& RegionsPAI() const;
108   const std::vector<G4String>& TypesPAI() const;
109 
110   void AddPhysics(const G4String& region, const G4String& type);
111   const std::vector<G4String>& RegionsPhysics() const;
112   const std::vector<G4String>& TypesPhysics() const;
113 
114   void SetSubCutRegion(const G4String& region);
115 
116   void SetProcessBiasingFactor(const G4String& procname, 
117                                G4double val, G4bool wflag);
118 
119   void ActivateForcedInteraction(const G4String& procname, 
120                                  const G4String& region,
121                                  G4double length, 
122                                  G4bool wflag);
123 
124   void ActivateSecondaryBiasing(const G4String& name,
125         const G4String& region, 
126         G4double factor,
127         G4double energyLimit);
128 
129   // initialisation methods
130   void DefineRegParamForLoss(G4VEnergyLossProcess*) const;
131   void DefineRegParamForEM(G4VEmProcess*) const;
132 
133   G4EmExtraParameters(G4EmExtraParameters &) = delete;
134   G4EmExtraParameters & operator=
135   (const G4EmExtraParameters &right) = delete;  
136 
137 private:
138 
139   G4String CheckRegion(const G4String&) const;
140 
141   void PrintWarning(G4ExceptionDescription& ed) const;
142 
143   G4EmExtraParametersMessenger* theMessenger;
144 
145   G4bool directionalSplitting;
146   G4bool quantumEntanglement;
147 
148   G4double dRoverRange;
149   G4double finalRange;
150   G4double dRoverRangeMuHad;
151   G4double finalRangeMuHad;
152   G4double dRoverRangeLIons;
153   G4double finalRangeLIons;
154   G4double dRoverRangeIons;
155   G4double finalRangeIons;
156 
157   G4double directionalSplittingRadius;
158   G4ThreeVector directionalSplittingTarget;
159 
160   std::vector<G4String>  m_particlesPAI;
161   std::vector<G4String>  m_regnamesPAI;
162   std::vector<G4String>  m_typesPAI;
163 
164   std::vector<G4String>  m_regnamesPhys;
165   std::vector<G4String>  m_typesPhys;
166 
167   std::vector<G4String>  m_regnamesSubCut;
168 
169   std::vector<G4String>  m_procBiasedXS;
170   std::vector<G4double>  m_factBiasedXS;
171   std::vector<G4bool>    m_weightBiasedXS;
172 
173   std::vector<G4String>  m_procForced;
174   std::vector<G4String>  m_regnamesForced;
175   std::vector<G4double>  m_lengthForced;
176   std::vector<G4bool>    m_weightForced;
177 
178   std::vector<G4String>  m_procBiasedSec;
179   std::vector<G4String>  m_regnamesBiasedSec;
180   std::vector<G4double>  m_factBiasedSec;
181   std::vector<G4double>  m_elimBiasedSec;
182 };
183 
184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
185 
186 #endif
187