Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Member Functions
exrdmAnalysisManager Class Reference

#include <exrdmAnalysisManager.hh>

Public Member Functions

void BookHisto ()
 
void BeginOfRun ()
 
void EndOfRun (G4int)
 
void BeginOfEvent ()
 
void EndOfEvent ()
 
void AddParticle (G4double, G4double, G4double, G4double)
 
void AddIsotope (G4double, G4double, G4double)
 
void AddEnergy (G4double, G4double, G4double)
 
void AddDecayProduct (G4double pid, G4int Z, G4int A, G4double energy, G4double time, G4double weight)
 
void SetVerbose (G4int val)
 
G4int GetVerbose () const
 
void SetFirstEventToDebug (G4int val)
 
G4int FirstEventToDebug () const
 
void SetLastEventToDebug (G4int val)
 
G4int LastEventToDebug () const
 
void SetMaxEnergyforHisto (G4double val)
 
G4double GetMaxEnergyforHisto () const
 
void SetMinEnergyforHisto (G4double val)
 
G4double GetMinEnergyforHisto () const
 
void SetNumBinforHisto (G4int val)
 
G4int GeNumBinforHisto () const
 
void SetThresholdEnergyforTarget (G4double val)
 
G4double GetThresholdEnergyforTarget () const
 
void SetThresholdEnergyforDetector (G4double val)
 
G4double GetThresholdEnergyforDetector () const
 
void SetPulseWidth (G4double val)
 
G4double GetPulseWidth () const
 

Static Public Member Functions

static exrdmAnalysisManagerGetInstance ()
 
static void Dispose ()
 

Detailed Description

Definition at line 61 of file exrdmAnalysisManager.hh.

Member Function Documentation

void exrdmAnalysisManager::AddDecayProduct ( G4double  pid,
G4int  Z,
G4int  A,
G4double  energy,
G4double  time,
G4double  weight 
)

Definition at line 329 of file exrdmAnalysisManager.cc.

References exrdmHisto::AddRow(), and exrdmHisto::FillTuple().

Referenced by exrdmSteppingAction::UserSteppingAction().

332 {
333  fHisto->FillTuple(3,0,pid);
334  fHisto->FillTuple(3,1,double(Z));
335  fHisto->FillTuple(3,2,double(A));
336  fHisto->FillTuple(3,3,energy);
337  fHisto->FillTuple(3,4,time/s);
338  fHisto->FillTuple(3,5,weight);
339  fHisto->AddRow(3);
340 }
void FillTuple(G4int, const G4String &, G4double)
Definition: exrdmHisto.cc:375
const XML_Char * s
void AddRow(G4int)
Definition: exrdmHisto.cc:419
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
void exrdmAnalysisManager::AddEnergy ( G4double  edep,
G4double  weight,
G4double  time 
)

Definition at line 279 of file exrdmAnalysisManager.cc.

References exrdmHisto::AddRow(), exrdmHisto::FillTuple(), G4cout, G4endl, python.hepunit::keV, python.hepunit::MeV, and python.hepunit::second.

Referenced by EndOfEvent(), and exrdmSteppingAction::UserSteppingAction().

281 {
282  if(1 < fVerbose) {
283  G4cout << "exrdmAnalysisManager::AddEnergy: e(keV)= " << edep/keV
284  << " weight = " << weight << " time (s) = " << time/second
285  << G4endl;
286  }
287  fHisto->FillTuple(2, 0, edep/MeV);
288  fHisto->FillTuple(2,1,time/second);
289  fHisto->FillTuple(2,2,weight);
290  fHisto->AddRow(2);
291  //
292  exrdmEnergyDeposition A(edep,time,weight);
293  fEdepo.push_back(A);
294 }
void FillTuple(G4int, const G4String &, G4double)
Definition: exrdmHisto.cc:375
void AddRow(G4int)
Definition: exrdmHisto.cc:419
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void exrdmAnalysisManager::AddIsotope ( G4double  pid,
G4double  weight,
G4double  time 
)

Definition at line 315 of file exrdmAnalysisManager.cc.

References exrdmHisto::AddRow(), exrdmHisto::FillTuple(), G4cout, G4endl, and python.hepunit::second.

