Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/medical/dna/dnadamage2/src/DetectorConstruction.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 // This example is provided by the Geant4-DNA collaboration
 27 // dnadamage2 example is derived from the chem6 example
 28 // chem6 example authors: W. G. Shin and S. Incerti (CENBG, France)
 29 //
 30 // Any report or published results obtained using the Geant4-DNA software
 31 // shall cite the following Geant4-DNA collaboration publication:
 32 // J. Appl. Phys. 125 (2019) 104301
 33 // Med. Phys. 45 (2018) e722-e739
 34 // J. Comput. Phys. 274 (2014) 841-882
 35 // Med. Phys. 37 (2010) 4692-4708
 36 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157-178
 37 // The Geant4-DNA web site is available at http://geant4-dna.org
 38 //
 39 // Authors: J. Naoki D. Kondo (UCSF, US)
 40 //          J. Ramos-Mendez and B. Faddegon (UCSF, US)
 41 //
 42 /// \file DetectorConstruction.cc
 43 /// \brief Implementation of the DetectorConstruction class
 44 
 45 #include "DetectorConstruction.hh"
 46 
 47 #include "ScoreLET.hh"
 48 #include "ScoreSpecies.hh"
 49 #include "ScoreStrandBreaks.hh"
 50 
 51 #include "G4Box.hh"
 52 #include "G4GeometryManager.hh"
 53 #include "G4LogicalVolume.hh"
 54 #include "G4LogicalVolumeStore.hh"
 55 #include "G4MultiFunctionalDetector.hh"
 56 #include "G4NistManager.hh"
 57 #include "G4PVPlacement.hh"
 58 #include "G4PhysicalConstants.hh"
 59 #include "G4RunManager.hh"
 60 #include "G4SDManager.hh"
 61 #include "G4SystemOfUnits.hh"
 62 #include "G4VPrimitiveScorer.hh"
 63 #include "G4VisAttributes.hh"
 64 
 65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
 66 
 67 DetectorConstruction::DetectorConstruction() : G4VUserDetectorConstruction()
 68 {
 69   fDetDir = new G4UIdirectory("/det/");
 70   fDetDir->SetGuidance("Detector control.");
 71 
 72   fSizeCmd = new G4UIcmdWithADoubleAndUnit("/det/setSize", this);
 73   fSizeCmd->SetDefaultUnit("um");
 74 
 75   fpOffSetFileUI = new G4UIcmdWithAString("/det/OffSetFile", this);
 76   fpPlasmidNbUI = new G4UIcmdWithAnInteger("/det/NbOfPlasmids", this);
 77 
 78   fpPlasmidFile = new G4UIcmdWithAString("/det/PlasmidFile", this);
 79   fpUseDNA = new G4UIcmdWithABool("/det/UseDNAVolumes", this);
 80 }
 81 
 82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
 83 
 84 DetectorConstruction::~DetectorConstruction()
 85 {
 86   delete fSizeCmd;
 87   delete fpOffSetFileUI;
 88   delete fpPlasmidFile;
 89   delete fpUseDNA;
 90 }
 91 
 92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
 93 
 94 void DetectorConstruction::SetNewValue(G4UIcommand* command, G4String newValue)
 95 {
 96   if (command == fSizeCmd) {
 97     G4double size = fSizeCmd->GetNewDoubleValue(newValue);
 98     SetSize(size);
 99   }
100 
101   if (command == fpPlasmidNbUI) {
102     fNbOfPlasmids = fpPlasmidNbUI->GetNewIntValue(newValue);
103   }
104 
105   if (command == fpOffSetFileUI) {
106     G4String File = newValue;
107     ReadOffsetFile(File);
108   }
109 
110   if (command == fpPlasmidFile) {
111     fPlasmidFile = newValue;
112   }
113 
114   if (command == fpUseDNA) {
115     fUseDNAVolumes = fpUseDNA->GetNewBoolValue(newValue);
116   }
117 }
118 
119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
120 
121 void DetectorConstruction::SetSize(G4double size)
122 {
123   G4GeometryManager::GetInstance()->OpenGeometry();
124   fPlasmidEnvelope->SetRadius(size / 2);
125   G4RunManager::GetRunManager()->GeometryHasBeenModified();
126   G4cout << "#### the geometry has been modified " << G4endl;
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
130 
131 G4VPhysicalVolume* DetectorConstruction::Construct()
132 {
133   // Clear DNA information Containers
134   fSampleDNANames.clear();
135   fSampleDNAPositions.clear();
136   fSampleDNADetails.clear();
137   fDNANames.clear();
138   fDNAPositions.clear();
139   fDNADetails.clear();
140 
141   // Water is defined from NIST material database
142   G4NistManager* man = G4NistManager::Instance();
143   G4Material* normalWater = man->FindOrBuildMaterial("G4_WATER");
144   G4Material* waterWorld =
145     man->BuildMaterialWithNewDensity("G4_WATER_WORLD", "G4_WATER", 1.0 * g / cm / cm / cm);
146 
147   //
148   // World
149   G4Box* solidWenvelope = new G4Box("SWorld", 1 * um, 1 * um, 1 * um);
150   G4LogicalVolume* logicWorld = new G4LogicalVolume(solidWenvelope, waterWorld, "LWorld");
151   G4VPhysicalVolume* physWorld =
152     new G4PVPlacement(0, G4ThreeVector(), "PWorld", logicWorld, 0, false, 0);
153 
154   // Molecule
155   fPlasmidEnvelope = new G4Orb("plasmidEnvelope", 0.5 * fWorldSize);
156 
157   G4LogicalVolume* tlogicPlasmid =
158     new G4LogicalVolume(fPlasmidEnvelope, normalWater, "PlasmidVolume");
159   new G4PVPlacement(0, G4ThreeVector(), tlogicPlasmid, "plasmid", logicWorld, false, 0);
160 
161   // World Vis Attributes
162   G4VisAttributes worldVis(G4Colour(0.0, 0.0, 1.0));
163   logicWorld->SetVisAttributes(worldVis);
164 
165   PhysGeoImport GeoImport = PhysGeoImport();
166 
167   G4LogicalVolume* logicStraightVoxel;
168 
169   logicStraightVoxel = GeoImport.CreateLogicVolumeXYZ(fPlasmidFile);
170 
171   fSampleDNANames = GeoImport.GetMoleculesNames();
172   fSampleDNAPositions = GeoImport.GetMoleculesPositions();
173   fSampleDNADetails = GeoImport.GetMoleculesDetails();
174 
175   for (int i = 0; i < fNbOfPlasmids; i++) {
176     if (fUseDNAVolumes)
177       new G4PVPlacement(0, fVOffset[i], logicStraightVoxel, "VoxelStraight", logicWorld, true, i);
178     AddDNAInformation(i, fVOffset[i]);
179   }
180 
181   // always return the physical World
182   return physWorld;
183 }
184 
185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
186 
187 void DetectorConstruction::ConstructSDandField()
188 {
189   G4SDManager::GetSDMpointer()->SetVerboseLevel(1);
190 
191   // declare World as a MultiFunctionalDetector scorer
192   //
193   G4MultiFunctionalDetector* mfDetector = new G4MultiFunctionalDetector("mfDetector");
194 
195   // LET scorer:
196   //  - scores restricted or unrestricted LET
197 
198   ScoreLET* LET = new ScoreLET("LET");
199   mfDetector->RegisterPrimitive(LET);
200 
201   // Species scorer:
202   //  - scores number of species over time
203   //  - compute the radiochemical yields (G values)
204 
205   G4VPrimitiveScorer* primitivSpecies = new ScoreSpecies("Species");
206   mfDetector->RegisterPrimitive(primitivSpecies);
207 
208   // SB Scorer: Requires access to Geometry
209   //  - scores number of Direct SB
210   //  - scores number of Indirect SB
211 
212   G4VPrimitiveScorer* primitiveSB = new ScoreStrandBreaks("StrandBreaks", this, &fWorldSize);
213   mfDetector->RegisterPrimitive(primitiveSB);
214 
215   // Attach Detectors to Plasmid Volumes
216   G4LogicalVolumeStore* theLogicalVolumes = G4LogicalVolumeStore::GetInstance();
217   for (size_t t = 0; t < theLogicalVolumes->size(); t++) {
218     G4String lVolumeName = (*theLogicalVolumes)[t]->GetName();
219     if (lVolumeName != "LWorld") {
220       (*theLogicalVolumes)[t]->SetSensitiveDetector(mfDetector);
221     }
222   }
223   G4SDManager::GetSDMpointer()->AddNewDetector(mfDetector);
224 }
225 
226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
227 
228 void DetectorConstruction::ReadOffsetFile(G4String file)
229 {
230   if (fVOffset.size() > 0) fVOffset.clear();
231 
232   std::ifstream OffSetFile;
233   OffSetFile.open(file);
234   G4double x, y, z;
235 
236   if (!OffSetFile) {
237     G4cout << "Plasmid Offset positions file not found!!!" << G4endl;
238     exit(1);
239   }
240 
241   else {
242     while (OffSetFile >> x >> y >> z) {
243       fVOffset.push_back(G4ThreeVector(x * nm, y * nm, z * nm));
244     }
245   }
246 }
247 
248 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
249 
250 void DetectorConstruction::AddDNAInformation(G4int copy, G4ThreeVector offset)
251 {
252   for (size_t i = 0; i < fSampleDNANames.size(); i++) {
253     fDNANames.push_back(fSampleDNANames[i]);
254     fDNAPositions.push_back(fSampleDNAPositions[i] + offset);
255 
256     std::vector<G4int> Details = fSampleDNADetails[i];
257     Details[0] = copy;
258     fDNADetails.push_back(Details);
259   }
260 }
261 
262 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
263