Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ParticleHPChannel.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 //
28  // Hadronic Process: Very Low Energy Neutron X-Sections
29  // original by H.P. Wellisch, TRIUMF, 14-Feb-97
30  // Builds and has the Cross-section data for one element and channel.
31 //
32 // Bug fixes and workarounds in the destructor, F.W.Jones 06-Jul-1999
33 // 070612 Fix memory leaking by T. Koi
34 //
35 // 080520 Delete unnecessary dependencies by T. Koi
36 
37 // P. Arce, June-2014 Conversion neutron_hp to particle_hp
38 //
39 #ifndef G4ParticleHPChannel_h
40 #define G4ParticleHPChannel_h 1
41 #include "globals.hh"
42 #include "G4ParticleHPIsoData.hh"
43 #include "G4ParticleHPVector.hh"
44 #include "G4Material.hh"
45 #include "G4HadProjectile.hh"
46 //#include "G4NeutronInelasticProcess.hh"
47 //#include "G4HadronFissionProcess.hh"
48 //#include "G4HadronElasticProcess.hh"
49 //#include "G4HadronCaptureProcess.hh"
50 #include "G4StableIsotopes.hh"
51 #include "G4ParticleHPCaptureFS.hh"
53 #include "G4Element.hh"
56 
58 {
59 public:
60 
62  : wendtFissionGenerator(getenv("G4NEUTRON_HP_USE_WENDT_FISSION_MODEL") == NULL ? NULL : G4WendtFissionFragmentGenerator::GetInstance())
63  {
64  theProjectile = const_cast<G4ParticleDefinition*>(projectile);
65  theChannelData = new G4ParticleHPVector;
66  theBuffer = 0;
67  theIsotopeWiseData = 0;
68  theFinalStates = 0;
69  active = 0;
70  registerCount = -1;
71  niso = -1;
72  theElement = NULL;
73  }
74 
76  : wendtFissionGenerator(getenv("G4NEUTRON_HP_USE_WENDT_FISSION_MODEL") == NULL ? NULL : G4WendtFissionFragmentGenerator::GetInstance())
77  {
78  theProjectile = G4Neutron::Neutron();
79  theChannelData = new G4ParticleHPVector;
80  theBuffer = 0;
81  theIsotopeWiseData = 0;
82  theFinalStates = 0;
83  active = 0;
84  registerCount = -1;
85  niso = -1;
86  theElement = NULL;
87  }
88 
90  {
91  delete theChannelData;
92  // Following statement disabled to avoid SEGV
93  // theBuffer is also deleted as "theChannelData" in
94  // ~G4ParticleHPIsoData. FWJ 06-Jul-1999
95  //if(theBuffer != 0) delete theBuffer;
96  if(theIsotopeWiseData != 0) delete [] theIsotopeWiseData;
97  // Deletion of FinalStates disabled to avoid endless looping
98  // in the destructor heirarchy. FWJ 06-Jul-1999
99  //if(theFinalStates != 0)
100  //{
101  // for(i=0; i<niso; i++)
102  // {
103  // delete theFinalStates[i];
104  // }
105  // delete [] theFinalStates;
106  //}
107  // FWJ experiment
108  //if(active!=0) delete [] active;
109 // T.K.
110  if ( theFinalStates != 0 )
111  {
112  for ( G4int i = 0 ; i < niso ; i++ )
113  {
114  delete theFinalStates[i];
115  }
116  delete [] theFinalStates;
117  }
118  if ( active != 0 ) delete [] active;
119 
120  }
121 
123 
125 
127 
128  inline G4bool IsActive(G4int isoNumber) { return active[isoNumber]; }
129 
130  inline G4bool HasFSData(G4int isoNumber) { return theFinalStates[isoNumber]->HasFSData(); }
131 
132  inline G4bool HasAnyData(G4int isoNumber) { return theFinalStates[isoNumber]->HasAnyData(); }
133 
135 
136  void Init(G4Element * theElement, const G4String dirName);
137 
138  void Init(G4Element * theElement, const G4String dirName, const G4String fsType);
139 
140  //void UpdateData(G4int A, G4int Z, G4int index, G4double abundance);
141  void UpdateData(G4int A, G4int Z, G4int index, G4double abundance, G4ParticleDefinition* projectile) { G4int M = 0; UpdateData( A, Z, M, index, abundance, projectile); };
142  void UpdateData(G4int A, G4int Z, G4int M, G4int index, G4double abundance, G4ParticleDefinition* projectile);
143 
144  void Harmonise(G4ParticleHPVector *& theStore, G4ParticleHPVector * theNew);
145 
146  G4HadFinalState * ApplyYourself(const G4HadProjectile & theTrack, G4int isoNumber=-1);
147 
148  inline G4int GetNiso() {return niso;}
149 
150  inline G4double GetN(G4int i) {return theFinalStates[i]->GetN();}
151  inline G4double GetZ(G4int i) {return theFinalStates[i]->GetZ();}
152  inline G4double GetM(G4int i) {return theFinalStates[i]->GetM();}
153 
155  {
156  G4bool result = false;
157  G4int i;
158  for(i=0; i<niso; i++)
159  {
160  if(theFinalStates[i]->HasAnyData()) result = true;
161  }
162  return result;
163  }
164 
165  void DumpInfo();
166 
167  G4String GetFSType() const {
168  return theFSType;
169  }
170 
172  return theFinalStates;
173  }
174 
175 private:
176  G4ParticleDefinition* theProjectile;
177 
178  G4ParticleHPVector * theChannelData; // total (element) cross-section for this channel
179  G4ParticleHPVector * theBuffer;
180 
181  G4ParticleHPIsoData * theIsotopeWiseData; // these are the isotope-wise cross-sections for each final state.
182  G4ParticleHPFinalState ** theFinalStates; // also these are isotope-wise pionters, parallel to the above.
183  G4bool * active;
184  G4int niso;
185 
186  G4StableIsotopes theStableOnes;
187 
188  G4String theDir;
189  G4String theFSType;
190  G4Element * theElement;
191 
192  G4int registerCount;
193 
194  G4WendtFissionFragmentGenerator* const wendtFissionGenerator;
195 
196 };
197 
198 #endif
G4double G4ParticleHPJENDLHEData::G4double result
void UpdateData(G4int A, G4int Z, G4int index, G4double abundance, G4ParticleDefinition *projectile)
G4double GetM(G4int i)
G4double GetN(G4int i)
int G4int
Definition: G4Types.hh:78
G4bool Register(G4ParticleHPFinalState *theFS)
G4ParticleHPFinalState ** GetFinalStates() const
G4String GetFSType() const
void Harmonise(G4ParticleHPVector *&theStore, G4ParticleHPVector *theNew)
double A(double temperature)
bool G4bool
Definition: G4Types.hh:79
G4double GetZ(G4int i)
void Init(G4Element *theElement, const G4String dirName)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4bool IsActive(G4int isoNumber)
G4bool HasAnyData(G4int isoNumber)
G4ParticleHPChannel(G4ParticleDefinition *projectile)
G4double GetFSCrossSection(G4double energy, G4int isoNumber)
G4double energy(const ThreeVector &p, const G4double m)
G4double GetWeightedXsec(G4double energy, G4int isoNumber)
G4bool HasFSData(G4int isoNumber)
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
double G4double
Definition: G4Types.hh:76
G4double GetXsec(G4double energy)
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)