Referenced by exrdmSteppingAction::UserSteppingAction().

317 {
318  if(1 < fVerbose) {
319  G4cout << "exrdmAnalysisManager::AddIsotope: " << pid
320  << G4endl;
321  }
322  fHisto->FillTuple(1,0,pid);
323  fHisto->FillTuple(1,1,time/second);
324  fHisto->FillTuple(1,2,weight);
325  fHisto->AddRow(1);
326 }
void FillTuple(G4int, const G4String &, G4double)
Definition: exrdmHisto.cc:375
void AddRow(G4int)
Definition: exrdmHisto.cc:419
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void exrdmAnalysisManager::AddParticle ( G4double  pid,
G4double  energy,
G4double  weight,
G4double  time 
)

Definition at line 298 of file exrdmAnalysisManager.cc.

References exrdmHisto::AddRow(), exrdmHisto::FillHisto(), exrdmHisto::FillTuple(), G4cout, G4endl, python.hepunit::MeV, and python.hepunit::second.

Referenced by exrdmSteppingAction::UserSteppingAction().

300 {
301  if(1 < fVerbose) {
302  G4cout << "exrdmAnalysisManager::AddParticle: " << pid
303  << G4endl;
304  }
305  fHisto->FillTuple(0,0, pid);
306  fHisto->FillTuple(0,1,energy/MeV);
307  fHisto->FillTuple(0,2,time/second);
308  fHisto->FillTuple(0,3,weight);
309  fHisto->AddRow(0);
310  // now fill th emission spectrum
311  if (energy>0.0) fHisto->FillHisto(6,energy/MeV,weight);
312 }
void FillTuple(G4int, const G4String &, G4double)
Definition: exrdmHisto.cc:375
void FillHisto(G4int, G4double, G4double)
Definition: exrdmHisto.cc:274
void AddRow(G4int)
Definition: exrdmHisto.cc:419
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void exrdmAnalysisManager::BeginOfEvent ( )

Definition at line 196 of file exrdmAnalysisManager.cc.

Referenced by exrdmEventAction::BeginOfEventAction().

197 {
198  fEdepo.clear();
199 }
void exrdmAnalysisManager::BeginOfRun ( )

Definition at line 130 of file exrdmAnalysisManager.cc.

References exrdmHisto::Book(), G4ProcessTable::FindProcess(), G4cout, G4endl, G4ProcessTable::GetProcessTable(), G4RadioactiveDecay::GetTheRadioactivityTables(), and G4RadioactiveDecay::IsAnalogueMonteCarlo().

Referenced by exrdmRunAction::BeginOfRunAction().

131 {
132  fHisto->Book();
133  G4cout
134  << "exrdmAnalysisManager: Histograms are booked and the run has been started"
135  << G4endl;
138  pTable->FindProcess("RadioactiveDecay", "GenericIon");
139  if (rDecay != NULL) {
140  if (!rDecay->IsAnalogueMonteCarlo()) {
141  std::vector<G4RadioactivityTable*> theTables =
142  rDecay->GetTheRadioactivityTables();
143  for (size_t i = 0 ; i < theTables.size(); i++) {
144  theTables[i]->GetTheMap()->clear();
145  G4cout << " Clear the radioactivity map: 0, new size = "
146  << theTables[i]->GetTheMap()->size() << G4endl;
147  }
148  }
149  }
150 
151 }
G4GLOB_DLL std::ostream G4cout
std::vector< G4RadioactivityTable * > GetTheRadioactivityTables()
void Book()
Definition: exrdmHisto.cc:113
#define G4endl
Definition: G4ios.hh:61
static G4ProcessTable * GetProcessTable()
G4VProcess * FindProcess(const G4String &processName, const G4String &particleName) const
void exrdmAnalysisManager::BookHisto ( )

Definition at line 94 of file exrdmAnalysisManager.cc.

References exrdmHisto::Add1D(), exrdmHisto::AddTuple(), and python.hepunit::MeV.

