Geant4  10.02.p03
ElectronRun Class Reference

#include <ElectronRun.hh>

Inheritance diagram for ElectronRun:
Collaboration diagram for ElectronRun:

Public Member Functions

 ElectronRun ()
 
virtual ~ElectronRun ()
 
virtual void RecordEvent (const G4Event *)
 
virtual void Merge (const G4Run *)
 
void DumpData (G4String &) const
 
- Public Member Functions inherited from G4Run
 G4Run ()
 
virtual ~G4Run ()
 
G4int GetRunID () const
 
G4int GetNumberOfEvent () const
 
G4int GetNumberOfEventToBeProcessed () const
 
const G4HCtableGetHCtable () const
 
const G4DCtableGetDCtable () const
 
const G4StringGetRandomNumberStatus () const
 
void SetRunID (G4int id)
 
void SetNumberOfEventToBeProcessed (G4int n_ev)
 
void SetHCtable (G4HCtable *HCtbl)
 
void SetDCtable (G4DCtable *DCtbl)
 
void SetRandomNumberStatus (G4String &st)
 
void StoreEvent (G4Event *evt)
 
const std::vector< const G4Event * > * GetEventVector () const
 

Private Member Functions

void Print (const std::vector< G4String > &title, const std::map< G4int, std::vector< G4double > > &out, G4String &) const
 

Private Attributes

std::map< G4int, G4THitsMap< G4double > *> fMap
 

Additional Inherited Members

- Protected Attributes inherited from G4Run
G4int runID
 
G4int numberOfEvent
 
G4int numberOfEventToBeProcessed
 
G4HCtableHCtable
 
G4DCtableDCtable
 
G4String randomNumberStatus
 
std::vector< const G4Event * > * eventVector
 

Detailed Description

Definition at line 43 of file ElectronRun.hh.

Constructor & Destructor Documentation

◆ ElectronRun()

ElectronRun::ElectronRun ( )

Definition at line 40 of file ElectronRun.cc.

41 : G4Run(), fMap()
42 {
43  fMap[0] = new G4THitsMap<G4double>("MyDetector", "cell flux");
44  fMap[1] = new G4THitsMap<G4double>("MyDetector", "e cell flux");
45  fMap[2] = new G4THitsMap<G4double>("MyDetector", "population");
46  fMap[3] = new G4THitsMap<G4double>("MyDetector", "e population");
47 }
G4Run()
Definition: G4Run.cc:34
std::map< G4int, G4THitsMap< G4double > *> fMap
Definition: ElectronRun.hh:58

◆ ~ElectronRun()

ElectronRun::~ElectronRun ( )
virtual

Definition at line 51 of file ElectronRun.cc.

52 {
53  // Important to clean up the map
54  std::map<G4int, G4THitsMap<G4double>* >::iterator iter = fMap.begin();
55 
56  while (iter != fMap.end()) {
57  delete iter->second;
58  iter++;
59  }
60 }
std::map< G4int, G4THitsMap< G4double > *> fMap
Definition: ElectronRun.hh:58

Member Function Documentation

◆ DumpData()

void ElectronRun::DumpData ( G4String outputFileSpec) const

Definition at line 95 of file ElectronRun.cc.

96 {
97  // Titles
98  std::vector<G4String> title;
99  title.push_back("Radius");
100 
101  // Output map - energy binning on x axis, theta on y
102  std::map< G4int, std::vector<G4double> > output;
103 
104  G4int nThetaBins = 233;
105 
106  // Energy bins depends on the number of scorers
107  G4int nEnergyBins = fMap.size();
108 
109  G4int i(0), j(0);
110 
111  // Initialise current to 0 in all bins
112  for (i=0; i<nThetaBins; i++) {
113  for (j=0; j<nEnergyBins; j++) {
114  output[i].push_back(0);
115  }
116  }
117 
118  i=0;
119 
120  // Fill the output map
121  std::map< G4int, G4THitsMap<G4double>* >::const_iterator iter = fMap.begin();
122 
123  while (iter != fMap.end()) {
124  G4THitsMap<G4double>* hitMap = iter->second;
125 
126  title.push_back(hitMap->GetName());
127 
128  std::map<G4int,G4double*>* myMap = hitMap->GetMap();
129 
130  for (j=0; j<nThetaBins; j++) {
131  G4double* current = (*myMap)[j];
132  if (0 != current) output[j][i] = (*current);
133  }
134 
135  i++;
136  iter++;
137  }
138 
139  Print(title, output, outputFileSpec);
140 }
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:68
int G4int
Definition: G4Types.hh:78
void Print(const std::vector< G4String > &title, const std::map< G4int, std::vector< G4double > > &out, G4String &) const
Definition: ElectronRun.cc:144
std::map< G4int, G4THitsMap< G4double > *> fMap
Definition: ElectronRun.hh:58
double G4double
Definition: G4Types.hh:76
subroutine title
Definition: hijing1.383.f:5981
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Merge()

