Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ElectronRun.cc
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 // $Id: ElectronRun.cc 78175 2013-12-04 14:11:29Z gunter $
27 //
28 /// \file medical/electronScattering2/src/ElectronRun.cc
29 /// \brief Implementation of the ElectronRun class
30 
31 #include "ElectronRun.hh"
33 #include "G4SDManager.hh"
34 #include "G4VPrimitiveScorer.hh"
35 #include "G4SystemOfUnits.hh"
36 #include <assert.h>
37 
38 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
39 
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 }
48 
49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50 
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 }
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63 
64 void ElectronRun::RecordEvent(const G4Event* anEvent)
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 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94 
95 void ElectronRun::DumpData(G4String &outputFileSpec) const
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 }
141 
142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
143 
144 void ElectronRun::Print(const std::vector<G4String>& title,
145  const std::map< G4int, std::vector<G4double> >&myMap,
146  G4String &outputFileSpec) const
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 }
195 
196 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
197 
198 void ElectronRun::Merge(const G4Run* aRun)
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 }
211 
212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Definition of the ElectronRun class.
virtual void Merge(const G4Run *)
Definition: G4Run.cc:54
std::ofstream outFile
Definition: GammaRayTel.cc:68
#define assert(x)
Definition: mymalloc.cc:1309
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
Definition: G4Run.hh:46
subroutine title(NA, NB, NCA, NCB)
Definition: dpm25nuc7.f:1744
virtual void RecordEvent(const G4Event *)
Definition: G4Run.cc:51
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:68
const XML_Char int const XML_Char * value
#define G4endl
Definition: G4ios.hh:61
void DumpData(G4String &) const
Definition: ElectronRun.cc:95
virtual ~ElectronRun()
Definition: ElectronRun.cc:51
G4HCofThisEvent * GetHCofThisEvent() const
Definition: G4Event.hh:174
double G4double
Definition: G4Types.hh:76
virtual void RecordEvent(const G4Event *)
Definition: ElectronRun.cc:64
virtual void Merge(const G4Run *)
Definition: ElectronRun.cc:198