Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/hadronic/Hadr09/include/HadronicGenerator.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 /// \file HadronicGenerator.hh
 27 /// \brief Definition of the HadronicGenerator class
 28 //
 29 //------------------------------------------------------------------------
 30 // Class: HadronicGenerator
 31 // Author: Alberto Ribon (CERN EP/SFT)
 32 // Date: May 2020
 33 //
 34 // This class shows how to use Geant4 as a generator for simulating
 35 // inelastic hadron-nuclear interactions.
 36 // Some of the most used hadronic models are currently supported in
 37 // this class:
 38 // - the hadronic string models Fritiof (FTF) and Quark-Gluon-String (QGS)
 39 //   coupled with Precompound/de-excitation
 40 // - the intranuclear cascade models: Bertini (BERT), Binary Cascade (BIC),
 41 //                                    and Liege (INCL)
 42 // Combinations of two models - in a transition energy interval, with a
 43 // linear probability as a function of the energy - are also available to
 44 // "mimic" the transition between hadronic models as in the most common
 45 // Geant4 reference physics lists.
 46 //
 47 // The current version of this class does NOT support:
 48 // -  hadron elastic interactions
 49 // -  neutron capture and fission
 50 // -  precise low-energy inelastic interactions of neutrons and
 51 //    charged particles (i.e. ParticleHP)
 52 // -  gamma/lepton-nuclear inelastic interactions
 53 //
 54 // This class does NOT use the Geant4 run-manager, and therefore should
 55 // be usable in a multi-threaded application, with one instance of this
 56 // class in each thread.
 57 //
 58 // This class has been inspired by test30 (whose author is Vladimir
 59 // Ivanchenko), with various simplifications and restricted to hadronic
 60 // inelastic interactions.
 61 //------------------------------------------------------------------------
 62 
 63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 65 
 66 #ifndef HadronicGenerator_h
 67 #define HadronicGenerator_h 1
 68 
 69 #include "G4HadronicProcess.hh"
 70 #include "G4ThreeVector.hh"
 71 #include "G4ios.hh"
 72 #include "globals.hh"
 73 
 74 #include <iomanip>
 75 #include <map>
 76 
 77 class G4ParticleDefinition;
 78 class G4VParticleChange;
 79 class G4ParticleTable;
 80 class G4Material;
 81 class G4HadronicInteraction;
 82 
 83 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 84 
 85 class HadronicGenerator
 86 {
 87     // This class provides the functionality of a "hadronic generator"
 88     // for Geant4 final-state inelastic hadronic collisions.
 89     // Only a few of the available Geant4 final-state hadronic inelastic
 90     // "physics cases" are currently available in this class - but it can
 91     // be extended to other cases if needed.
 92     // It is important to notice that this class does NOT use the Geant4
 93     // run-manager, so it should work fine in a multi-threaded environment,
 94     // with a separate instance of this class in each thread.
 95   public:
 96     explicit HadronicGenerator(const G4String physicsCase = "FTFP_BERT_ATL");
 97     // Currently supported final-state hadronic inelastic "physics cases":
 98     // -  Hadronic models :        BERT, BIC, IonBIC, INCL, FTFP, QGSP
 99     // -  "Physics-list proxies" : FTFP_BERT_ATL (default), FTFP_BERT,
100     //                             QGSP_BERT, QGSP_BIC, FTFP_INCLXX
101     //    (i.e. they are not real, complete physics lists - for instance
102     //     they do not have: transportation, electromagnetic physics,
103     //     hadron elastic scattering, neutron fission and capture, etc. -
104     //     however, they cover all hadron types and all energies by
105     //     combining different hadronic models, i.e. there are transitions
106     //     between two hadronic models in well-defined energy intervals,
107     //     e.g. "FTFP_BERT" has the transition between BERT and FTFP
108     //     hadronic models; moreover, the transition intervals used in
109     //     our "physics cases"might not be the same as in the corresponding
110     //     physics lists).
111 
112     ~HadronicGenerator();
113 
114     inline G4bool IsPhysicsCaseSupported() const;
115     // Returns "true" if the physicsCase is supported; "false" otherwise.
116 
117     G4bool IsApplicable(const G4String& nameProjectile, const G4double projectileEnergy) const;
118     G4bool IsApplicable(G4ParticleDefinition* projectileDefinition,
119                         const G4double projectileEnergy) const;
120     // Returns "true" if the specified projectile (either by name or particle definition)
121     // of given energy is applicable, "false" otherwise.
122 
123     G4VParticleChange* GenerateInteraction(const G4String& nameProjectile,
124                                            const G4double projectileEnergy,
125                                            const G4ThreeVector& projectileDirection,
126                                            G4Material* targetMaterial);
127     G4VParticleChange* GenerateInteraction(G4ParticleDefinition* projectileDefinition,
128                                            const G4double projectileEnergy,
129                                            const G4ThreeVector& projectileDirection,
130                                            G4Material* targetMaterial);
131     // This is the main method provided by the class:
132     // in input it receives the projectile (either by name or particle definition),
133     // its energy, its direction and the target material, and it returns one sampled
134     // final-state of the inelastic hadron-nuclear collision as modelled by the
135     // final-state hadronic inelastic "physics case" specified in the constructor.
136     // If the required hadronic collision is not possible, then the method returns
137     // immediately an empty "G4VParticleChange", i.e. without secondaries produced.
138 
139     inline G4HadronicProcess* GetHadronicProcess() const;
140     inline G4HadronicInteraction* GetHadronicInteraction() const;
141     // Returns the hadronic process and the hadronic interaction, respectively,
142     // that handled the last call of "GenerateInteraction".
143 
144     G4double GetImpactParameter() const;
145     G4int GetNumberOfTargetSpectatorNucleons() const;
146     G4int GetNumberOfProjectileSpectatorNucleons() const;
147     G4int GetNumberOfNNcollisions() const;
148     // In the case of hadronic interactions handled by the FTF model, returns,
149     // respectively, the impact parameter, the number of target/projectile
150     // spectator nucleons, and the number of nucleon-nucleon collisions,
151     // else, returns a negative value (-999).
152 
153   private:
154     G4String fPhysicsCase;
155     G4bool fPhysicsCaseIsSupported;
156     G4HadronicProcess* fLastHadronicProcess;
157     G4ParticleTable* fPartTable;
158     std::map<G4ParticleDefinition*, G4HadronicProcess*> fProcessMap;
159 };
160 
161 inline G4bool HadronicGenerator::IsPhysicsCaseSupported() const
162 {
163   return fPhysicsCaseIsSupported;
164 }
165 
166 inline G4HadronicProcess* HadronicGenerator::GetHadronicProcess() const
167 {
168   return fLastHadronicProcess;
169 }
170 
171 inline G4HadronicInteraction* HadronicGenerator::GetHadronicInteraction() const
172 {
173   return fLastHadronicProcess == nullptr ? nullptr : fLastHadronicProcess->GetHadronicInteraction();
174 }
175 
176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177 
178 #endif
179