Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
HadrontherapyInteractionParameters.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 // This is the *BASIC* version of Hadrontherapy, a Geant4-based application
27 // See more at: http://g4advancedexamples.lngs.infn.it/Examples/hadrontherapy
28 //
29 // Visit the Hadrontherapy web site (http://www.lns.infn.it/link/Hadrontherapy) to request
30 // the *COMPLETE* version of this program, together with its documentation;
31 // Hadrontherapy (both basic and full version) are supported by the Italian INFN
32 // Institute in the framework of the MC-INFN Group
33 //
34 
35 #include <fstream>
36 #include <iostream>
37 #include <sstream>
38 #include <cmath>
39 #include <vector>
40 
44 
45 #include "globals.hh"
46 #include "G4SystemOfUnits.hh"
47 #include "G4UnitsTable.hh"
48 #include "G4UImanager.hh"
49 #include "G4RunManager.hh"
50 #include "G4LossTableManager.hh"
51 #include "G4Material.hh"
52 #include "G4MaterialCutsCouple.hh"
53 #include "G4ParticleDefinition.hh"
54 #include "G4ParticleTable.hh"
55 #include "G4NistManager.hh"
56 #include "G4Element.hh"
57 #include "G4StateManager.hh"
58 
60  nistEle(new G4NistElementBuilder(0)),
61  nistMat(new G4NistMaterialBuilder(nistEle, 0)),
62  data(G4cout.rdbuf()),
63  pMessenger(0),
64  beamFlag(false)
65 #ifdef G4ANALYSIS_USE_ROOT
66  ,theRootCanvas(0),
67  theRootGraph(0)
68 #endif
69 {
70  if (wantMessenger) pMessenger = new HadrontherapyParameterMessenger(this);
71 }
72 
74 {
75  if (pMessenger) delete pMessenger;
76  delete nistMat;
77  delete nistEle;
78 }
79 
81  const G4ParticleDefinition* pDef,
82  const G4Material* pMat,
83  G4double dens)
84 {
85  if (dens) return ComputeTotalDEDX(ene, pDef, pMat)/dens;
86  return ComputeTotalDEDX(ene, pDef, pMat);
87 }
89 {
90  // Check arguments
91  if ( !ParseArg(vararg)) return false;
92  // Clear previous energy & mass sp vectors
93  energy.clear();
94  massDedx.clear();
95  // log scale
96  if (kinEmin != kinEmax && npoints >1)
97  {
98  G4double logmin = std::log10(kinEmin);
99  G4double logmax = std::log10(kinEmax);
100  G4double en;
101  // uniform log space
102  for (G4double c = 0.; c < npoints; c++)
103  {
104  en = std::pow(10., logmin + ( c*(logmax-logmin) / (npoints - 1.)) );
105  energy.push_back(en/MeV);
106  dedxtot = ComputeTotalDEDX (en, particle, material);
107  massDedx.push_back ( (dedxtot / density)/(MeV*cm2/g) );
108  }
109  }
110  else // one point only
111  {
112  energy.push_back(kinEmin/MeV);
113  dedxtot = ComputeTotalDEDX (kinEmin, particle, material);
114  massDedx.push_back ( (dedxtot / density)/(MeV*cm2/g) );
115  }
116 
117  G4cout.precision(6);
118  data << "MeV " << "MeV*cm2/g " << particle << " (into " <<
119  material << ", density = " << G4BestUnit(density,"Volumic Mass") << ")" << G4endl;
120  data << G4endl;
121  data << std::left << std::setfill(' ');
122  for (size_t i=0; i<energy.size(); i++){
123  data << std::setw(16) << energy[i] << massDedx[i] << G4endl;
124  }
125  outfile.close();
126  // This will plot
127 #ifdef G4ANALYSIS_USE_ROOT
128  PlotStopping("pdf");
129 #endif
130 
131 // Info to user
132  G4String ofName = (filename == "") ? "User terminal": filename;
133  G4cout << "User choice:\n";
134  G4cout << "Kinetic energy lower limit= "<< G4BestUnit(kinEmin,"Energy") <<
135  ", Kinetic energy upper limit= " << G4BestUnit(kinEmax,"Energy") <<
136  ", npoints= "<< npoints << ", particle= \"" << particle <<
137  "\", material= \"" << material << "\", filename= \""<<
138  ofName << "\"" << G4endl;
139  return true;
140 }
142 // Save Plot
143 #ifdef G4ANALYSIS_USE_ROOT
144 void HadrontherapyInteractionParameters::PlotStopping(const G4String& filetype)
145 {
146  if (!theRootCanvas)
147  {
148  gROOT->Reset();
149  gROOT->SetStyle("Plain");
150  theRootCanvas = new TCanvas("theRootCanvas","Interaction Parameters",200, 10, 600,400);
151  theRootCanvas -> SetFillColor(20);
152  theRootCanvas -> SetBorderMode(1);
153  theRootCanvas -> SetBorderSize(1);
154  theRootCanvas -> SetFrameBorderMode(0);
155  theRootCanvas -> SetGrid();
156  // Use global pad: root manual pgg 109,...
157  }
158 
159  if (theRootGraph) delete theRootGraph;
160  theRootGraph = new TGraph(energy.size(), &energy[0], &massDedx[0]);
161  //theRootGraph = new TGraph();
162  axisX = theRootGraph -> GetXaxis(),
163  axisY = theRootGraph -> GetYaxis();
164  axisX -> SetTitle("MeV");
165  axisY -> SetTitle("Stopping Power (MeV cm2/g)");
166  //axisX -> SetNdivisions(500,kTRUE);
167  //axisX -> SetTickLength(0.03);
168  //axisX -> SetLabelOffset(2.005);
169  axisX -> SetAxisColor(2);
170  axisY -> SetAxisColor(2);
171  gPad -> SetLogx(1);
172  gPad -> SetLogy(1);
173  theRootGraph -> SetMarkerColor(4);
174  theRootGraph -> SetMarkerStyle(20);// circle
175  theRootGraph -> SetMarkerSize(.5);
176 
177  G4String gName = particle.substr(0, particle.find("[") ); // cut excitation energy
178  gName = gName + "_" + material;
179  G4String fName = "./referenceData/interaction/" + gName + "." + filetype;
180  theRootGraph -> SetTitle(gName);
181  theRootGraph -> Draw("AP");
182  //theRootCanvas -> Update();
183  //theRootCanvas -> Draw();
184  theRootCanvas -> SaveAs(fName);
185 }
186 #endif
187 
188 // Search for user material choice inside G4NistManager database
189 G4Material* HadrontherapyInteractionParameters::GetNistMaterial(G4String mat)
190 {
191  Pmaterial = G4NistManager::Instance()->FindOrBuildMaterial(mat);
192  if (Pmaterial) density = Pmaterial -> GetDensity();
193  return Pmaterial;
194 }
195 // Parse arguments line
197 {
198  kinEmin = kinEmax = npoints = 0.;
199  particle = material = filename = "";
200  // set internal variables
201  std::istringstream strParam(vararg);
202  // TODO here check for number and parameters consistency
203  strParam >> std::skipws >> material >> kinEmin >> kinEmax >> npoints >> particle >> filename;
204  // npoints must be an integer!
205  npoints = std::floor(npoints);
206 
207 // Check that kinEmax >= kinEmin > 0 && npoints >= 1
208 // TODO NIST points and linear scale
209  if (kinEmax == 0. && kinEmin > 0. ) kinEmax = kinEmin;
210  if (kinEmax == 0. && kinEmin == 0. ) kinEmax = kinEmin = 1.*MeV;
211  if (kinEmax < kinEmin)
212  {
213  G4cout << "WARNING: kinEmin must not exceed kinEmax!" << G4endl;
214  G4cout << "Usage: /parameter/command material kinetic Emin kinetic Emax nPoints [particle] [output filename]" << G4endl;
215  return false;
216  }
217  if (npoints < 1) npoints = 1;
218 
219  // check if element/material is into database
220  if (!GetNistMaterial(material) )
221  {
222  G4cout << "WARNING: material \"" << material << "\" doesn't exist in NIST elements/materials"
223  " table [$G4INSTALL/source/materials/src/G4NistMaterialBuilder.cc]" << G4endl;
224  G4cout << "Use command \"/parameter/nist\" to see full materials list" << G4endl;
225  return false;
226  }
227  // Check for particle
228  if (particle == "") particle = "proton"; // default to "proton"
229  else if ( !FindParticle(particle) )
230  {
231  G4cout << "WARNING: Particle \"" << particle << "\" isn't supported." << G4endl;
232  G4cout << "Try the command \"/particle/list\" to get full supported particles list." << G4endl;
233  G4cout << "If you are interested in an ion that isn't in this list you must give it to the particle gun."
234  "\nTry the commands:\n/gun/particle ion"
235  "\n/gun/ion <atomic number> <mass number> <[charge]>" << G4endl << G4endl;
236  return false;
237  }
238  // start physics by forcing a G4RunManager::BeamOn():
239  BeamOn();
240  // Set output file
241  if( filename != "" )
242  {
243  outfile.open(filename,std::ios_base::trunc); // overwrite existing file
244  data.rdbuf(outfile.rdbuf());
245  }
246  else data.rdbuf(G4cout.rdbuf()); // output is G4cout!
247  return true;
248 }
249 // Force physics tables build
251 {
252  // first check if RunManager is above G4State_Idle
254  G4ApplicationState aState = mState -> GetCurrentState();
255  if ( aState <= G4State_Idle && beamFlag == false)
256  {
257  G4cout << "Issuing a G4RunManager::beamOn()... ";
258  G4cout << "Current Run State is " << mState -> GetStateString( aState ) << G4endl;
260  beamFlag = true;
261  }
262 
263 }
264 // print a list of Nist elements and materials
266 {
267 /*
268  $G4INSTALL/source/materials/src/G4NistElementBuilder.cc
269  You can also construct a new material by the ConstructNewMaterial method:
270  see $G4INSTALL/source/materials/src/G4NistMaterialBuilder.cc
271 */
272  // Get simplest full list
273  if (vararg =="list")
274  {
275  const std::vector<G4String>& vec = nistMat -> GetMaterialNames();
276  for (size_t i=0; i<vec.size(); i++)
277  {
278  G4cout << std::setw(12) << std::left << i+1 << vec[i] << G4endl;
279  }
280  G4cout << G4endl;
281  }
282  else if (vararg =="all" || vararg =="simple" || vararg =="compound" || vararg =="hep" )
283  {
284  nistMat -> ListMaterials(vararg);
285  }
286 }
287 
288