File indexing completed on 2025-02-23 09:20:10
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030 #include <map>
0031
0032 #include "globals.hh"
0033
0034 #include "G4HumanPhantomConstruction.hh"
0035
0036 #include "G4SystemOfUnits.hh"
0037 #include "G4HumanPhantomSD.hh"
0038 #include "G4SDManager.hh"
0039
0040
0041
0042
0043
0044 #include "G4PhantomBuilder.hh"
0045 #include "G4FemaleBuilder.hh"
0046 #include "G4MaleBuilder.hh"
0047 #include "G4PhantomHeadBuilder.hh"
0048 #include "G4RunManager.hh"
0049 #include "G4HumanPhantomMaterial.hh"
0050 #include "G4Box.hh"
0051 #include "G4LogicalVolume.hh"
0052 #include "G4VPhysicalVolume.hh"
0053 #include "G4VisAttributes.hh"
0054 #include "G4Colour.hh"
0055 #include "G4PVPlacement.hh"
0056
0057 G4HumanPhantomConstruction::G4HumanPhantomConstruction()
0058 {
0059 fMessenger = new G4HumanPhantomMessenger(this);
0060 fMaterial = new G4HumanPhantomMaterial();
0061 }
0062
0063 G4HumanPhantomConstruction::~G4HumanPhantomConstruction()
0064 {
0065 delete fMaterial;
0066 delete fMessenger;
0067 }
0068
0069 G4VPhysicalVolume* G4HumanPhantomConstruction::Construct()
0070 {
0071 fMaterial -> DefineMaterials();
0072
0073 G4BasePhantomBuilder* builder = nullptr;
0074
0075 if (fModel == "MIRDHead" || fModel == "ORNLHead")
0076 {
0077 G4cout << "HeadBuilder instantiated" << G4endl;
0078 builder = new G4PhantomHeadBuilder;
0079 if (fModel == "MIRDHead") builder->SetModel("MIRD");
0080 else if (fModel == "ORNLHead") builder->SetModel("ORNLMale");
0081 }
0082 else
0083 {
0084 if (fSex =="Female")
0085 {
0086 builder = new G4FemaleBuilder;
0087 builder->SetModel(fModel);
0088 G4cout << fModel << " "<< fSex << G4endl;
0089 }
0090 else if (fSex == "Male")
0091 {
0092 builder = new G4MaleBuilder;
0093 builder->SetModel(fModel);
0094 }
0095 }
0096
0097 builder->SetMotherVolume(ConstructWorld());
0098
0099
0100
0101 builder->BuildHead("black", false, fSensitivities["Head"]);
0102 builder->BuildSkull("orange", false, fSensitivities["Skull"]);
0103 builder->BuildBrain("yellow", true, fSensitivities["Brain"]);
0104
0105 if (fModel != "MIRDHead" && fModel != "ORNLHead")
0106 {
0107
0108 builder->BuildTrunk("yellow", false, fSensitivities["Trunk"]);
0109
0110 builder->BuildLeftLeg("yellow", false, fSensitivities["LeftLeg"]);
0111 builder->BuildRightLeg("yellow", false, fSensitivities["RightLeg"]);
0112
0113 builder->BuildLeftArmBone("grey", true, fSensitivities["LeftArmBone"]);
0114 builder->BuildRightArmBone("grey", true, fSensitivities["RightArmBone"]);
0115
0116 builder->BuildLeftLegBone("grey", true, fSensitivities["LeftLegBone"]);
0117 builder ->BuildRightLegBone("grey", true, fSensitivities["RightLegBone"]);
0118
0119 builder->BuildUpperSpine("yellow", true, fSensitivities["UpperSpine"]);
0120
0121 if (fModel == "MIRD")
0122 {
0123 builder->BuildLeftScapula("grey", true, fSensitivities["LeftScapula"]);
0124 builder->BuildRightScapula("grey", true, fSensitivities["RightScapula"]);
0125 builder->BuildLeftAdrenal("yellow", true, fSensitivities["LeftAdrenal"]);
0126 builder->BuildRightAdrenal("yellow", true, fSensitivities["RightAdrenal"]);
0127 builder->BuildThymus("orange", true, fSensitivities["Thymus"]);
0128 builder->BuildLeftClavicle("grey", true, fSensitivities["LeftClavicle"]);
0129 builder->BuildRightClavicle("grey", true, fSensitivities["RightClavicle"]);
0130 builder->BuildSmallIntestine("orange", true, fSensitivities["SmallIntestine"]);
0131 builder->BuildRibCage("grey", true, fSensitivities["RibCage"]);
0132 }
0133
0134 builder->BuildMiddleLowerSpine("yellow", true, fSensitivities["MiddleLowerSpine"]);
0135
0136 builder->BuildPelvis("grey", true, fSensitivities["Pelvis"]);
0137
0138 builder->BuildStomach("orange", true, fSensitivities["Stomach"]);
0139 builder->BuildUpperLargeIntestine("lightBlue", true, fSensitivities["UpperLargeIntestine"]);
0140 builder->BuildLowerLargeIntestine("lightBlue", true, fSensitivities["LowerLargeIntestine"]);
0141 builder->BuildSpleen("green", true, fSensitivities["Spleen"]);
0142 builder->BuildPancreas("purple", true, fSensitivities["Pancreas"]);
0143
0144
0145 builder->BuildLeftKidney("green", true, fSensitivities["LeftKidney"]);
0146 builder->BuildRightKidney("green", true, fSensitivities["RightKidney"]);
0147 builder->BuildUrinaryBladder("green", true, fSensitivities["UrinaryBladder"]);
0148
0149
0150
0151
0152
0153
0154 if(fSex=="Female"){
0155
0156 builder->BuildLeftOvary("purple", true, fSensitivities["LeftOvary"]);
0157 builder->BuildRightOvary("purple", true, fSensitivities["RightOvary"]);
0158 builder->BuildUterus("purple", true, fSensitivities["Uterus"]);
0159
0160 if (fModel == "ORNLFemale" || fModel == "MIRD")
0161 {
0162 builder->BuildLeftBreast("purple", true,fSensitivities["LeftBreast"]);
0163 builder->BuildRightBreast("purple", true,fSensitivities["RightBreast"]);
0164 }
0165 }
0166
0167 if(fSex=="Male"){
0168
0169 if (fModel == "MIRD"){
0170 builder -> BuildMaleGenitalia("yellow",false, fSensitivities["MaleGenitalia"]);
0171 builder -> BuildLeftTeste("purple",true, fSensitivities["LeftTeste"]);
0172 builder -> BuildRightTeste("purple",true, fSensitivities["RightTeste"]);
0173 }
0174 else G4cout << "ORNL does not have model for male genitalia and testes yet" << G4endl;
0175 }
0176
0177 }
0178 G4VPhysicalVolume* result=builder->GetPhantom();
0179 delete builder;
0180 return result;
0181 }
0182
0183 void G4HumanPhantomConstruction::SetBodyPartSensitivity(G4String, G4bool)
0184 {
0185 G4cout << "This method is not currently working !!!!" << G4endl;
0186 }
0187
0188 G4VPhysicalVolume* G4HumanPhantomConstruction::ConstructWorld()
0189 {
0190 G4Material* air = fMaterial -> GetMaterial("Air");
0191
0192
0193
0194 G4double worldSize = 1.5 *m ;
0195 G4Box* world = new G4Box("world", worldSize, worldSize, worldSize);
0196
0197 auto* logicWorld = new G4LogicalVolume(world,
0198 air,
0199 "logicalWorld", nullptr, nullptr, nullptr);
0200
0201 G4VPhysicalVolume* motherVolume = new G4PVPlacement(nullptr,G4ThreeVector(),
0202 "physicalWorld",
0203 logicWorld,
0204 nullptr,
0205 false,
0206 0);
0207
0208
0209 auto* WorldVisAtt = new G4VisAttributes(G4Colour(0.94,0.5,0.5));
0210
0211 WorldVisAtt->SetForceSolid(false);
0212 logicWorld->SetVisAttributes(G4VisAttributes::GetInvisible());
0213
0214 return motherVolume;
0215 }
0216
0217 void G4HumanPhantomConstruction::SetPhantomSex(G4String newSex)
0218 {
0219 fSex=newSex;
0220
0221 if (fSex == "Male")
0222 {
0223 G4cout << ">> Male Phantom will be built." << G4endl;
0224 }
0225 if (fSex == "Female")
0226 {
0227 G4cout << ">> Female Phantom will be built." << G4endl;
0228 }
0229 if ((fSex != "Female") && (fSex != "Male"))
0230 G4cout << fSex << " can not be defined!" << G4endl;
0231 }
0232
0233 void G4HumanPhantomConstruction::SetPhantomModel(G4String newModel)
0234 {
0235 fModel = newModel;
0236
0237 if (fModel == "MIRD")
0238 {
0239 G4cout<<" >> Phantom " << fModel << " will be built."<<G4endl;
0240 }
0241 if (fModel == "ORNLFemale")
0242 {
0243 G4cout<<" >> Phantom " << fModel << " will be built."<<G4endl;
0244 }
0245
0246 if (fModel == "ORNLMale")
0247 {
0248 G4cout<<" >> Phantom " << fModel << " will be built."<<G4endl;
0249 }
0250
0251 if (fModel == "MIRDHead")
0252 {
0253 G4cout<<" >> Phantom " << fModel << " will be built."<<G4endl;
0254 }
0255
0256 if (fModel == "ORNLHead")
0257 {
0258 G4cout<<" >> Phantom " << fModel << " will be built."<<G4endl;
0259 }
0260 }
0261
0262 void G4HumanPhantomConstruction::ConstructSDandField()
0263 {
0264 auto* SD = new G4HumanPhantomSD("SD", "HumanPhantomCollection");
0265 G4SDManager::GetSDMpointer()->AddNewDetector(SD);
0266
0267 if (fModel != "ORNLMale" && fModel != "ORNLFemale" && fModel!= "ORNLHead")
0268 {
0269 SetSensitiveDetector("logicalHead",SD);
0270 SetSensitiveDetector("logicalSkull",SD);
0271 SetSensitiveDetector("logicalBrain",SD);
0272 if (fModel != "MIRDHead")
0273 {
0274 SetSensitiveDetector("logicalTrunk",SD);
0275 SetSensitiveDetector("logicalLeftLeg",SD);
0276 SetSensitiveDetector("logicalRightLeg",SD);
0277 SetSensitiveDetector("logicalLeftArmBone",SD);
0278 SetSensitiveDetector("logicalRightArmBone",SD);
0279 SetSensitiveDetector("logicalLeftLegBone",SD);
0280 SetSensitiveDetector("logicalRightLegBone",SD);
0281 SetSensitiveDetector("logicalUpperSpine",SD);
0282 SetSensitiveDetector("logicalLeftScapula",SD);
0283 SetSensitiveDetector("logicalRightScapula",SD);
0284 SetSensitiveDetector("logicalLeftAdrenal",SD);
0285 SetSensitiveDetector("logicalRightAdrenal",SD);
0286 SetSensitiveDetector("logicalThymus",SD);
0287 SetSensitiveDetector("logicalLeftClavicle",SD);
0288 SetSensitiveDetector("logicalRightClavicle",SD);
0289 SetSensitiveDetector("logicalSmallIntestine",SD);
0290 SetSensitiveDetector("logicalRibCage",SD);
0291 SetSensitiveDetector("logicalMiddleLowerSpine",SD);
0292 SetSensitiveDetector("logicalStomach",SD);
0293 SetSensitiveDetector("logicalUpperLargeIntestine",SD);
0294 SetSensitiveDetector("logicalLowerLargeIntestine",SD);
0295 SetSensitiveDetector("logicalSpleen",SD);
0296 SetSensitiveDetector("logicalPancreas",SD);
0297 SetSensitiveDetector("logicalLeftKidney",SD);
0298 SetSensitiveDetector("logicalRightKidney",SD);
0299 SetSensitiveDetector("logicalUrinaryBladder",SD);
0300
0301 if(fSex=="Female"){
0302 SetSensitiveDetector("logicalLeftOvary",SD);
0303 SetSensitiveDetector("logicalRightOvary",SD);
0304 SetSensitiveDetector("logicalUterus",SD);
0305 SetSensitiveDetector("logicalLeftBreast",SD);
0306 SetSensitiveDetector("logicalRightBreast",SD);
0307 }
0308 else if(fSex=="Male"){
0309 SetSensitiveDetector("logicalMaleGenitalia",SD);
0310 SetSensitiveDetector("logicalLeftTeste",SD);
0311 SetSensitiveDetector("logicalRightTeste",SD);
0312 }
0313 }
0314 }else
0315 {
0316 SetSensitiveDetector("HeadVolume",SD);
0317 SetSensitiveDetector("SkullVolume",SD);
0318 SetSensitiveDetector("BrainVolume",SD);
0319 G4cout <<"ORNL model!!!! Head is sensitive only!!!" << G4endl;
0320 }
0321 }