95 {
96  fHistEMax = 15.0*MeV;
97  fHistEMin = .0*MeV;
98  fHistNBin = 100;
99 
100  fHisto->Add1D("H10",
101  "Energy deposit (MeV) in the traget",fHistNBin,fHistEMin,fHistEMax,MeV);
102  fHisto->Add1D("H11",
103  "Energy deposit (MeV) in the detector",fHistNBin,fHistEMin,fHistEMax,MeV);
104  fHisto->Add1D("H12",
105  "Total energy spectrum (MeV) of the traget and detector",fHistNBin,
106  fHistEMin,fHistEMax,MeV);
107  fHisto->Add1D("H13",
108  "Coincidence spectrum (MeV) between the traget and detector",fHistNBin,
109  fHistEMin,fHistEMax,MeV);
110  fHisto->Add1D("H14",
111  "Anti-coincidence spectrum (MeV) in the traget",fHistNBin,
112  fHistEMin,fHistEMax,MeV);
113  fHisto->Add1D("H15",
114  "Anti-coincidence spectrum (MeV) in the detector",fHistNBin,
115  fHistEMin,fHistEMax,MeV);
116  fHisto->Add1D("H16",
117  "Decay emission spectrum (MeV)",fHistNBin,fHistEMin,fHistEMax,MeV);
118  // in aida these histos are indiced from 0-6
119  //
120  fHisto->AddTuple( "T1", "Emitted Particles",
121  "double PID, Energy, Time, Weight" );
122  fHisto->AddTuple( "T2", "RadioIsotopes","double PID, Time, Weight" );
123  fHisto->AddTuple( "T3", "Energy Depositions","double Energy, Time, Weight" );
124  fHisto->AddTuple( "RDecayProducts", "All Products of RDecay",
125  "double PID, Z, A, Energy, Time, Weight" );
126 }
void AddTuple(const G4String &, const G4String &, const G4String &)
Definition: exrdmHisto.cc:333
void Add1D(const G4String &, const G4String &, G4int nb=100, G4double x1=0., G4double x2=1., G4double u=1.)
Definition: exrdmHisto.cc:226
void exrdmAnalysisManager::Dispose ( )
static

Definition at line 59 of file exrdmAnalysisManager.cc.

Referenced by main().

60 {
61  delete fManager;
62  fManager = 0;
63 }
void exrdmAnalysisManager::EndOfEvent ( )

Definition at line 203 of file exrdmAnalysisManager.cc.

References AddEnergy(), exrdmHisto::FillHisto(), and sort().

Referenced by exrdmEventAction::EndOfEventAction().

