QwGeant4
HistoManager Class Reference

#include <HistoManager.hh>

+ Collaboration diagram for HistoManager:

Public Member Functions

 ~HistoManager ()
 
void BeginOfRun ()
 
void EndOfRun ()
 
void SetVerbose (G4int val)
 
void SetParticleName (const G4String &)
 
void SetElementName (const G4String &)
 
void SetNumberOfBinsE (G4int val)
 
void SetNumberOfBinsP (G4int val)
 
void SetMinKinEnergy (G4double val)
 
void SetMaxKinEnergy (G4double val)
 
void SetMinMomentum (G4double val)
 
void SetMaxMomentum (G4double val)
 
G4int GetVerbose () const
 

Static Public Member Functions

static HistoManagerGetPointer ()
 

Private Member Functions

 HistoManager ()
 

Private Attributes

HistofHisto
 
const G4ParticleDefinition * fNeutron
 
G4String fParticleName
 
G4String fElementName
 
G4double fMinKinEnergy
 
G4double fMaxKinEnergy
 
G4double fMinMomentum
 
G4double fMaxMomentum
 
G4int fVerbose
 
G4int fBinsE
 
G4int fBinsP
 
G4bool fIsInitialised
 

Static Private Attributes

static HistoManagerfManager = 0
 

Detailed Description

Definition at line 60 of file HistoManager.hh.

Constructor & Destructor Documentation

HistoManager::HistoManager ( )
private

Definition at line 81 of file HistoManager.cc.

References fBinsE, fBinsP, fElementName, fHisto, fIsInitialised, fMaxKinEnergy, fMaxMomentum, fMinKinEnergy, fMinMomentum, fNeutron, fParticleName, fVerbose, and Histo::SetVerbose().

82 {
83  fHisto = new Histo();
84  fNeutron = G4Neutron::Neutron();
85  fVerbose = 1;
86 
87  fParticleName = "proton";
88  fElementName = "Al";
89 
90  fMinKinEnergy = 0.1*MeV;
91  fMaxKinEnergy = 10*TeV;
92  fMinMomentum = 1*MeV;
93  fMaxMomentum = 10*TeV;
94 
95  fBinsE = 800;
96  fBinsP = 700;
97 
99 
100  fIsInitialised = false;
101 }
G4bool fIsInitialised
G4String fElementName
Definition: Histo.hh:56
G4double fMinKinEnergy
const G4ParticleDefinition * fNeutron
Definition: HistoManager.hh:99
G4double fMinMomentum
void SetVerbose(G4int val)
Definition: Histo.hh:106
G4double fMaxMomentum
G4double fMaxKinEnergy
G4String fParticleName
Histo * fHisto
Definition: HistoManager.hh:97

+ Here is the call graph for this function:

HistoManager::~HistoManager ( )

Definition at line 105 of file HistoManager.cc.

References fHisto.

106 {
107  delete fHisto;
108 }
Histo * fHisto
Definition: HistoManager.hh:97

Member Function Documentation

void HistoManager::BeginOfRun ( )

Definition at line 112 of file HistoManager.cc.

References Histo::Add1D(), Histo::Book(), fBinsE, fBinsP, fHisto, fIsInitialised, fMaxKinEnergy, fMaxMomentum, fMinKinEnergy, fMinMomentum, and Histo::SetHisto1D().

