Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/medical/electronScattering2/src/ElectronRun.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 27 /// \file medical/electronScattering2/src/ElectronRun.cc
 28 /// \brief Implementation of the ElectronRun class
 29 
 30 #include "ElectronRun.hh"
 31 
 32 #include "G4MultiFunctionalDetector.hh"
 33 #include "G4SDManager.hh"
 34 #include "G4SystemOfUnits.hh"
 35 #include "G4VPrimitiveScorer.hh"
 36 
 37 #include <assert.h>
 38 
 39 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 40 
 41 ElectronRun::ElectronRun() : 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 
 51 ElectronRun::~ElectronRun()
 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. / ((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......
212