Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/physics_lists/constructors/electromagnetic/include/G4GammaGeneralProcess.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 //
 29 // GEANT4 Class header file
 30 //
 31 //
 32 // File name:     G4GammaGeneralProcess
 33 //
 34 // Author:        Vladimir Ivanchenko
 35 //
 36 // Creation date: 19.07.2018
 37 //
 38 // Modifications:
 39 //
 40 // Class Description:
 41 //
 42 // It is the gamma super process
 43 
 44 // -------------------------------------------------------------------
 45 //
 46 
 47 #ifndef G4GammaGeneralProcess_h
 48 #define G4GammaGeneralProcess_h 1
 49 
 50 #include "G4VEmProcess.hh"
 51 #include "globals.hh"
 52 #include "G4EmDataHandler.hh"
 53 
 54 class G4Step;
 55 class G4Track;
 56 class G4ParticleDefinition;
 57 class G4VParticleChange;
 58 class G4GammaConversionToMuons;
 59 class G4HadronicProcess;
 60 class G4MaterialCutsCouple;
 61 
 62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 63 
 64 class G4GammaGeneralProcess : public G4VEmProcess
 65 {
 66 public:
 67 
 68   explicit G4GammaGeneralProcess(const G4String& pname="GammaGeneralProc");
 69 
 70   ~G4GammaGeneralProcess() override;
 71 
 72   G4bool IsApplicable(const G4ParticleDefinition&) override;
 73 
 74   void AddEmProcess(G4VEmProcess*);
 75 
 76   void AddMMProcess(G4GammaConversionToMuons*);
 77 
 78   void AddHadProcess(G4HadronicProcess*);
 79 
 80   void ProcessDescription(std::ostream& outFile) const override;
 81 
 82 protected:
 83 
 84   void InitialiseProcess(const G4ParticleDefinition*) override;
 85 
 86 public:
 87 
 88   // Initialise for build of tables
 89   void PreparePhysicsTable(const G4ParticleDefinition&) override;
 90 
 91   // Build physics table during initialisation
 92   void BuildPhysicsTable(const G4ParticleDefinition&) override;
 93 
 94   // Called before tracking of each new G4Track
 95   void StartTracking(G4Track*) override;
 96 
 97   // implementation of virtual method, specific for G4GammaGeneralProcess
 98   G4double PostStepGetPhysicalInteractionLength(
 99                              const G4Track& track,
100                              G4double   previousStepSize,
101                              G4ForceCondition* condition) override;
102 
103   // implementation of virtual method, specific for G4GammaGeneralProcess
104   G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&) override;
105 
106   // Store PhysicsTable in a file.
107   // Return false in case of failure at I/O
108   G4bool StorePhysicsTable(const G4ParticleDefinition*,
109                            const G4String& directory,
110                            G4bool ascii = false) override;
111 
112   // Retrieve Physics from a file.
113   // (return true if the Physics Table can be build by using file)
114   // (return false if the process has no functionality or in case of failure)
115   // File name should is constructed as processName+particleName and the
116   // should be placed under the directory specifed by the argument.
117   G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
118                               const G4String& directory,
119                               G4bool ascii) override;
120 
121   // Return sub-process limiting current step
122   const G4VProcess* GetCreatorProcess() const override;
123   inline const G4VProcess* GetSelectedProcess() const;
124 
125   // Temporary methods
126   const G4String& GetSubProcessName() const;
127   G4int GetSubProcessSubType() const;
128 
129   G4VEmProcess* GetEmProcess(const G4String& name) override;
130 
131   inline G4HadronicProcess* GetGammaNuclear() const;
132 
133   // hide copy constructor and assignment operator
134   G4GammaGeneralProcess(G4GammaGeneralProcess &) = delete;
135   G4GammaGeneralProcess & operator=
136   (const G4GammaGeneralProcess &right) = delete;
137 
138 protected:
139 
140   G4double GetMeanFreePath(const G4Track& track, G4double previousStepSize,
141                            G4ForceCondition* condition) override;
142 
143   inline G4double ComputeGeneralLambda(std::size_t idxe, std::size_t idxt);
144 
145   inline G4double GetProbability(std::size_t idxt);
146 
147   inline void SelectedProcess(const G4Step& step, G4VProcess* ptr);
148 
149   inline void SelectEmProcess(const G4Step&, G4VEmProcess*);
150 
151   void SelectHadProcess(const G4Track&, const G4Step&, G4HadronicProcess*);
152 
153   // It returns the cross section per volume for energy/ material
154   G4double TotalCrossSectionPerVolume();
155 
156 private:
157 
158   G4bool RetrieveTable(G4VEmProcess*, const G4String& directory,
159                        G4bool ascii);
160 
161 protected:
162 
163   G4HadronicProcess*           theGammaNuclear = nullptr;
164   G4VProcess*                  selectedProc = nullptr;
165 
166   G4double                     preStepLogE = 1.0;
167   G4double                     factor = 1.0;
168 
169 
170 private:
171   static G4EmDataHandler*      theHandler;
172   static const std::size_t     nTables = 15;
173   static G4bool                theT[nTables];
174   static G4String              nameT[nTables];
175 
176   G4VEmProcess*                thePhotoElectric = nullptr;
177   G4VEmProcess*                theCompton = nullptr;
178   G4VEmProcess*                theConversionEE = nullptr;
179   G4VEmProcess*                theRayleigh = nullptr;
180   G4GammaConversionToMuons*    theConversionMM = nullptr;
181 
182   G4double                     minPEEnergy;
183   G4double                     minEEEnergy;
184   G4double                     minMMEnergy;
185   G4double                     peLambda = 0.0;
186 
187   std::size_t                  nLowE = 40;
188   std::size_t                  nHighE = 50;
189   std::size_t                  idxEnergy = 0;
190 };
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
193 
194 inline G4double
195 G4GammaGeneralProcess::ComputeGeneralLambda(std::size_t idxe, std::size_t idxt)
196 {
197   idxEnergy = idxe;
198   return factor*theHandler->GetVector(idxt, basedCoupleIndex)
199     ->LogVectorValue(preStepKinEnergy, preStepLogE);
200 }
201 
202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
203 
204 inline G4double G4GammaGeneralProcess::GetProbability(std::size_t idxt)
205 {
206   return theHandler->GetVector(idxt, basedCoupleIndex)
207     ->LogVectorValue(preStepKinEnergy, preStepLogE);
208 }
209 
210 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
211 
212 inline void
213 G4GammaGeneralProcess::SelectedProcess(const G4Step& step, G4VProcess* ptr)
214 {
215   selectedProc = ptr;
216   step.GetPostStepPoint()->SetProcessDefinedStep(ptr);
217 }
218 
219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
220 
221 inline void
222 G4GammaGeneralProcess::SelectEmProcess(const G4Step& step, G4VEmProcess* proc)
223 {
224   proc->CurrentSetup(currentCouple, preStepKinEnergy);
225   SelectedProcess(step, proc);
226 }
227 
228 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
229 
230 inline const G4VProcess* G4GammaGeneralProcess::GetSelectedProcess() const
231 {
232   return selectedProc;
233 }
234 
235 
236 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
237 
238 inline G4HadronicProcess* G4GammaGeneralProcess::GetGammaNuclear() const
239 {
240   return theGammaNuclear;
241 }
242 
243 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
244 
245 #endif
246