113 {
114  G4double p1 = std::log10(fMinMomentum/GeV);
115  G4double p2 = std::log10(fMaxMomentum/GeV);
116  G4double e1 = std::log10(fMinKinEnergy/MeV);
117  G4double e2 = std::log10(fMaxKinEnergy/MeV);
118 
119  //G4cout<<"e1= "<<e1<<" e2= "<<e2<<" p1= "<<p1<<" p2= "<<p2<<G4endl;
120 
121  if(fIsInitialised) {
122  fHisto->SetHisto1D(0,fBinsP,p1,p2,1.0);
123  fHisto->SetHisto1D(1,fBinsE,e1,e2,1.0);
124  fHisto->SetHisto1D(2,fBinsP,p1,p2,1.0);
125  fHisto->SetHisto1D(3,fBinsE,e1,e2,1.0);
126  fHisto->SetHisto1D(4,fBinsE,e1,e2,1.0);
127  fHisto->SetHisto1D(5,fBinsE,e1,e2,1.0);
128  fHisto->SetHisto1D(6,fBinsE,e1,e2,1.0);
129  fHisto->SetHisto1D(7,fBinsE,e1,e2,1.0);
130 
131  } else {
132  fHisto->Add1D("h1","Elastic cross section (barn,1.0) as a functions of log10(p/GeV)",
133  fBinsP,p1,p2,1.0);
134  fHisto->Add1D("h2","Elastic cross section (barn) as a functions of log10(E/MeV)",
135  fBinsE,e1,e2,1.0);
136  fHisto->Add1D("h3","Inelastic cross section (barn) as a functions of log10(p/GeV)",
137  fBinsP,p1,p2,1.0);
138  fHisto->Add1D("h4","Inelastic cross section (barn) as a functions of log10(E/MeV)",
139  fBinsE,e1,e2,1.0);
140  fHisto->Add1D("h5","Capture cross section (barn) as a functions of log10(E/MeV)",
141  fBinsE,e1,e2,1.0);
142  fHisto->Add1D("h6","Fission cross section (barn) as a functions of log10(E/MeV)",
143  fBinsE,e1,e2,1.0);
144  fHisto->Add1D("h7","Charge exchange cross section (barn) as a functions of log10(E/MeV)",
145  fBinsE,e1,e2,1.0);
146  fHisto->Add1D("h8","Total cross section (barn) as a functions of log10(E/MeV)",
147  fBinsE,e1,e2,1.0);
148  }
149 
150  fIsInitialised = true;
151  fHisto->Book();
152 }
G4bool fIsInitialised
G4double fMinKinEnergy
G4double fMinMomentum
G4double fMaxMomentum
void Book()
Definition: Histo.cc:77
G4double fMaxKinEnergy
void Add1D(const G4String &, const G4String &, G4int nb, G4double x1, G4double x2, G4double u=1.)
Definition: Histo.cc:153
void SetHisto1D(G4int, G4int, G4double, G4double, G4double)
Definition: Histo.cc:176
Histo * fHisto
Definition: HistoManager.hh:97

+ Here is the call graph for this function:

void HistoManager::EndOfRun ( )

Definition at line 156 of file HistoManager.cc.

References fBinsE, fBinsP, fElementName, fHisto, Histo::Fill(), fMaxKinEnergy, fMaxMomentum, fMinKinEnergy, fMinMomentum, fNeutron, fParticleName, fVerbose, and Histo::Save().

