Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/advanced/human_phantom/src/G4HumanPhantomConstruction.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 // Previous authors: G. Guerrieri, S. Guatelli, and M. G. Pia, INFN Genova, Italy
 27 // Authors (since 2007): S. Guatelli, University of Wollongong, Australia
 28 //
 29 
 30 #include <map>
 31 
 32 #include "globals.hh"
 33 
 34 #include "G4HumanPhantomConstruction.hh"
 35 
 36 #include "G4SystemOfUnits.hh"
 37 #include "G4HumanPhantomSD.hh"
 38 #include "G4SDManager.hh"
 39 
 40 //#include "G4VBodyFactory.hh"
 41 //#include "G4MIRDBodyFactory.hh"
 42 //#include "G4ORNLBodyFactory.hh"
 43 
 44 #include "G4PhantomBuilder.hh"
 45 #include "G4FemaleBuilder.hh"
 46 #include "G4MaleBuilder.hh"
 47 #include "G4PhantomHeadBuilder.hh"
 48 #include "G4RunManager.hh"
 49 #include "G4HumanPhantomMaterial.hh"
 50 #include "G4Box.hh"
 51 #include "G4LogicalVolume.hh"
 52 #include "G4VPhysicalVolume.hh"
 53 #include "G4VisAttributes.hh"
 54 #include "G4Colour.hh"
 55 #include "G4PVPlacement.hh"
 56 
 57 G4HumanPhantomConstruction::G4HumanPhantomConstruction()
 58 {
 59   fMessenger = new G4HumanPhantomMessenger(this);
 60   fMaterial = new G4HumanPhantomMaterial();
 61 }
 62 
 63 G4HumanPhantomConstruction::~G4HumanPhantomConstruction()
 64 {
 65   delete fMaterial;
 66   delete fMessenger;
 67 }
 68 
 69 G4VPhysicalVolume* G4HumanPhantomConstruction::Construct()
 70 {
 71   fMaterial -> DefineMaterials();
 72   
 73   G4BasePhantomBuilder*  builder = nullptr;
 74 
 75   if (fModel == "MIRDHead" || fModel == "ORNLHead")
 76     { 
 77       G4cout << "HeadBuilder instantiated" << G4endl;
 78       builder = new G4PhantomHeadBuilder;
 79       if (fModel ==  "MIRDHead") builder->SetModel("MIRD");
 80       else if (fModel ==  "ORNLHead") builder->SetModel("ORNLMale");
 81     }
 82   else
 83     {  
 84       if (fSex =="Female") 
 85   { 
 86     builder = new G4FemaleBuilder;
 87     builder->SetModel(fModel);
 88     G4cout << fModel << " "<< fSex << G4endl;
 89   }
 90       else if (fSex == "Male") 
 91   {
 92     builder = new G4MaleBuilder;
 93     builder->SetModel(fModel);
 94   }
 95     }
 96   
 97   builder->SetMotherVolume(ConstructWorld());
 98   
 99   // the argument indicates the sensitivity of the volume
100   
101   builder->BuildHead("black", false, fSensitivities["Head"]);
102   builder->BuildSkull("orange", false, fSensitivities["Skull"]); 
103   builder->BuildBrain("yellow", true, fSensitivities["Brain"]); 
104 
105   if (fModel != "MIRDHead" && fModel != "ORNLHead")
106     { 
107       //  builder->SetModel(model);
108       builder->BuildTrunk("yellow", false, fSensitivities["Trunk"]);
109       
110       builder->BuildLeftLeg("yellow", false, fSensitivities["LeftLeg"]);
111       builder->BuildRightLeg("yellow", false, fSensitivities["RightLeg"]);
112       
113       builder->BuildLeftArmBone("grey", true, fSensitivities["LeftArmBone"]);
114       builder->BuildRightArmBone("grey", true, fSensitivities["RightArmBone"]);  
115     
116       builder->BuildLeftLegBone("grey", true, fSensitivities["LeftLegBone"]);
117       builder ->BuildRightLegBone("grey", true, fSensitivities["RightLegBone"]);
118  
119       builder->BuildUpperSpine("yellow", true, fSensitivities["UpperSpine"]); 
120       
121       if (fModel == "MIRD") 
122   {
123     builder->BuildLeftScapula("grey", true, fSensitivities["LeftScapula"]); 
124     builder->BuildRightScapula("grey", true, fSensitivities["RightScapula"]);
125     builder->BuildLeftAdrenal("yellow", true, fSensitivities["LeftAdrenal"]);
126     builder->BuildRightAdrenal("yellow", true, fSensitivities["RightAdrenal"]);
127     builder->BuildThymus("orange", true, fSensitivities["Thymus"]); 
128     builder->BuildLeftClavicle("grey", true, fSensitivities["LeftClavicle"]);
129     builder->BuildRightClavicle("grey", true, fSensitivities["RightClavicle"]);
130     builder->BuildSmallIntestine("orange", true, fSensitivities["SmallIntestine"]);
131     builder->BuildRibCage("grey", true, fSensitivities["RibCage"]); 
132   }
133   
134       builder->BuildMiddleLowerSpine("yellow", true, fSensitivities["MiddleLowerSpine"]);
135   
136       builder->BuildPelvis("grey", true, fSensitivities["Pelvis"]); 
137   
138       builder->BuildStomach("orange", true, fSensitivities["Stomach"]); 
139       builder->BuildUpperLargeIntestine("lightBlue", true, fSensitivities["UpperLargeIntestine"]);
140       builder->BuildLowerLargeIntestine("lightBlue", true, fSensitivities["LowerLargeIntestine"]);
141       builder->BuildSpleen("green", true, fSensitivities["Spleen"]);
142       builder->BuildPancreas("purple", true, fSensitivities["Pancreas"]); 
143       //builder->BuildLiver("orange", true,sensitivities["Liver"]); 
144 
145       builder->BuildLeftKidney("green", true, fSensitivities["LeftKidney"]);
146       builder->BuildRightKidney("green", true, fSensitivities["RightKidney"]);
147       builder->BuildUrinaryBladder("green", true, fSensitivities["UrinaryBladder"]);
148  
149       //builder->BuildHeart("red", true,fSensitivities["Hearth"]);// to do MIRD
150      // builder->BuildLeftLung("blue", true, fSensitivities["LeftLung"]);
151       //builder->BuildRightLung("blue", true, fSensitivities["RightLung"]);
152      // builder->BuildThyroid("orange", true, fSensitivities["Thyroid"]); 
153 
154       if(fSex=="Female"){
155 
156   builder->BuildLeftOvary("purple", true, fSensitivities["LeftOvary"]);
157   builder->BuildRightOvary("purple", true, fSensitivities["RightOvary"]);
158   builder->BuildUterus("purple", true, fSensitivities["Uterus"]);
159 
160   if (fModel == "ORNLFemale" || fModel == "MIRD")
161     {
162       builder->BuildLeftBreast("purple", true,fSensitivities["LeftBreast"]); 
163       builder->BuildRightBreast("purple", true,fSensitivities["RightBreast"]);
164     }
165       }
166       
167       if(fSex=="Male"){
168   
169   if (fModel == "MIRD"){ 
170     builder -> BuildMaleGenitalia("yellow",false, fSensitivities["MaleGenitalia"]);
171     builder -> BuildLeftTeste("purple",true, fSensitivities["LeftTeste"]);
172     builder -> BuildRightTeste("purple",true, fSensitivities["RightTeste"]);
173   }
174   else G4cout <<  "ORNL does not have model for male genitalia and testes yet" << G4endl;
175       }
176       
177     }
178   G4VPhysicalVolume* result=builder->GetPhantom(); 
179   delete builder;
180   return result; 
181 }
182 
183 void  G4HumanPhantomConstruction::SetBodyPartSensitivity(G4String, G4bool)
184 {
185   G4cout << "This method is not currently working !!!!" << G4endl;
186 }
187 
188 G4VPhysicalVolume* G4HumanPhantomConstruction::ConstructWorld()
189 {
190   G4Material* air = fMaterial -> GetMaterial("Air");
191 
192   // World Volume
193   // G4double worldSize = 1.*m ;
194   G4double worldSize = 1.5 *m ;
195   G4Box* world = new G4Box("world", worldSize, worldSize, worldSize);
196 
197   auto* logicWorld = new G4LogicalVolume(world, 
198                   air, 
199                   "logicalWorld", nullptr, nullptr, nullptr);
200 
201   G4VPhysicalVolume* motherVolume = new G4PVPlacement(nullptr,G4ThreeVector(),
202                   "physicalWorld",
203                   logicWorld,
204                   nullptr,
205                   false,
206                   0);
207 
208   // Visualization Attributes
209   auto* WorldVisAtt = new G4VisAttributes(G4Colour(0.94,0.5,0.5));
210     
211   WorldVisAtt->SetForceSolid(false);
212   logicWorld->SetVisAttributes(G4VisAttributes::GetInvisible());
213  
214   return motherVolume;
215 }
216 
217 void G4HumanPhantomConstruction::SetPhantomSex(G4String newSex)
218 {
219   fSex=newSex;
220 
221   if (fSex == "Male")
222     {
223       G4cout << ">> Male Phantom will be built." << G4endl;
224     }
225   if (fSex == "Female")
226     {
227       G4cout << ">> Female Phantom will be built." << G4endl;
228     }
229   if ((fSex != "Female") && (fSex != "Male"))
230     G4cout << fSex << " can not be defined!" << G4endl;
231 }
232 
233 void G4HumanPhantomConstruction::SetPhantomModel(G4String newModel)
234 {
235   fModel = newModel;
236 
237   if (fModel == "MIRD")
238     {
239       G4cout<<" >> Phantom " << fModel << " will be built."<<G4endl;
240     }
241   if (fModel == "ORNLFemale")
242     {
243       G4cout<<" >> Phantom " << fModel << " will be built."<<G4endl;
244     }
245 
246   if (fModel == "ORNLMale")
247     {
248       G4cout<<" >> Phantom " << fModel << " will be built."<<G4endl;
249     }
250 
251   if (fModel == "MIRDHead")
252     {
253       G4cout<<" >> Phantom " << fModel << " will be built."<<G4endl;
254     }
255 
256   if (fModel == "ORNLHead")
257     {
258       G4cout<<" >> Phantom " << fModel << " will be built."<<G4endl;
259     }
260 }
261 
262 void G4HumanPhantomConstruction::ConstructSDandField()
263 {
264    auto* SD = new G4HumanPhantomSD("SD", "HumanPhantomCollection");
265    G4SDManager::GetSDMpointer()->AddNewDetector(SD);
266 
267 if (fModel != "ORNLMale" && fModel != "ORNLFemale" && fModel!= "ORNLHead")  
268 {
269   SetSensitiveDetector("logicalHead",SD);
270   SetSensitiveDetector("logicalSkull",SD);
271   SetSensitiveDetector("logicalBrain",SD);
272   if (fModel != "MIRDHead")
273     { 
274       SetSensitiveDetector("logicalTrunk",SD);
275       SetSensitiveDetector("logicalLeftLeg",SD); 
276       SetSensitiveDetector("logicalRightLeg",SD);
277       SetSensitiveDetector("logicalLeftArmBone",SD); 
278       SetSensitiveDetector("logicalRightArmBone",SD);
279       SetSensitiveDetector("logicalLeftLegBone",SD); 
280       SetSensitiveDetector("logicalRightLegBone",SD);
281       SetSensitiveDetector("logicalUpperSpine",SD);
282       SetSensitiveDetector("logicalLeftScapula",SD);
283       SetSensitiveDetector("logicalRightScapula",SD);
284       SetSensitiveDetector("logicalLeftAdrenal",SD);
285       SetSensitiveDetector("logicalRightAdrenal",SD);      
286       SetSensitiveDetector("logicalThymus",SD);      
287       SetSensitiveDetector("logicalLeftClavicle",SD);
288       SetSensitiveDetector("logicalRightClavicle",SD);
289       SetSensitiveDetector("logicalSmallIntestine",SD); 
290       SetSensitiveDetector("logicalRibCage",SD);       
291       SetSensitiveDetector("logicalMiddleLowerSpine",SD); 
292       SetSensitiveDetector("logicalStomach",SD); 
293       SetSensitiveDetector("logicalUpperLargeIntestine",SD);
294       SetSensitiveDetector("logicalLowerLargeIntestine",SD);
295       SetSensitiveDetector("logicalSpleen",SD);
296       SetSensitiveDetector("logicalPancreas",SD);
297       SetSensitiveDetector("logicalLeftKidney",SD);
298       SetSensitiveDetector("logicalRightKidney",SD);       
299       SetSensitiveDetector("logicalUrinaryBladder",SD);
300 
301       if(fSex=="Female"){
302   SetSensitiveDetector("logicalLeftOvary",SD);
303         SetSensitiveDetector("logicalRightOvary",SD); 
304         SetSensitiveDetector("logicalUterus",SD);
305         SetSensitiveDetector("logicalLeftBreast",SD);
306         SetSensitiveDetector("logicalRightBreast",SD); 
307   }
308       else if(fSex=="Male"){
309     SetSensitiveDetector("logicalMaleGenitalia",SD);
310           SetSensitiveDetector("logicalLeftTeste",SD);
311     SetSensitiveDetector("logicalRightTeste",SD);
312   }
313       }
314   }else 
315  { 
316   SetSensitiveDetector("HeadVolume",SD);
317   SetSensitiveDetector("SkullVolume",SD);
318   SetSensitiveDetector("BrainVolume",SD);
319   G4cout <<"ORNL model!!!! Head is sensitive only!!!" << G4endl;
320 }   
321 }
322