Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4QEnvironment.hh
Go to the documentation of this file.
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// $Id$
28//
29// ---------------- G4QEnvironment ----------------
30// by Mikhail Kossov, August 2000.
31// header for Multy Quasmon Environment in the CHIPS Model
32// ------------------------------------------------------------
33// Short description: The G4QEnvironment class corresponds to the nuclear
34// environment, containing excited nucleons or nuclear clusters
35// (G4Quasmons). In the constructer the nucleus (G4QNucleus) is clusterized
36// and then the projectile hadron (or hadrons) can create one (or a few)
37// Quasmons, which are fragmented by the Fragment member function. As a
38// result a vector of G4QHadrons is created, which can include the residual
39// nucleus in the ground state.
40//---------------------------------------------------------------------
41
42#ifndef G4QEnvironment_h
43#define G4QEnvironment_h 1
44
45#include "G4RandomDirection.hh"
46#include "G4QuasmonVector.hh"
47#include "G4QFreeScattering.hh"
48#include "G4QThd.hh"
49
51{
52public:
53 G4QEnvironment(const G4QNucleus& theEnv); // Create Env and add Quasmons
54 // later
55 G4QEnvironment(const G4QHadronVector& projHadrons, const G4int targPDG);
56 G4QEnvironment(const G4QEnvironment& right); // copy QEnvironment by value
57 G4QEnvironment(G4QEnvironment* right); // copy QEnvironment by pointer
58 ~G4QEnvironment(); // Public Destructor
59
60 // Overloaded operators
61 const G4QEnvironment& operator=(const G4QEnvironment& right);
62 G4bool operator==(const G4QEnvironment &right) const;
63 G4bool operator!=(const G4QEnvironment &right) const;
64
65 //Selectors
67 G4QuasmonVector* GetQuasmons(); // User is responsible for Destroy/Clear/Delete
68 G4QHadronVector* GetQHadrons(); // User is responsible for Destroy/Clear/Delete
69 G4QHadronVector* GetProjectiles(); // User is responsible for Destroy/Clear/Delete
70
71 // Modifiers
72 void AddQuasmon(G4Quasmon* Q); // Add aQuasmon to theEnvironment
73 G4QHadronVector* Fragment(); // User must clear and destroy the G4QHadronVec
74
75 // External managment
76 void DecayBaryon(G4QHadron* dB, G4QHadronVector* HV); // Decay baryon
77 void DecayMeson(G4QHadron* dB, G4QHadronVector* HV); // Decay meson
78 void DecayAntistrange(G4QHadron* aS, G4QHadronVector* HV);// Decay AntiStrNuc
79 void CheckMassShell(G4QHadronVector* HV); // All hadrons on mass shell
80
81 // Static functions
82 static void SetParameters(G4double solAn=0.4,G4bool efFlag=false,G4double piThresh=141.4,
83 G4double mpisq=20000., G4double dinum=1880.);
84 static void OpenElectromagneticDecays();
85 static void CloseElectromagneticDecays();
86
87protected:
88 void CleanUpQHadrons(); // Only another G4QEnvironment issue can use it (can make mess)
89 void FillQHadrons(G4QHadronVector* input); //Only another G4QEnvironment issue can use it
90
91private:
92 G4QHadronVector* FSInteraction(); // Final State Interaction after Hadronization
93 G4QHadronVector HadronizeQEnvironment(); // Main HadronizationFunction used in Fragment
94 void CopyAndDeleteHadronVector(G4QHadronVector* HV);//Copy HadrVect to Output
95 void CreateQuasmon(const G4QContent& pQC, const G4LorentzVector& p4M, G4bool f=false);
96 void InitClustersVector(G4int maxC, G4int maxA);//Init.NucClust's for 1st int
97 void CleanUp(); // Makes theEnvironment=vacuum & kill Quasmons
98 void PrepareInteractionProbabilities(const G4QContent& projQC, G4double AP);
99 void EvaporateResidual(G4QHadron* h, G4bool f=true);// Evaporate NuclearFragm
100 G4bool CheckGroundState(G4Quasmon* quasm,G4bool corFlag=false);//as G4Q for QHV
101 G4bool DecayInEnvQ(G4Quasmon* quasm); // Use befor evaporation in PANIC case
102
103// Body
104private:
105 // Static Parameters
106 static G4bool WeakDecays; // Flag for opening WeakDecays (notUsed: allAreClosed)
107 static G4bool ElMaDecays; // Flag for opening ElectroMagDecays (true by default)
108 static G4bool EnergyFlux; // Flag for Energy Flux use instead of Multy Quasmon
109 static G4double SolidAngle; // Part of Solid Angle to capture secondaries(@@A-dep)
110 static G4double PiPrThresh; // Pion Production Threshold for gammas
111 static G4double M2ShiftVir; // Shift for M2=-Q2=m_pi^2 of the virtual gamma
112 static G4double DiNuclMass; // Double Nucleon Mass for virtual normalization
113 // Output hadrons
114 G4QHadronVector theQHadrons; // Vector of generated secondary hadrons
115 // Internal working objects
116 G4QHadronVector intQHadrons; // Vector of postponed secondary hadrons
117 G4QCHIPSWorld* theWorld; // the CHIPS World
118 G4int nBarClust; // Maximum BarionNumber of clusters (To optimize calc)
119 G4double f2all; // Ratio of freeNucleons to free+denseNucleons
120 G4QuasmonVector theQuasmons; // Intermediate vectorOfQuasmons before fragmentation
121 G4QCandidateVector theQCandidates; // Vector of possible candidates to clusters
122 G4QNucleus theEnvironment; // InitialNucleus (later ResidualNuclearEnvironment)
123 G4LorentzVector tot4Mom; // Total 4-momentum in the reaction
124 G4int totCharge; // Total charge in the reaction (for current control)
125 G4int totBaryoN; // Total baryon number in the reaction (for cur.cont)
126 G4QHadronVector theProjectiles; // Vector of projectiles in the interaction
127 G4int theTargetPDG; // PDG of the target nucleus in the interaction
128 G4QFreeScattering* theQFScat; // Pointer to the CHIPS Quasi-Free Scatterer
129};
130
131// Inline functions
133 {return this == &rhs;}
135 {return this != &rhs;}
136inline G4QNucleus G4QEnvironment::GetEnvironment() const {return theEnvironment;}
137
138#endif
std::vector< G4QCandidate * > G4QCandidateVector
std::vector< G4QHadron * > G4QHadronVector
std::vector< G4Quasmon * > G4QuasmonVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
void DecayMeson(G4QHadron *dB, G4QHadronVector *HV)
G4QHadronVector * GetProjectiles()
void FillQHadrons(G4QHadronVector *input)
void DecayAntistrange(G4QHadron *aS, G4QHadronVector *HV)
G4bool operator!=(const G4QEnvironment &right) const
static void CloseElectromagneticDecays()
G4QHadronVector * GetQHadrons()
G4QuasmonVector * GetQuasmons()
void CheckMassShell(G4QHadronVector *HV)
G4QNucleus GetEnvironment() const
void AddQuasmon(G4Quasmon *Q)
G4bool operator==(const G4QEnvironment &right) const
static void SetParameters(G4double solAn=0.4, G4bool efFlag=false, G4double piThresh=141.4, G4double mpisq=20000., G4double dinum=1880.)
void DecayBaryon(G4QHadron *dB, G4QHadronVector *HV)
const G4QEnvironment & operator=(const G4QEnvironment &right)
G4QHadronVector * Fragment()
static void OpenElectromagneticDecays()