Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:21:53

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 // This example is provided by the Geant4-DNA collaboration
0027 // dnadamage2 example is derived from the chem6 example
0028 // chem6 example authors: W. G. Shin and S. Incerti (CENBG, France)
0029 //
0030 // Any report or published results obtained using the Geant4-DNA software
0031 // shall cite the following Geant4-DNA collaboration publication:
0032 // J. Appl. Phys. 125 (2019) 104301
0033 // Med. Phys. 45 (2018) e722-e739
0034 // J. Comput. Phys. 274 (2014) 841-882
0035 // Med. Phys. 37 (2010) 4692-4708
0036 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157-178
0037 // The Geant4-DNA web site is available at http://geant4-dna.org
0038 //
0039 // Authors: J. Naoki D. Kondo (UCSF, US)
0040 //          J. Ramos-Mendez and B. Faddegon (UCSF, US)
0041 //
0042 /// \file DetectorConstruction.cc
0043 /// \brief Implementation of the DetectorConstruction class
0044 
0045 #include "DetectorConstruction.hh"
0046 
0047 #include "ScoreLET.hh"
0048 #include "ScoreSpecies.hh"
0049 #include "ScoreStrandBreaks.hh"
0050 
0051 #include "G4Box.hh"
0052 #include "G4GeometryManager.hh"
0053 #include "G4LogicalVolume.hh"
0054 #include "G4LogicalVolumeStore.hh"
0055 #include "G4MultiFunctionalDetector.hh"
0056 #include "G4NistManager.hh"
0057 #include "G4PVPlacement.hh"
0058 #include "G4PhysicalConstants.hh"
0059 #include "G4RunManager.hh"
0060 #include "G4SDManager.hh"
0061 #include "G4SystemOfUnits.hh"
0062 #include "G4VPrimitiveScorer.hh"
0063 #include "G4VisAttributes.hh"
0064 
0065 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0066 
0067 DetectorConstruction::DetectorConstruction() : G4VUserDetectorConstruction()
0068 {
0069   fDetDir = new G4UIdirectory("/det/");
0070   fDetDir->SetGuidance("Detector control.");
0071 
0072   fSizeCmd = new G4UIcmdWithADoubleAndUnit("/det/setSize", this);
0073   fSizeCmd->SetDefaultUnit("um");
0074 
0075   fpOffSetFileUI = new G4UIcmdWithAString("/det/OffSetFile", this);
0076   fpPlasmidNbUI = new G4UIcmdWithAnInteger("/det/NbOfPlasmids", this);
0077 
0078   fpPlasmidFile = new G4UIcmdWithAString("/det/PlasmidFile", this);
0079   fpUseDNA = new G4UIcmdWithABool("/det/UseDNAVolumes", this);
0080 }
0081 
0082 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0083 
0084 DetectorConstruction::~DetectorConstruction()
0085 {
0086   delete fSizeCmd;
0087   delete fpOffSetFileUI;
0088   delete fpPlasmidFile;
0089   delete fpUseDNA;
0090 }
0091 
0092 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0093 
0094 void DetectorConstruction::SetNewValue(G4UIcommand* command, G4String newValue)
0095 {
0096   if (command == fSizeCmd) {
0097     G4double size = fSizeCmd->GetNewDoubleValue(newValue);
0098     SetSize(size);
0099   }
0100 
0101   if (command == fpPlasmidNbUI) {
0102     fNbOfPlasmids = fpPlasmidNbUI->GetNewIntValue(newValue);
0103   }
0104 
0105   if (command == fpOffSetFileUI) {
0106     G4String File = newValue;
0107     ReadOffsetFile(File);
0108   }
0109 
0110   if (command == fpPlasmidFile) {
0111     fPlasmidFile = newValue;
0112   }
0113 
0114   if (command == fpUseDNA) {
0115     fUseDNAVolumes = fpUseDNA->GetNewBoolValue(newValue);
0116   }
0117 }
0118 
0119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0120 
0121 void DetectorConstruction::SetSize(G4double size)
0122 {
0123   G4GeometryManager::GetInstance()->OpenGeometry();
0124   fPlasmidEnvelope->SetRadius(size / 2);
0125   G4RunManager::GetRunManager()->GeometryHasBeenModified();
0126   G4cout << "#### the geometry has been modified " << G4endl;
0127 }
0128 
0129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0130 
0131 G4VPhysicalVolume* DetectorConstruction::Construct()
0132 {
0133   // Clear DNA information Containers
0134   fSampleDNANames.clear();
0135   fSampleDNAPositions.clear();
0136   fSampleDNADetails.clear();
0137   fDNANames.clear();
0138   fDNAPositions.clear();
0139   fDNADetails.clear();
0140 
0141   // Water is defined from NIST material database
0142   G4NistManager* man = G4NistManager::Instance();
0143   G4Material* normalWater = man->FindOrBuildMaterial("G4_WATER");
0144   G4Material* waterWorld =
0145     man->BuildMaterialWithNewDensity("G4_WATER_WORLD", "G4_WATER", 1.0 * g / cm / cm / cm);
0146 
0147   //
0148   // World
0149   G4Box* solidWenvelope = new G4Box("SWorld", 1 * um, 1 * um, 1 * um);
0150   G4LogicalVolume* logicWorld = new G4LogicalVolume(solidWenvelope, waterWorld, "LWorld");
0151   G4VPhysicalVolume* physWorld =
0152     new G4PVPlacement(0, G4ThreeVector(), "PWorld", logicWorld, 0, false, 0);
0153 
0154   // Molecule
0155   fPlasmidEnvelope = new G4Orb("plasmidEnvelope", 0.5 * fWorldSize);
0156 
0157   G4LogicalVolume* tlogicPlasmid =
0158     new G4LogicalVolume(fPlasmidEnvelope, normalWater, "PlasmidVolume");
0159   new G4PVPlacement(0, G4ThreeVector(), tlogicPlasmid, "plasmid", logicWorld, false, 0);
0160 
0161   // World Vis Attributes
0162   G4VisAttributes worldVis(G4Colour(0.0, 0.0, 1.0));
0163   logicWorld->SetVisAttributes(worldVis);
0164 
0165   PhysGeoImport GeoImport = PhysGeoImport();
0166 
0167   G4LogicalVolume* logicStraightVoxel;
0168 
0169   logicStraightVoxel = GeoImport.CreateLogicVolumeXYZ(fPlasmidFile);
0170 
0171   fSampleDNANames = GeoImport.GetMoleculesNames();
0172   fSampleDNAPositions = GeoImport.GetMoleculesPositions();
0173   fSampleDNADetails = GeoImport.GetMoleculesDetails();
0174 
0175   for (int i = 0; i < fNbOfPlasmids; i++) {
0176     if (fUseDNAVolumes)
0177       new G4PVPlacement(0, fVOffset[i], logicStraightVoxel, "VoxelStraight", logicWorld, true, i);
0178     AddDNAInformation(i, fVOffset[i]);
0179   }
0180 
0181   // always return the physical World
0182   return physWorld;
0183 }
0184 
0185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0186 
0187 void DetectorConstruction::ConstructSDandField()
0188 {
0189   G4SDManager::GetSDMpointer()->SetVerboseLevel(1);
0190 
0191   // declare World as a MultiFunctionalDetector scorer
0192   //
0193   G4MultiFunctionalDetector* mfDetector = new G4MultiFunctionalDetector("mfDetector");
0194 
0195   // LET scorer:
0196   //  - scores restricted or unrestricted LET
0197 
0198   ScoreLET* LET = new ScoreLET("LET");
0199   mfDetector->RegisterPrimitive(LET);
0200 
0201   // Species scorer:
0202   //  - scores number of species over time
0203   //  - compute the radiochemical yields (G values)
0204 
0205   G4VPrimitiveScorer* primitivSpecies = new ScoreSpecies("Species");
0206   mfDetector->RegisterPrimitive(primitivSpecies);
0207 
0208   // SB Scorer: Requires access to Geometry
0209   //  - scores number of Direct SB
0210   //  - scores number of Indirect SB
0211 
0212   G4VPrimitiveScorer* primitiveSB = new ScoreStrandBreaks("StrandBreaks", this, &fWorldSize);
0213   mfDetector->RegisterPrimitive(primitiveSB);
0214 
0215   // Attach Detectors to Plasmid Volumes
0216   G4LogicalVolumeStore* theLogicalVolumes = G4LogicalVolumeStore::GetInstance();
0217   for (size_t t = 0; t < theLogicalVolumes->size(); t++) {
0218     G4String lVolumeName = (*theLogicalVolumes)[t]->GetName();
0219     if (lVolumeName != "LWorld") {
0220       (*theLogicalVolumes)[t]->SetSensitiveDetector(mfDetector);
0221     }
0222   }
0223   G4SDManager::GetSDMpointer()->AddNewDetector(mfDetector);
0224 }
0225 
0226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0227 
0228 void DetectorConstruction::ReadOffsetFile(G4String file)
0229 {
0230   if (fVOffset.size() > 0) fVOffset.clear();
0231 
0232   std::ifstream OffSetFile;
0233   OffSetFile.open(file);
0234   G4double x, y, z;
0235 
0236   if (!OffSetFile) {
0237     G4cout << "Plasmid Offset positions file not found!!!" << G4endl;
0238     exit(1);
0239   }
0240 
0241   else {
0242     while (OffSetFile >> x >> y >> z) {
0243       fVOffset.push_back(G4ThreeVector(x * nm, y * nm, z * nm));
0244     }
0245   }
0246 }
0247 
0248 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
0249 
0250 void DetectorConstruction::AddDNAInformation(G4int copy, G4ThreeVector offset)
0251 {
0252   for (size_t i = 0; i < fSampleDNANames.size(); i++) {
0253     fDNANames.push_back(fSampleDNANames[i]);
0254     fDNAPositions.push_back(fSampleDNAPositions[i] + offset);
0255 
0256     std::vector<G4int> Details = fSampleDNADetails[i];
0257     Details[0] = copy;
0258     fDNADetails.push_back(Details);
0259   }
0260 }
0261 
0262 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....