157 {
158  if(fVerbose > 0) {
159  G4cout << "HistoManager: End of run actions are started" << G4endl;
160  }
161 
162  const G4Element* elm =
163  G4NistManager::Instance()->FindOrBuildElement(fElementName);
164  const G4Material* mat =
165  G4NistManager::Instance()->FindOrBuildMaterial("G4_"+fElementName);
166  const G4ParticleDefinition* particle =
167  G4ParticleTable::GetParticleTable()->FindParticle(fParticleName);
168 
169  G4cout << "### Fill Cross Sections for " << fParticleName
170  << " off " << fElementName
171  << G4endl;
172  if(fVerbose > 0) {
173  G4cout << "------------------------------------------------------------------------"
174  << G4endl;
175  G4cout << " N E(MeV) Elastic(b) Inelastic(b)";
176  if(particle == fNeutron) { G4cout << " Capture(b) Fission(b)"; }
177  G4cout << " Total(b)" << G4endl;
178  G4cout << "------------------------------------------------------------------------"
179  << G4endl;
180  }
181  if(!particle || !elm) {
182  G4cout << "HistoManager WARNING Particle or element undefined" << G4endl;
183  return;
184  }
185 
186  G4int prec = G4cout.precision();
187  G4cout.precision(4);
188 
189  G4HadronicProcessStore* store = G4HadronicProcessStore::Instance();
190  G4double mass = particle->GetPDGMass();
191 
192  // Build histograms
193 
194  G4double p1 = std::log10(fMinMomentum/GeV);
195  G4double p2 = std::log10(fMaxMomentum/GeV);
196  G4double e1 = std::log10(fMinKinEnergy/MeV);
197  G4double e2 = std::log10(fMaxKinEnergy/MeV);
198  G4double de = (e2 - e1)/G4double(fBinsE);
199  G4double dp = (p2 - p1)/G4double(fBinsP);
200 
201  G4double x = e1 - de*0.5;
202  G4double e, p, xs, xtot;
203  G4int i;
204  for(i=0; i<fBinsE; i++) {
205  x += de;
206  e = std::pow(10.,x)*MeV;
207  if(fVerbose>0) G4cout << std::setw(5) << i << std::setw(12) << e;
208  xs = store->GetElasticCrossSectionPerAtom(particle,e,elm,mat);
209  xtot = xs;
210  if(fVerbose>0) G4cout << std::setw(12) << xs/barn;
211  fHisto->Fill(1, x, xs/barn);
212  xs = store->GetInelasticCrossSectionPerAtom(particle,e,elm,mat);
213  xtot += xs;
214  if(fVerbose>0) G4cout << " " << std::setw(12) << xs/barn;
215  fHisto->Fill(3, x, xs/barn);
216  if(particle == fNeutron) {
217  xs = store->GetCaptureCrossSectionPerAtom(particle,e,elm,mat);
218  xtot += xs;
219  if(fVerbose>0) G4cout << " " << std::setw(12) << xs/barn;
220  fHisto->Fill(4, x, xs/barn);
221  xs = store->GetFissionCrossSectionPerAtom(particle,e,elm,mat);
222  xtot += xs;
223  if(fVerbose>0) G4cout << " " << std::setw(12) << xs/barn;
224  fHisto->Fill(5, x, xs/barn);
225  }
226  xs = store->GetChargeExchangeCrossSectionPerAtom(particle,e,elm,mat);
227  if(fVerbose>0) G4cout << " " << std::setw(12) << xtot/barn << G4endl;
228  fHisto->Fill(6, x, xs/barn);
229  fHisto->Fill(7, x, xtot/barn);
230  }
231 
232  x = p1 - dp*0.5;
233  for(i=0; i<fBinsP; i++) {
234  x += dp;
235  p = std::pow(10.,x)*GeV;
236  e = std::sqrt(p*p + mass*mass) - mass;
237  xs = store->GetElasticCrossSectionPerAtom(particle,e,elm,mat);
238  fHisto->Fill(0, x, xs/barn);
239  xs = store->GetInelasticCrossSectionPerAtom(particle,e,elm,mat);
240  fHisto->Fill(2, x, xs/barn);
241  }
242  if(fVerbose > 0) {
243  G4cout << "-----------------------------------------------------------------"
244  << G4endl;
245  }
246  G4cout.precision(prec);
247  fHisto->Save();
248 }
G4String fElementName
G4double fMinKinEnergy
void Fill(G4int, G4double, G4double)
Definition: Histo.cc:211
void Save()
Definition: Histo.cc:129
const G4ParticleDefinition * fNeutron
Definition: HistoManager.hh:99
G4double fMinMomentum
G4double fMaxMomentum
G4double fMaxKinEnergy
G4String fParticleName
Histo * fHisto
Definition: HistoManager.hh:97

+ Here is the call graph for this function:

HistoManager * HistoManager::GetPointer ( )
static

Definition at line 70 of file HistoManager.cc.

References fManager.

Referenced by RunAction::BeginOfRunAction(), RunAction::EndOfRunAction(), and DetectorMessenger::SetNewValue().

71 {
72  if(!fManager) {
73  static HistoManager manager;
74  fManager = &manager;
75  }
76  return fManager;
77 }
static HistoManager * fManager
Definition: HistoManager.hh:95

+ Here is the caller graph for this function:

G4int HistoManager::GetVerbose ( ) const
inline

Definition at line 156 of file HistoManager.hh.

References fVerbose.

157 {
158  return fVerbose;
159 }
void HistoManager::SetElementName ( const G4String &  name)
inline

Definition at line 121 of file HistoManager.hh.

References fElementName.

Referenced by DetectorMessenger::SetNewValue().

122 {
123  fElementName = name;
124 }
G4String fElementName

+ Here is the caller graph for this function:

void HistoManager::SetMaxKinEnergy ( G4double  val)
inline

Definition at line 141 of file HistoManager.hh.

References fMaxKinEnergy, and fMinKinEnergy.

Referenced by DetectorMessenger::SetNewValue().

142 {
143  if(val>fMinKinEnergy) { fMaxKinEnergy = val; }
144 }
G4double fMinKinEnergy
G4double fMaxKinEnergy

+ Here is the caller graph for this function:

void HistoManager::SetMaxMomentum ( G4double  val)
inline

Definition at line 151 of file HistoManager.hh.

References fMaxMomentum, and fMinMomentum.

Referenced by DetectorMessenger::SetNewValue().

152 {
153  if(val>fMinMomentum) { fMaxMomentum = val; }
154 }
G4double fMinMomentum
G4double fMaxMomentum

+ Here is the caller graph for this function:

void HistoManager::SetMinKinEnergy ( G4double  val)
inline

Definition at line 136 of file HistoManager.hh.

References fMaxKinEnergy, and fMinKinEnergy.

Referenced by DetectorMessenger::SetNewValue().