204 {
205  if (fEdepo.size()) {
206  std::sort(fEdepo.begin(),fEdepo.end());
207  G4double TarE = fEdepo[0].GetEnergy();
208  G4double Time = fEdepo[0].GetTime();
209  G4double TarW = fEdepo[0].GetEnergy()*fEdepo[0].GetWeight();
210  G4double DetE = 0.;
211  G4double DetW = 0.;
212  G4double ComW = 0.;
213  if (TarE< 0.) {
214  DetE = -TarE;
215  DetW = -TarW;
216  TarE = 0.;
217  TarW = 0.;
218  }
219  for (size_t i = 1; i < fEdepo.size(); i++) {
220  if (std::fabs(fEdepo[i].GetTime() - Time) <= fPulseWidth) {
221  if ( fEdepo[i].GetEnergy() > 0. ) {
222  TarE += fEdepo[i].GetEnergy();
223  TarW += fEdepo[i].GetEnergy()*fEdepo[i].GetWeight();
224  } else {
225  DetE -= fEdepo[i].GetEnergy();
226  DetW -= fEdepo[i].GetEnergy()*fEdepo[i].GetWeight();
227  }
228  } else {
229  // tallying
230  if (TarE || DetE) ComW = (TarW+DetW)/(TarE+DetE);
231  if (TarE) TarW /= TarE;
232  if (DetE) DetW /= DetE;
233  //
234  if (TarE) fHisto->FillHisto(0,TarE,TarW); // target histogram
235  if (DetE) fHisto->FillHisto(1,DetE,DetW); // Detector histogram
236  if (TarE+DetE) fHisto->FillHisto(2,DetE+TarE,ComW); // Total histogram
237  if (DetE >= fDetectorThresE && TarE >= fTargetThresE )
238  fHisto->FillHisto(3,DetE,DetW); // coincidence histogram
239  if (TarE >= fTargetThresE && DetE < fDetectorThresE)
240  fHisto->FillHisto(4,TarE,TarW); // target anti-coincidence histogram
241  if (TarE < fTargetThresE && DetE >= fDetectorThresE)
242  fHisto->FillHisto(5,DetE,DetW); // detector anti-coincidence histogram
243  //
244  //reset
245  TarE = fEdepo[i].GetEnergy();
246  Time = fEdepo[i].GetTime();
247  TarW = fEdepo[i].GetEnergy()*fEdepo[i].GetWeight();
248  DetE = 0.;
249  DetW = 0.;
250  if (TarE < 0) {
251  DetE = -TarE;
252  DetW = -TarW;
253  TarE = 0.;
254  TarW = 0.;
255  }
256  }
257  }
258  //tally the last hit
259  if (TarE || DetE) ComW = (TarW+DetW)/(TarE+DetE);
260  if (TarE) TarW /= TarE;
261  if (DetE) DetW /= DetE;
262  //
263  if (TarE) fHisto->FillHisto(0,TarE,TarW); // target histogram
264  if (DetE) fHisto->FillHisto(1,DetE,DetW); // Detector histogram
265  if (TarE+DetE) fHisto->FillHisto(2,DetE+TarE,ComW); // Total histogram
266  if (DetE >= fDetectorThresE && TarE >= fTargetThresE )
267  fHisto->FillHisto(3,DetE,DetW); // coincidence histogram
268  if (TarE >= fTargetThresE && DetE < fDetectorThresE)
269  fHisto->FillHisto(4,TarE,TarW); // target anti-coincidence histogram
270  if (TarE < fTargetThresE && DetE >= fDetectorThresE)
271  fHisto->FillHisto(5,DetE,DetW); // detector anti-coincidence histogram
272  // now add zero energy to separat events
273  AddEnergy(0.,0.,0.);
274  }
275 }
void AddEnergy(G4double, G4double, G4double)
void FillHisto(G4int, G4double, G4double)
Definition: exrdmHisto.cc:274
subroutine sort(A, N)
Definition: dpm25nuc7.f:4670
double G4double
Definition: G4Types.hh:76
void exrdmAnalysisManager::EndOfRun ( G4int  nevent)

Definition at line 155 of file exrdmAnalysisManager.cc.

References G4ProcessTable::FindProcess(), G4cout, G4endl, exrdmHisto::GetFileName(), G4ProcessTable::GetProcessTable(), G4RadioactiveDecay::GetTheRadioactivityTables(), G4RadioactiveDecay::IsAnalogueMonteCarlo(), and exrdmHisto::Save().

Referenced by exrdmRunAction::EndOfRunAction().

156 {
157  fHisto->Save();
158  G4cout << "exrdmAnalysisManager: Histograms have been saved!" << G4endl;
159 
160  // Output the induced radioactivities
161  // in VR mode only
162  //
165  pTable->FindProcess("RadioactiveDecay", "GenericIon");
166  if (rDecay != NULL) {
167  if (!rDecay->IsAnalogueMonteCarlo()) {
168  G4String fileName = fHisto->GetFileName() + ".activity" ;
169  std::ofstream outfile (fileName, std::ios::out );
170  std::vector<G4RadioactivityTable*> theTables =
171  rDecay->GetTheRadioactivityTables();
172  for (size_t i = 0 ; i < theTables.size(); i++) {
173  G4double rate, error;
174  outfile << "Radioactivities in decay window no. " << i << G4endl;
175  outfile
176  << "Z \tA \tE \tActivity (decays/window) \tError (decays/window) "
177  << G4endl;
178  map<G4ThreeVector,G4TwoVector> *aMap = theTables[i]->GetTheMap();
179  map<G4ThreeVector,G4TwoVector>::iterator iter;
180  for(iter=aMap->begin(); iter != aMap->end(); iter++) {
181  rate = iter->second.x()/nevent;
182  error = std::sqrt(iter->second.y())/nevent;
183  if ( rate < 0.) rate = 0.; // statically it can be < 0.
184  outfile << iter->first.x() <<"\t"<< iter->first.y() <<"\t"
185  << iter->first.z() << "\t" << rate <<"\t" << error
186  << G4endl;
187  }
188  outfile << G4endl;
189  }
190  outfile.close();
191  }
192  }
193 }
const G4String & GetFileName() const
Definition: exrdmHisto.cc:447
G4GLOB_DLL std::ostream G4cout
std::vector< G4RadioactivityTable * > GetTheRadioactivityTables()
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static G4ProcessTable * GetProcessTable()
G4VProcess * FindProcess(const G4String &processName, const G4String &particleName) const
void Save()
Definition: exrdmHisto.cc:196
G4int exrdmAnalysisManager::FirstEventToDebug ( ) const
inline