void ElectronRun::Merge ( const G4Run aRun)
virtual

Reimplemented from G4Run.

Definition at line 198 of file ElectronRun.cc.

199 {
200  // This method is called at the end of the run for each worker thread.
201  // It accumulates the worker's results into global results.
202  const ElectronRun* localRun = static_cast<const ElectronRun*>(aRun);
203  const std::map< G4int, G4THitsMap<G4double>* >& localMap = localRun->fMap;
204  std::map< G4int, G4THitsMap<G4double>* >::const_iterator iter = localMap.begin();
205  for ( ; iter != localMap.end() ; ++iter)
206  (*(fMap[iter->first])) += (*(iter->second));
207 
208  // This call lets Geant4 maintain overall summary information.
209  G4Run::Merge(aRun);
210 }
virtual void Merge(const G4Run *)
Definition: G4Run.cc:54
std::map< G4int, G4THitsMap< G4double > *> fMap
Definition: ElectronRun.hh:58
Here is the call graph for this function:

◆ Print()

void ElectronRun::Print ( const std::vector< G4String > &  title,
const std::map< G4int, std::vector< G4double > > &  out,
G4String outputFileSpec 
) const
private

Definition at line 144 of file ElectronRun.cc.

147 {
148  // Print to G4cout and an output file
149  std::ofstream outFile(outputFileSpec);
150 
151  // Print title vector
152  std::vector<G4String>::const_iterator titleIter = title.begin();
153 
154  while (titleIter != title.end()) {
155  G4cout << std::setw(8)<<*titleIter<<" ";
156  titleIter++;
157  }
158 
159  G4cout<<G4endl;
160 
161  // Print current data
162  std::map< G4int, std::vector<G4double> >::const_iterator iter = myMap.begin();
163 
164  while (iter != myMap.end()) {
165  G4cout << std::setw(8)<<std::setprecision(3)<< iter->first<<" ";
166 
167  std::vector<G4double>::const_iterator energyBinIter = iter->second.begin();
168 
169  // The built-in scorers do not automatically account for the area of the
170  // cylinder replica rings. We must account for this now by multiplying our result
171  // by the ratio of the area of the full cylinder end over the area of the actual
172  // scoring ring.
173  // In this ratio, PI divides out, as does the width of the scoring rings.
174  // Left with just the number of rings squared divided by the ring index plus
175  // 1 squared minus ring index squared.
176  G4int ringNum = iter->first;
177  G4double areaCorrection = 233.*233. /
178  ( (ringNum+1)*(ringNum+1) - ringNum*ringNum );
179  G4int counter = 0;
180 
181  while (energyBinIter != iter->second.end()) {
182  G4double value = *energyBinIter;
183  if (counter < 2) value = value*areaCorrection;
184  G4cout << std::setw(10)<<std::setprecision(5)<< value*mm*mm<<" ";
185  outFile << value*mm*mm;
186  if (counter < 3) outFile <<",";
187  counter++;
188  energyBinIter++;
189  }
190  outFile<<G4endl;
191  G4cout<<G4endl;
192  iter++;
193  }
194 }
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
subroutine title
Definition: hijing1.383.f:5981
static const double mm
Definition: G4SIunits.hh:114
Here is the caller graph for this function:

◆ RecordEvent()

void ElectronRun::RecordEvent ( const G4Event anEvent)
virtual

Reimplemented from G4Run.

Definition at line 64 of file ElectronRun.cc.

65 {
66  // Get the hits collection
67  G4HCofThisEvent* eventHitCollection = anEvent->GetHCofThisEvent();
68 
69  if (!eventHitCollection) return;
70 
71  // Update our private fMap
72  std::map< G4int, G4THitsMap<G4double>* >::iterator iter = fMap.begin();
73 
74  while (iter != fMap.end()) {
75  G4int id = iter->first;
76 
77  // Get the hit collection corresponding to "id"
78  G4THitsMap<G4double>* eventHitsMap
79  = dynamic_cast< G4THitsMap<G4double>* >(eventHitCollection->GetHC(id));
80 
81  // Expect this to exist
82  assert (0 != eventHitsMap);
83 
84  // Accumulate event data into our G4THitsMap<G4double> map
85  *(iter->second) += *eventHitsMap;
86 
87  iter++;
88  }
89 
90  G4Run::RecordEvent(anEvent);
91 }
int G4int
Definition: G4Types.hh:78
G4HCofThisEvent * GetHCofThisEvent() const
Definition: G4Event.hh:185
virtual void RecordEvent(const G4Event *)
Definition: G4Run.cc:51
std::map< G4int, G4THitsMap< G4double > *> fMap
Definition: ElectronRun.hh:58
Here is the call graph for this function:

Member Data Documentation

◆ fMap

std::map<G4int, G4THitsMap<G4double>* > ElectronRun::fMap
private

Definition at line 58 of file ElectronRun.hh.


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