ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ElectronRun.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file ElectronRun.cc
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 //
29 
30 #include "ElectronRun.hh"
32 #include "G4SDManager.hh"
33 #include "G4VPrimitiveScorer.hh"
34 #include "G4SystemOfUnits.hh"
35 #include <assert.h>
36 
37 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
38 
40 : G4Run(), fMap()
41 {
42  fMap[0] = new G4THitsMap<G4double>("MyDetector", "cell flux");
43  fMap[1] = new G4THitsMap<G4double>("MyDetector", "e cell flux");
44  fMap[2] = new G4THitsMap<G4double>("MyDetector", "population");
45  fMap[3] = new G4THitsMap<G4double>("MyDetector", "e population");
46 }
47 
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49 
51 {
52  // Important to clean up the map
53  std::map<G4int, G4THitsMap<G4double>* >::iterator iter = fMap.begin();
54 
55  while (iter != fMap.end()) {
56  delete iter->second;
57  iter++;
58  }
59 }
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62 
63 void ElectronRun::RecordEvent(const G4Event* anEvent)
64 {
65  // Get the hits collection
66  G4HCofThisEvent* eventHitCollection = anEvent->GetHCofThisEvent();
67 
68  if (!eventHitCollection) return;
69 
70  // Update our private fMap
71  std::map< G4int, G4THitsMap<G4double>* >::iterator iter = fMap.begin();
72 
73  while (iter != fMap.end()) {
74  G4int id = iter->first;
75 
76  // Get the hit collection corresponding to "id"
77  G4THitsMap<G4double>* eventHitsMap
78  = dynamic_cast< G4THitsMap<G4double>* >(eventHitCollection->GetHC(id));
79 
80  // Expect this to exist
81  assert (0 != eventHitsMap);
82 
83  // Accumulate event data into our G4THitsMap<G4double> map
84  *(iter->second) += *eventHitsMap;
85 
86  iter++;
87  }
88 
89  G4Run::RecordEvent(anEvent);
90 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93 
94 void ElectronRun::DumpData(G4String &outputFileSpec) const
95 {
96  // Titles
97  std::vector<G4String> title;
98  title.push_back("Radius");
99 
100  // Output map - energy binning on x axis, theta on y
101  std::map< G4int, std::vector<G4double> > output;
102 
103  G4int nThetaBins = 233;
104 
105  // Energy bins depends on the number of scorers
106  G4int nEnergyBins = fMap.size();
107 
108  G4int i(0), j(0);
109 
110  // Initialise current to 0 in all bins
111  for (i=0; i<nThetaBins; i++) {
112  for (j=0; j<nEnergyBins; j++) {
113  output[i].push_back(0);
114  }
115  }
116 
117  i=0;
118 
119  // Fill the output map
120  std::map< G4int, G4THitsMap<G4double>* >::const_iterator iter = fMap.begin();
121 
122  while (iter != fMap.end()) {
123  G4THitsMap<G4double>* hitMap = iter->second;
124 
125  title.push_back(hitMap->GetName());
126 
127  std::map<G4int,G4double*>* myMap = hitMap->GetMap();
128 
129  for (j=0; j<nThetaBins; j++) {
130  G4double* current = (*myMap)[j];
131  if (0 != current) output[j][i] = (*current);
132  }
133 
134  i++;
135  iter++;
136  }
137 
138  Print(title, output, outputFileSpec);
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142 
143 void ElectronRun::Print(const std::vector<G4String>& title,
144  const std::map< G4int, std::vector<G4double> >&myMap,
145  G4String &outputFileSpec) const
146 {
147  // Print to G4cout and an output file
148  std::ofstream outFile(outputFileSpec);
149 
150  // Print title vector
151  std::vector<G4String>::const_iterator titleIter = title.begin();
152 
153  while (titleIter != title.end()) {
154  G4cout << std::setw(8)<<*titleIter<<" ";
155  titleIter++;
156  }
157 
158  G4cout<<G4endl;
159 
160  // Print current data
161  std::map< G4int, std::vector<G4double> >::const_iterator iter = myMap.begin();
162 
163  while (iter != myMap.end()) {
164  G4cout << std::setw(8)<<std::setprecision(3)<< iter->first<<" ";
165 
166  std::vector<G4double>::const_iterator energyBinIter = iter->second.begin();
167 
168  // The built-in scorers do not automatically account for the area of the
169  // cylinder replica rings. We must account for this now by multiplying our result
170  // by the ratio of the area of the full cylinder end over the area of the actual
171  // scoring ring.
172  // In this ratio, PI divides out, as does the width of the scoring rings.
173  // Left with just the number of rings squared divided by the ring index plus
174  // 1 squared minus ring index squared.
175  G4int ringNum = iter->first;
176  G4double areaCorrection = 233.*233. /
177  ( (ringNum+1)*(ringNum+1) - ringNum*ringNum );
178  G4int counter = 0;
179 
180  while (energyBinIter != iter->second.end()) {
181  G4double value = *energyBinIter;
182  if (counter < 2) value = value*areaCorrection;
183  G4cout << std::setw(10)<<std::setprecision(5)<< value*mm*mm<<" ";
184  outFile << value*mm*mm;
185  if (counter < 3) outFile <<",";
186  counter++;
187  energyBinIter++;
188  }
189  outFile<<G4endl;
190  G4cout<<G4endl;
191  iter++;
192  }
193 }
194 
195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
196 
197 void ElectronRun::Merge(const G4Run* aRun)
198 {
199  // This method is called at the end of the run for each worker thread.
200  // It accumulates the worker's results into global results.
201  const ElectronRun* localRun = static_cast<const ElectronRun*>(aRun);
202  const std::map< G4int, G4THitsMap<G4double>* >& localMap = localRun->fMap;
203  std::map< G4int, G4THitsMap<G4double>* >::const_iterator iter = localMap.begin();
204  for ( ; iter != localMap.end() ; ++iter)
205  (*(fMap[iter->first])) += (*(iter->second));
206 
207  // This call lets Geant4 maintain overall summary information.
208  G4Run::Merge(aRun);
209 }
210 
211 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......