Definition at line 95 of file exrdmAnalysisManager.hh.

95 {return fNEvt1;};
G4int exrdmAnalysisManager::GeNumBinforHisto ( ) const
inline

Definition at line 104 of file exrdmAnalysisManager.hh.

104 {return fHistNBin;};
exrdmAnalysisManager * exrdmAnalysisManager::GetInstance ( void  )
static
G4double exrdmAnalysisManager::GetMaxEnergyforHisto ( ) const
inline

Definition at line 100 of file exrdmAnalysisManager.hh.

100 {return fHistEMax;};
G4double exrdmAnalysisManager::GetMinEnergyforHisto ( ) const
inline

Definition at line 102 of file exrdmAnalysisManager.hh.

102 {return fHistEMin;};
G4double exrdmAnalysisManager::GetPulseWidth ( ) const
inline

Definition at line 111 of file exrdmAnalysisManager.hh.

111 {return fPulseWidth;};
G4double exrdmAnalysisManager::GetThresholdEnergyforDetector ( ) const
inline

Definition at line 109 of file exrdmAnalysisManager.hh.

109 {return fDetectorThresE;};
G4double exrdmAnalysisManager::GetThresholdEnergyforTarget ( ) const
inline

Definition at line 107 of file exrdmAnalysisManager.hh.

107 {return fTargetThresE;};
G4int exrdmAnalysisManager::GetVerbose ( ) const
inline

Definition at line 92 of file exrdmAnalysisManager.hh.

92 {return fVerbose;};
G4int exrdmAnalysisManager::LastEventToDebug ( ) const
inline

Definition at line 97 of file exrdmAnalysisManager.hh.

97 {return fNEvt2;};
void exrdmAnalysisManager::SetFirstEventToDebug ( G4int  val)
inline

Definition at line 94 of file exrdmAnalysisManager.hh.

94 {fNEvt1 = val;};
void exrdmAnalysisManager::SetLastEventToDebug ( G4int  val)
inline

Definition at line 96 of file exrdmAnalysisManager.hh.

96 {fNEvt2 = val;};
void exrdmAnalysisManager::SetMaxEnergyforHisto ( G4double  val)
inline

Definition at line 99 of file exrdmAnalysisManager.hh.

99 {fHistEMax = val;};
void exrdmAnalysisManager::SetMinEnergyforHisto ( G4double  val)
inline

Definition at line 101 of file exrdmAnalysisManager.hh.

101 {fHistEMin = val;};
void exrdmAnalysisManager::SetNumBinforHisto ( G4int  val)
inline

Definition at line 103 of file exrdmAnalysisManager.hh.

103 {fHistNBin = val;};
void exrdmAnalysisManager::SetPulseWidth ( G4double  val)
inline

Definition at line 110 of file exrdmAnalysisManager.hh.

110 {fPulseWidth = val;};
void exrdmAnalysisManager::SetThresholdEnergyforDetector ( G4double  val)
inline

Definition at line 108 of file exrdmAnalysisManager.hh.

108 {fDetectorThresE = val;};
void exrdmAnalysisManager::SetThresholdEnergyforTarget ( G4double  val)
inline

Definition at line 106 of file exrdmAnalysisManager.hh.

106 {fTargetThresE = val;};
void exrdmAnalysisManager::SetVerbose ( G4int  val)
inline

Definition at line 91 of file exrdmAnalysisManager.hh.

91 {fVerbose = val;};

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