Geant4 Cross Reference |
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 /// \file medical/GammaTherapy/src/HistoManager.cc 27 /// \brief Implementation of the HistoManager class 28 // 29 // 30 // 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 33 34 #include "HistoManager.hh" 35 36 #include "DetectorConstruction.hh" 37 38 #include "G4SystemOfUnits.hh" 39 40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 41 42 HistoManager::HistoManager() 43 { 44 fNBinsZ = 60; 45 fNBinsR = 80; 46 fNBinsE = 200; 47 48 fAbsorberZ = 300. * mm; 49 fAbsorberR = 200. * mm; 50 fScoreZ = 100. * mm; 51 fMaxEnergy = 50. * MeV; 52 53 fStepZ = fStepR = fStepE = 0.0; 54 55 Book(); 56 } 57 58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 59 60 HistoManager::~HistoManager() {} 61 62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 63 64 void HistoManager::Book() 65 { 66 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance(); 67 68 // Create or get analysis manager 69 analysisManager->SetDefaultFileType("root"); 70 analysisManager->SetVerboseLevel(1); 71 analysisManager->SetActivation(true); // enable inactivation of histograms 72 73 // Creating an 1-dimensional histograms in the root directory of the tree 74 fHisto.assign(10, 0); 75 int iHisto = 0; 76 fHisto[iHisto] = analysisManager->CreateH1( 77 "10", "Energy deposit at radius (mm) normalised on 1st channel", fNBinsR, 0., fAbsorberR / mm); 78 79 iHisto++; 80 fHisto[iHisto] = analysisManager->CreateH1( 81 "11", "Energy deposit at radius (mm) normalised to integral", fNBinsR, 0., fAbsorberR / mm); 82 83 iHisto++; 84 fHisto[iHisto] = analysisManager->CreateH1( 85 "12", "Energy deposit (MeV/kg/electron) at radius (mm)", fNBinsR, 0., fAbsorberR / mm); 86 87 iHisto++; 88 fHisto[iHisto] = analysisManager->CreateH1("13", "Energy profile (MeV/kg/electron) over Z (mm)", 89 fNBinsZ, 0., fAbsorberZ / mm); 90 91 iHisto++; 92 fHisto[iHisto] = 93 analysisManager->CreateH1("14", "Energy profile (MeV/kg/electron) over Z (mm) at Central Voxel", 94 fNBinsZ, 0., fAbsorberZ / mm); 95 96 iHisto++; 97 fHisto[iHisto] = analysisManager->CreateH1("15", "Energy (MeV) of fGamma produced in the target", 98 fNBinsE, 0., fMaxEnergy / MeV); 99 100 iHisto++; 101 fHisto[iHisto] = analysisManager->CreateH1("16", "Energy (MeV) of fGamma before phantom", fNBinsE, 102 0., fMaxEnergy / MeV); 103 104 iHisto++; 105 fHisto[iHisto] = analysisManager->CreateH1("17", "Energy (MeV) of electrons produced in phantom", 106 fNBinsE, 0., fMaxEnergy / MeV); 107 108 iHisto++; 109 fHisto[iHisto] = analysisManager->CreateH1("18", "Energy (MeV) of electrons produced in target", 110 fNBinsE, 0., fMaxEnergy / MeV); 111 112 iHisto++; 113 fHisto[iHisto] = analysisManager->CreateH1( 114 "19", "Gamma Energy Fluence (MeV/cm2) at radius(mm) in front of phantom", fNBinsR, 0., 115 fAbsorberR / mm); 116 117 // Create all histograms as inactivated 118 // as we have not yet set nbins, vmin, vmax 119 for (int i = 0; i < iHisto + 1; i++) 120 analysisManager->SetH1Activation(i, false); 121 } 122 123 void HistoManager::Update(DetectorConstruction* det, bool bForceActivation) 124 { 125 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance(); 126 127 if (bForceActivation) { 128 for (int i = 0; i < (int)fHisto.size(); i++) 129 analysisManager->SetH1Activation(fHisto[i], true); 130 analysisManager->SetActivation(true); 131 } 132 133 if (analysisManager->IsActive()) { 134 // Check nBinsR / fAbsorberR histograms 135 if (det->GetNumberDivR() != fNBinsR || std::fabs(det->GetAbsorberR() - fAbsorberR) > 0.01 * mm) 136 { 137 fNBinsR = det->GetNumberDivR(); 138 fAbsorberR = det->GetAbsorberR(); 139 std::vector<G4int> histoId{0, 1, 2, 9}; 140 for (auto v : histoId) { 141 analysisManager->SetH1(fHisto[v], fNBinsR, 0., fAbsorberR / mm); 142 } 143 } 144 145 // Check nBinsZ / fAbsorberZ histograms 146 if (det->GetNumberDivZ() != fNBinsZ || std::fabs(det->GetAbsorberZ() - fAbsorberZ) > 0.01 * mm) 147 { 148 fNBinsZ = det->GetNumberDivZ(); 149 fAbsorberZ = det->GetAbsorberZ(); 150 std::vector<G4int> histoId{3, 4}; 151 for (auto v : histoId) { 152 analysisManager->SetH1(fHisto[v], fNBinsZ, 0., fAbsorberZ / mm); 153 } 154 } 155 156 // Check nBinsE / fAbsorberE histograms 157 if (det->GetNumberDivE() != fNBinsE || std::fabs(det->GetMaxEnergy() - fMaxEnergy) > 0.01) { 158 fNBinsE = det->GetNumberDivE(); 159 fMaxEnergy = det->GetMaxEnergy(); 160 std::vector<G4int> histoId{5, 6, 7, 8}; 161 for (auto v : histoId) { 162 analysisManager->SetH1(fHisto[v], fNBinsE, 0., fMaxEnergy / MeV); 163 } 164 } 165 } 166 } 167 168 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 169 170 void HistoManager::DumpHistoParameters() 171 { 172 if (!G4AnalysisManager::Instance()->IsActive()) return; 173 174 for (int i = 0; i < (int)fHisto.size(); i++) { 175 G4int histoId = fHisto[i]; 176 G4String title = G4AnalysisManager::Instance()->GetH1Title(histoId); 177 G4int nbins = G4AnalysisManager::Instance()->GetH1Nbins(histoId); 178 G4double xmin = G4AnalysisManager::Instance()->GetH1Xmin(histoId); 179 G4double xmax = G4AnalysisManager::Instance()->GetH1Xmax(histoId); 180 G4double width = G4AnalysisManager::Instance()->GetH1Width(histoId); 181 G4cout << "Histogram parameters : " << i << " " << histoId << " : " << nbins << " "; 182 G4cout << xmin << "/" << xmax << " " << width << G4endl; 183 } 184 } 185