File indexing completed on 2025-02-23 09:21:53
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
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
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
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
0083
0084 DetectorConstruction::~DetectorConstruction()
0085 {
0086 delete fSizeCmd;
0087 delete fpOffSetFileUI;
0088 delete fpPlasmidFile;
0089 delete fpUseDNA;
0090 }
0091
0092
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
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
0130
0131 G4VPhysicalVolume* DetectorConstruction::Construct()
0132 {
0133
0134 fSampleDNANames.clear();
0135 fSampleDNAPositions.clear();
0136 fSampleDNADetails.clear();
0137 fDNANames.clear();
0138 fDNAPositions.clear();
0139 fDNADetails.clear();
0140
0141
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
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
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
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
0182 return physWorld;
0183 }
0184
0185
0186
0187 void DetectorConstruction::ConstructSDandField()
0188 {
0189 G4SDManager::GetSDMpointer()->SetVerboseLevel(1);
0190
0191
0192
0193 G4MultiFunctionalDetector* mfDetector = new G4MultiFunctionalDetector("mfDetector");
0194
0195
0196
0197
0198 ScoreLET* LET = new ScoreLET("LET");
0199 mfDetector->RegisterPrimitive(LET);
0200
0201
0202
0203
0204
0205 G4VPrimitiveScorer* primitivSpecies = new ScoreSpecies("Species");
0206 mfDetector->RegisterPrimitive(primitivSpecies);
0207
0208
0209
0210
0211
0212 G4VPrimitiveScorer* primitiveSB = new ScoreStrandBreaks("StrandBreaks", this, &fWorldSize);
0213 mfDetector->RegisterPrimitive(primitiveSB);
0214
0215
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
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
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