Geant4  10.01.p01
G4NeutronHPChannel.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 #ifndef G4NeutronHPChannel_h
38 #define G4NeutronHPChannel_h 1
39 #include "globals.hh"
40 #include "G4NeutronHPIsoData.hh"
41 #include "G4NeutronHPVector.hh"
42 #include "G4Material.hh"
43 #include "G4HadProjectile.hh"
44 //#include "G4NeutronInelasticProcess.hh"
45 //#include "G4HadronFissionProcess.hh"
46 //#include "G4HadronElasticProcess.hh"
47 //#include "G4HadronCaptureProcess.hh"
48 #include "G4StableIsotopes.hh"
49 #include "G4NeutronHPCaptureFS.hh"
50 #include "G4NeutronHPFinalState.hh"
51 #include "G4Element.hh"
53 
55 {
56 public:
57 
59  : wendtFissionGenerator(getenv("G4NEUTRON_HP_USE_WENDT_FISSION_MODEL") == NULL ? NULL : G4WendtFissionFragmentGenerator::GetInstance())
60  {
62  theBuffer = 0;
64  theFinalStates = 0;
65  active = 0;
66  registerCount = -1;
67  }
68 
70  {
71  delete theChannelData;
72  // Following statement disabled to avoid SEGV
73  // theBuffer is also deleted as "theChannelData" in
74  // ~G4NeutronHPIsoData. FWJ 06-Jul-1999
75  //if(theBuffer != 0) delete theBuffer;
76  if(theIsotopeWiseData != 0) delete [] theIsotopeWiseData;
77  // Deletion of FinalStates disabled to avoid endless looping
78  // in the destructor heirarchy. FWJ 06-Jul-1999
79  //if(theFinalStates != 0)
80  //{
81  // for(i=0; i<niso; i++)
82  // {
83  // delete theFinalStates[i];
84  // }
85  // delete [] theFinalStates;
86  //}
87  // FWJ experiment
88  //if(active!=0) delete [] active;
89 // T.K.
90  if ( theFinalStates != 0 )
91  {
92  for ( G4int i = 0 ; i < niso ; i++ )
93  {
94  delete theFinalStates[i];
95  }
96  delete [] theFinalStates;
97  }
98  if ( active != 0 ) delete [] active;
99 
100  }
101 
103 
105 
107 
108  inline G4bool IsActive(G4int isoNumber) { return active[isoNumber]; }
109 
110  inline G4bool HasFSData(G4int isoNumber) { return theFinalStates[isoNumber]->HasFSData(); }
111 
112  inline G4bool HasAnyData(G4int isoNumber) { return theFinalStates[isoNumber]->HasAnyData(); }
113 
115 
116  void Init(G4Element * theElement, const G4String dirName);
117 
118  void Init(G4Element * theElement, const G4String dirName, const G4String fsType);
119 
120  //void UpdateData(G4int A, G4int Z, G4int index, G4double abundance);
121  void UpdateData(G4int A, G4int Z, G4int index, G4double abundance) { G4int M = 0; UpdateData( A, Z, M, index, abundance); };
122  void UpdateData(G4int A, G4int Z, G4int M, G4int index, G4double abundance);
123 
124  void Harmonise(G4NeutronHPVector *& theStore, G4NeutronHPVector * theNew);
125 
126  G4HadFinalState * ApplyYourself(const G4HadProjectile & theTrack, G4int isoNumber=-1);
127 
128  inline G4int GetNiso() {return niso;}
129 
130  inline G4double GetN(G4int i) {return theFinalStates[i]->GetN();}
131  inline G4double GetZ(G4int i) {return theFinalStates[i]->GetZ();}
132  inline G4double GetM(G4int i) {return theFinalStates[i]->GetM();}
133 
135  {
136  G4bool result = false;
137  G4int i;
138  for(i=0; i<niso; i++)
139  {
140  if(theFinalStates[i]->HasAnyData()) result = true;
141  }
142  return result;
143  }
144 
145 private:
146  G4NeutronHPVector * theChannelData; // total (element) cross-section for this channel
148 
149  G4NeutronHPIsoData * theIsotopeWiseData; // these are the isotope-wise cross-sections for each final state.
150  G4NeutronHPFinalState ** theFinalStates; // also these are isotope-wise pionters, parallel to the above.
153 
155 
159 
161 
163 
164 };
165 
166 #endif
void Harmonise(G4NeutronHPVector *&theStore, G4NeutronHPVector *theNew)
G4bool Register(G4NeutronHPFinalState *theFS)
G4double GetWeightedXsec(G4double energy, G4int isoNumber)
G4NeutronHPFinalState ** theFinalStates
G4NeutronHPVector * theChannelData
G4double GetZ(G4int i)
void Init(G4Element *theElement, const G4String dirName)
int G4int
Definition: G4Types.hh:78
G4NeutronHPIsoData * theIsotopeWiseData
G4NeutronHPVector * theBuffer
G4double GetN(G4int i)
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)
bool G4bool
Definition: G4Types.hh:79
G4StableIsotopes theStableOnes
void UpdateData(G4int A, G4int Z, G4int index, G4double abundance)
static const G4double A[nN]
G4double GetXsec(G4double energy)
G4double GetFSCrossSection(G4double energy, G4int isoNumber)
G4double energy(const ThreeVector &p, const G4double m)
G4bool IsActive(G4int isoNumber)
G4bool HasAnyData(G4int isoNumber)
G4bool HasFSData(G4int isoNumber)
G4double GetM(G4int i)
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
double G4double
Definition: G4Types.hh:76
G4WendtFissionFragmentGenerator *const wendtFissionGenerator