137 {
138  if(val>0 && val<fMaxKinEnergy) { fMinKinEnergy = val; }
139 }
G4double fMinKinEnergy
G4double fMaxKinEnergy

+ Here is the caller graph for this function:

void HistoManager::SetMinMomentum ( G4double  val)
inline

Definition at line 146 of file HistoManager.hh.

References fMaxMomentum, and fMinMomentum.

Referenced by DetectorMessenger::SetNewValue().

147 {
148  if(val>0 && val<fMaxMomentum) { fMinMomentum = val; }
149 }
G4double fMinMomentum
G4double fMaxMomentum

+ Here is the caller graph for this function:

void HistoManager::SetNumberOfBinsE ( G4int  val)
inline

Definition at line 126 of file HistoManager.hh.

References fBinsE.

Referenced by DetectorMessenger::SetNewValue().

127 {
128  if(val>0) { fBinsE = val; }
129 }

+ Here is the caller graph for this function:

void HistoManager::SetNumberOfBinsP ( G4int  val)
inline

Definition at line 131 of file HistoManager.hh.

References fBinsP.

Referenced by DetectorMessenger::SetNewValue().

132 {
133  if(val>0) { fBinsP = val; }
134 }

+ Here is the caller graph for this function:

void HistoManager::SetParticleName ( const G4String &  name)
inline

Definition at line 116 of file HistoManager.hh.

References fParticleName.

Referenced by DetectorMessenger::SetNewValue().

117 {
118  fParticleName = name;
119 }
G4String fParticleName

+ Here is the caller graph for this function:

void HistoManager::SetVerbose ( G4int  val)

Definition at line 252 of file HistoManager.cc.

References fHisto, fVerbose, and Histo::SetVerbose().

Referenced by DetectorMessenger::SetNewValue().

253 {
254  fVerbose = val;
255  fHisto->SetVerbose(val);
256 }
void SetVerbose(G4int val)
Definition: Histo.hh:106
Histo * fHisto
Definition: HistoManager.hh:97

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

Field Documentation

G4int HistoManager::fBinsE
private

Definition at line 110 of file HistoManager.hh.

Referenced by BeginOfRun(), EndOfRun(), HistoManager(), and SetNumberOfBinsE().

G4int HistoManager::fBinsP
private

Definition at line 111 of file HistoManager.hh.

Referenced by BeginOfRun(), EndOfRun(), HistoManager(), and SetNumberOfBinsP().

G4String HistoManager::fElementName
private

Definition at line 102 of file HistoManager.hh.

Referenced by EndOfRun(), HistoManager(), and SetElementName().

Histo* HistoManager::fHisto
private

Definition at line 97 of file HistoManager.hh.

Referenced by BeginOfRun(), EndOfRun(), HistoManager(), SetVerbose(), and ~HistoManager().

G4bool HistoManager::fIsInitialised
private

Definition at line 113 of file HistoManager.hh.

Referenced by BeginOfRun(), and HistoManager().

HistoManager * HistoManager::fManager = 0
staticprivate

Definition at line 95 of file HistoManager.hh.

Referenced by GetPointer().

G4double HistoManager::fMaxKinEnergy
private

Definition at line 105 of file HistoManager.hh.

Referenced by BeginOfRun(), EndOfRun(), HistoManager(), SetMaxKinEnergy(), and SetMinKinEnergy().

G4double HistoManager::fMaxMomentum
private

Definition at line 107 of file HistoManager.hh.

Referenced by BeginOfRun(), EndOfRun(), HistoManager(), SetMaxMomentum(), and SetMinMomentum().

G4double HistoManager::fMinKinEnergy
private

Definition at line 104 of file HistoManager.hh.

Referenced by BeginOfRun(), EndOfRun(), HistoManager(), SetMaxKinEnergy(), and SetMinKinEnergy().

G4double HistoManager::fMinMomentum
private

Definition at line 106 of file HistoManager.hh.

Referenced by BeginOfRun(), EndOfRun(), HistoManager(), SetMaxMomentum(), and SetMinMomentum().

const G4ParticleDefinition* HistoManager::fNeutron
private

Definition at line 99 of file HistoManager.hh.

Referenced by EndOfRun(), and HistoManager().

G4String HistoManager::fParticleName
private

Definition at line 101 of file HistoManager.hh.

Referenced by EndOfRun(), HistoManager(), and SetParticleName().

G4int HistoManager::fVerbose
private

Definition at line 109 of file HistoManager.hh.

Referenced by EndOfRun(), GetVerbose(), HistoManager(), and SetVerbose().


The documentation for this class was generated from the following files: