Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/management/include/G4HadronicProcess.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 // G4HadronicProcess
 32 //
 33 // This is the top level Hadronic Process class
 34 // The inelastic, elastic, capture, and fission processes
 35 // should derive from this class
 36 //
 37 // original by H.P.Wellisch
 38 // J.L. Chuma, TRIUMF, 10-Mar-1997
 39 // Last modified: 04-Apr-1997
 40 // 19-May-2008 V.Ivanchenko cleanup and added comments
 41 // 05-Jul-2010 V.Ivanchenko cleanup commented lines
 42 // 28-Jul-2012 M.Maire add function GetTargetDefinition() 
 43 // 14-Sep-2012 Inherit from RestDiscrete, use subtype code (now in ctor) to
 44 //    configure base-class
 45 // 28-Sep-2012 M. Kelsey -- Undo inheritance change, keep new ctor
 46 
 47 #ifndef G4HadronicProcess_h
 48 #define G4HadronicProcess_h 1
 49  
 50 #include "globals.hh"
 51 #include "G4VDiscreteProcess.hh"
 52 #include "G4EnergyRangeManager.hh"
 53 #include "G4Nucleus.hh" 
 54 #include "G4ReactionProduct.hh"
 55 #include "G4HadronicProcessType.hh"
 56 #include "G4CrossSectionDataStore.hh"
 57 #include "G4Material.hh"
 58 #include "G4DynamicParticle.hh"
 59 #include "G4ThreeVector.hh"
 60 #include "G4HadXSTypes.hh"
 61 #include <vector>
 62 
 63 class G4Track;
 64 class G4Step;
 65 class G4Element;
 66 class G4ParticleChange;
 67 class G4HadronicInteraction;
 68 class G4HadronicProcessStore;
 69 class G4VCrossSectionDataSet;
 70 class G4VLeadingParticleBiasing;
 71 class G4ParticleDefinition;
 72 
 73 class G4HadronicProcess : public G4VDiscreteProcess
 74 {
 75 public:
 76   G4HadronicProcess(const G4String& processName="Hadronic",
 77         G4ProcessType procType=fHadronic);    
 78 
 79   // Preferred signature for subclasses, specifying their subtype here
 80   G4HadronicProcess(const G4String& processName, 
 81         G4HadronicProcessType subType);    
 82 
 83   ~G4HadronicProcess() override;
 84 
 85   // register generator of secondaries
 86   void RegisterMe(G4HadronicInteraction* a);
 87 
 88   // get cross section per element
 89   G4double GetElementCrossSection(const G4DynamicParticle * part, 
 90           const G4Element * elm, 
 91           const G4Material* mat = nullptr);
 92 
 93   // obsolete method to get cross section per element
 94   inline
 95   G4double GetMicroscopicCrossSection(const G4DynamicParticle * part, 
 96               const G4Element * elm, 
 97               const G4Material* mat = nullptr);
 98 
 99   // initialisation for a new track
100   void StartTracking(G4Track* track) override;
101 
102   // compute step limit
103   G4double PostStepGetPhysicalInteractionLength(const G4Track& track,
104            G4double, G4ForceCondition*) override;
105 
106   // generic PostStepDoIt recommended for all derived classes
107   G4VParticleChange* PostStepDoIt(const G4Track& aTrack, 
108           const G4Step& aStep) override;
109 
110   // initialisation of physics tables and G4HadronicProcessStore
111   void PreparePhysicsTable(const G4ParticleDefinition&) override;
112 
113   // build physics tables and print out the configuration of the process
114   void BuildPhysicsTable(const G4ParticleDefinition&) override;
115 
116   // dump physics tables 
117   void DumpPhysicsTable(const G4ParticleDefinition& p);
118 
119   // add cross section data set
120   void AddDataSet(G4VCrossSectionDataSet * aDataSet);
121 
122   // access to the list of hadronic interactions
123   std::vector<G4HadronicInteraction*>& GetHadronicInteractionList();
124 
125   // access to an hadronic interaction by name
126   G4HadronicInteraction* GetHadronicModel(const G4String&);
127 
128   // access to the chosen generator
129   inline G4HadronicInteraction* GetHadronicInteraction() const;
130   
131   // get inverse cross section per volume
132   G4double GetMeanFreePath(const G4Track &aTrack, G4double, 
133          G4ForceCondition *) override;
134 
135   // access to the target nucleus
136   inline const G4Nucleus* GetTargetNucleus() const;
137 
138   inline G4Nucleus* GetTargetNucleusPointer();
139   
140   inline const G4Isotope* GetTargetIsotope();
141 
142   // methods needed for implementation of integral XS
143   G4double ComputeCrossSection(const G4ParticleDefinition*,
144                                const G4Material*,
145                                const G4double kinEnergy);
146 
147   inline G4HadXSType CrossSectionType() const;
148   inline void SetCrossSectionType(G4HadXSType val);
149 
150   void ProcessDescription(std::ostream& outFile) const override;
151  
152   // scale cross section
153   void BiasCrossSectionByFactor(G4double aScale);
154   void MultiplyCrossSectionBy(G4double factor);
155   inline G4double CrossSectionFactor() const;
156 
157   // Integral option 
158   inline void SetIntegral(G4bool val);
159 
160   // Energy-momentum non-conservation limits and reporting
161   inline void SetEpReportLevel(G4int level);
162   inline void SetEnergyMomentumCheckLevels(G4double relativeLevel,
163                                            G4double absoluteLevel);
164   inline std::pair<G4double, G4double> GetEnergyMomentumCheckLevels() const;
165 
166   // access to the cross section data store
167   inline G4CrossSectionDataStore* GetCrossSectionDataStore();
168 
169   // access to the data for integral XS method
170   inline std::vector<G4TwoPeaksHadXS*>* TwoPeaksXS() const;
171   inline std::vector<G4double>* EnergyOfCrossSectionMax() const;
172 
173   // hide assignment operator as private 
174   G4HadronicProcess& operator=(const G4HadronicProcess& right) = delete;
175   G4HadronicProcess(const G4HadronicProcess&) = delete;
176 
177 protected:
178 
179   // generic method to choose secondary generator 
180   // recommended for all derived classes
181   inline G4HadronicInteraction* ChooseHadronicInteraction(
182       const G4HadProjectile & aHadProjectile, G4Nucleus& aTargetNucleus,
183       const G4Material* aMaterial, const G4Element* anElement);
184               
185   // access to the cross section data set
186   inline G4double GetLastCrossSection();
187 
188   // fill result
189   void FillResult(G4HadFinalState* aR, const G4Track& aT);
190 
191   void DumpState(const G4Track&, const G4String&, G4ExceptionDescription&);
192 
193   // Check the result for catastrophic energy non-conservation
194   G4HadFinalState* CheckResult(const G4HadProjectile& thePro,
195              const G4Nucleus& targetNucleus, 
196              G4HadFinalState* result);
197 
198   // Check 4-momentum balance
199   void CheckEnergyMomentumConservation(const G4Track&, const G4Nucleus&);
200 
201 private:
202 
203   void InitialiseLocal();
204   void UpdateCrossSectionAndMFP(const G4double kinEnergy);
205   void RecomputeXSandMFP(const G4double kinEnergy);
206 
207   inline void DefineXSandMFP();
208   inline void ComputeXSandMFP();
209 
210   G4double XBiasSurvivalProbability();
211   G4double XBiasSecondaryWeight();
212 
213 protected:
214 
215   G4HadProjectile thePro;
216 
217   G4ParticleChange* theTotalResult;
218   G4CrossSectionDataStore* theCrossSectionDataStore;
219 
220   G4double fWeight = 1.0;
221   G4double aScaleFactor = 1.0;
222   G4double theLastCrossSection = 0.0;
223   G4double mfpKinEnergy = DBL_MAX;
224   G4int epReportLevel = 0;
225 
226   G4HadXSType fXSType = fHadNoIntegral;
227 
228 private:
229     
230   G4EnergyRangeManager theEnergyRangeManager;
231   G4Nucleus targetNucleus;
232     
233   G4HadronicInteraction* theInteraction = nullptr;
234   G4HadronicProcessStore* theProcessStore;
235   const G4HadronicProcess* masterProcess = nullptr;
236   const G4ParticleDefinition* firstParticle = nullptr;
237   const G4ParticleDefinition* currentParticle = nullptr;
238   const G4Material* currentMat = nullptr;
239   const G4DynamicParticle* fDynParticle = nullptr;
240 
241   std::vector<G4double>* theEnergyOfCrossSectionMax = nullptr;
242   std::vector<G4TwoPeaksHadXS*>* fXSpeaks = nullptr;
243 
244   G4double theMFP = DBL_MAX;
245   G4double minKinEnergy;
246 
247   // counters
248   G4int nMatWarn = 0;
249   G4int nKaonWarn = 0;
250   G4int nICelectrons = 0;
251   G4int matIdx = 0;
252 
253   // flags
254   G4bool levelsSetByProcess = false;
255   G4bool useIntegralXS = true;
256   G4bool isMaster = true;
257 
258   G4ThreeVector unitVector;
259 
260   // Energy-momentum checking
261   std::pair<G4double, G4double> epCheckLevels;
262   std::vector<G4VLeadingParticleBiasing*> theBias;
263 };
264 
265 inline G4double G4HadronicProcess::
266 GetMicroscopicCrossSection(const G4DynamicParticle * part, 
267          const G4Element * elm, 
268          const G4Material* mat)
269 { 
270   return GetElementCrossSection(part, elm, mat);
271 }
272 
273 inline G4HadronicInteraction* 
274 G4HadronicProcess::GetHadronicInteraction() const
275 { 
276   return theInteraction;
277 }
278 
279 inline const G4Nucleus*
280 G4HadronicProcess::GetTargetNucleus() const
281 {
282   return &targetNucleus;
283 }
284   
285 inline const G4Isotope* G4HadronicProcess::GetTargetIsotope()
286 { 
287   return targetNucleus.GetIsotope();
288 }
289 
290 inline G4HadXSType 
291 G4HadronicProcess::CrossSectionType() const
292 { 
293   return fXSType;
294 }
295 
296 inline void 
297 G4HadronicProcess::SetCrossSectionType(G4HadXSType val)
298 { 
299   fXSType = val;
300 }
301 
302 inline G4double G4HadronicProcess::CrossSectionFactor() const
303 { 
304   return aScaleFactor;
305 }
306 
307 inline void G4HadronicProcess::SetIntegral(G4bool val)
308 { 
309   useIntegralXS = val;
310 }
311 
312 inline void G4HadronicProcess::SetEpReportLevel(G4int level)
313 { 
314   epReportLevel = level;
315 }
316 
317 inline void 
318 G4HadronicProcess::SetEnergyMomentumCheckLevels(G4double relativeLevel, 
319                                                 G4double absoluteLevel)
320 { 
321   epCheckLevels.first = relativeLevel;
322   epCheckLevels.second = absoluteLevel;
323   levelsSetByProcess = true;
324 }
325 
326 inline std::pair<G4double, G4double>
327 G4HadronicProcess::GetEnergyMomentumCheckLevels() const
328 { 
329   return epCheckLevels;
330 }
331 
332 inline G4CrossSectionDataStore*
333 G4HadronicProcess::GetCrossSectionDataStore()
334 {
335   return theCrossSectionDataStore;
336 }
337 
338 inline std::vector<G4TwoPeaksHadXS*>* 
339 G4HadronicProcess::TwoPeaksXS() const
340 { 
341   return fXSpeaks;
342 }
343 
344 inline std::vector<G4double>*
345 G4HadronicProcess::EnergyOfCrossSectionMax() const
346 {
347   return theEnergyOfCrossSectionMax;
348 }
349 
350 inline G4HadronicInteraction* G4HadronicProcess::
351 ChooseHadronicInteraction(const G4HadProjectile& aHadProjectile,
352                           G4Nucleus& aTargetNucleus,
353                           const G4Material* aMaterial,
354                           const G4Element* anElement)
355 { 
356   return theEnergyRangeManager.GetHadronicInteraction(aHadProjectile, 
357                                                       aTargetNucleus,
358                                                       aMaterial,anElement);
359 }
360 
361 inline G4Nucleus* G4HadronicProcess::GetTargetNucleusPointer() 
362 { 
363   return &targetNucleus;
364 }
365 
366 inline G4double G4HadronicProcess::GetLastCrossSection() 
367 { 
368   return theLastCrossSection;
369 }
370 
371 inline void G4HadronicProcess::DefineXSandMFP()
372 {
373   theLastCrossSection = aScaleFactor*
374     theCrossSectionDataStore->GetCrossSection(fDynParticle, currentMat);
375   theMFP = (theLastCrossSection > 0.0) ? 1.0/theLastCrossSection : DBL_MAX;
376 }
377 
378 inline void G4HadronicProcess::ComputeXSandMFP()
379 {
380   theLastCrossSection = aScaleFactor*
381     theCrossSectionDataStore->ComputeCrossSection(fDynParticle, currentMat);
382   theMFP = (theLastCrossSection > 0.0) ? 1.0/theLastCrossSection : DBL_MAX;
383 }
384  
385 #endif
386  
387