Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-02-23 09:22:01

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 //
0027 #include "DNAGeometryMessenger.hh"
0028 
0029 #include "DNAGeometry.hh"
0030 #include "UtilityFunctions.hh"
0031 
0032 #include "G4UIcmdWith3VectorAndUnit.hh"
0033 #include "G4UIcmdWithABool.hh"
0034 #include "G4UIcmdWithADouble.hh"
0035 #include "G4UIcmdWithADoubleAndUnit.hh"
0036 #include "G4UIcmdWithAString.hh"
0037 #include "G4UIcmdWithAnInteger.hh"
0038 #include "G4UIcmdWithoutParameter.hh"
0039 #include "G4UIdirectory.hh"
0040 
0041 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0042 
0043 DNAGeometryMessenger::DNAGeometryMessenger(DNAGeometry* dnaGeometry)
0044   : fpDNAGeometry(dnaGeometry),
0045     fpDNADirectory(new G4UIdirectory("/dnageom/")),
0046     fpNewPlacementVolumeFile(new G4UIcmdWithAString("/dnageom/placementVolume", this)),
0047     fpVoxelPlacementsFile(new G4UIcmdWithAString("/dnageom/definitionFile", this)),
0048     fpVoxelSideLength(new G4UIcmdWith3VectorAndUnit("/dnageom/placementSize", this)),
0049     fpFractalScaling(new G4UIcmdWith3VectorAndUnit("/dnageom/fractalScaling", this)),
0050     fpAnglesAsPi(new G4UIcmdWithABool("/dnageom/setVoxelPlacementAnglesAsMultiplesOfPi", this)),
0051     fpCustomMoleculeSizes(new G4UIcmdWithABool("/dnageom/useCustomMoleculeSizes", this)),
0052     fpAddMoleculeSize(new G4UIcmdWithAString("/dnageom/moleculeSize", this)),
0053     fpCheckDNAOverlaps(new G4UIcmdWithABool("/dnageom/checkOverlaps", this)),
0054     fpVerbosity(new G4UIcmdWithAnInteger("/dnageom/verbose", this)),
0055     fpSmartless(new G4UIcmdWithAnInteger("/dnageom/setSmartVoxels", this)),
0056     fpIntRangeDirect(new G4UIcmdWithADoubleAndUnit("/dnageom/interactionDirectRange", this)),
0057     fpRadicalKillDistance(new G4UIcmdWithADoubleAndUnit("/dnageom/radicalKillDistance", this)),
0058     fpUseHistoneScav(new G4UIcmdWithABool("/dnageom/activateHistoneScavenging", this)),
0059     fpDrawCellVolumes(new G4UIcmdWithABool("/dnageom/drawCellVolumes", this)),
0060     fpTestDirectory(new G4UIdirectory("/dnatests/")),
0061     fpChromosomeTest(new G4UIcmdWithoutParameter("/dnatests/chromosome", this)),
0062     fpBasePairTest(new G4UIcmdWithoutParameter("/dnatests/basepairs", this)),
0063     fpUniqueIDTest(new G4UIcmdWithoutParameter("/dnatests/uniqueid", this))
0064 {
0065   // Fractals
0066   fpDNADirectory->SetGuidance("Commands to control the fractal.");
0067 
0068   fpNewPlacementVolumeFile->SetGuidance("Set a placement volume");
0069   fpNewPlacementVolumeFile->SetGuidance("format: name path");
0070   fpNewPlacementVolumeFile->SetParameterName("name path", false);
0071 
0072   fpVoxelPlacementsFile->SetGuidance("Path to file that defines placement locations");
0073   fpVoxelPlacementsFile->SetParameterName("path", false);
0074 
0075   fpVoxelSideLength->SetGuidance("Side length for each placement (x, y, z)");
0076   fpVoxelSideLength->SetParameterName("xlength", "ylength", "zlength", false);
0077 
0078   fpFractalScaling->SetGuidance("Scaling of XYZ in fractal definition file");
0079   fpFractalScaling->SetParameterName("xlength", "ylength", "zlength", false);
0080 
0081   fpAnglesAsPi->SetGuidance("Take the angles in the voxel placement file as multiples of pi");
0082   fpAnglesAsPi->SetGuidance("E.g. set to true if the angle 0.5 should mean 90 degrees");
0083   fpAnglesAsPi->SetParameterName("true/false angles as pi", false);
0084   fpAnglesAsPi->SetDefaultValue(false);
0085 
0086   fpAnglesAsPi->SetGuidance("Enable custom molecule sizes");
0087   fpAnglesAsPi->SetGuidance("These can now be set via /dnageom/moleculeSize");
0088   fpAnglesAsPi->SetParameterName("true/false custom sizes", false);
0089   fpAnglesAsPi->SetDefaultValue(false);
0090 
0091   fpAddMoleculeSize->SetGuidance("Set a molecule size in angstrom.");
0092   fpAddMoleculeSize->SetGuidance("format: molecule_name x y z");
0093   fpAddMoleculeSize->SetGuidance("E.G.: PHOSPHATE 3 4 5");
0094   fpAddMoleculeSize->SetGuidance("Note: molecule names are case insensitive");
0095   fpAddMoleculeSize->SetParameterName("name x y z units", false);
0096 
0097   // control
0098   fpCheckDNAOverlaps->SetGuidance("Check overlaps in DNA geometry region");
0099   fpCheckDNAOverlaps->SetParameterName("true/false check overlaps", false);
0100 
0101   fpVerbosity->SetGuidance("Verbosity for DNA geometry");
0102   fpVerbosity->SetParameterName("int verbose level", false);
0103 
0104   fpSmartless->SetGuidance("Optimisation value (int) for smart voxels");
0105   fpSmartless->SetGuidance("The G4 default is 2");
0106   fpSmartless->SetParameterName("Optimasation value", false);
0107 
0108   fpIntRangeDirect->SetGuidance("Critical range to start recording SSBs from direct effects");
0109   fpIntRangeDirect->SetParameterName("Range", false);
0110 
0111   fpRadicalKillDistance->SetGuidance("Distance from base pairs at which to kill radicals");
0112   fpRadicalKillDistance->SetParameterName("Range", false);
0113 
0114   fpUseHistoneScav->SetGuidance("Activate Histone scavenging function with default radius");
0115   fpUseHistoneScav->SetGuidance("Radius can be controlled with /dnageom/histoneScavengingRadius");
0116   fpUseHistoneScav->SetParameterName("true/false, set histone scavenging on", false);
0117 
0118   fpDrawCellVolumes->SetGuidance(
0119     "Draw cell/chromosome volumes rather than DNA (makes DNA invisible)");
0120   fpDrawCellVolumes->SetParameterName("true/false draw cell volumes", false);
0121 
0122   // Tests
0123   fpTestDirectory->SetGuidance("Tests of the DNA geometry");
0124 
0125   fpChromosomeTest->SetGuidance("Test Chromosome Placement classes");
0126 
0127   fpBasePairTest->SetGuidance("Test Base Pair indices are correct");
0128 
0129   fpUniqueIDTest->SetGuidance("Test Unique ID Algorithm");
0130 }
0131 
0132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0133 
0134 void DNAGeometryMessenger::SetNewValue(G4UIcommand* command, G4String newValue)
0135 {
0136   // Fractals
0137   if (command == fpNewPlacementVolumeFile.get()) {
0138     std::vector<G4String> cmd = utility::Split(newValue, ' ');
0139     // input validation
0140     if (!(cmd.size() == 2 || cmd.size() == 3)) {
0141       G4cout << "Invalid input. Command takes 2 string parameters only" << G4endl
0142              << "separated by a space, or three parameters" << G4endl
0143              << "where the last is a boolean." << G4endl;
0144     }
0145     else {
0146       G4String vname = (G4String)cmd[0];
0147       G4String path = (G4String)cmd[1];
0148       G4bool twist = false;
0149       if (cmd.size() == 3) {
0150         if ((G4String)cmd[2] == "true") {
0151           twist = true;
0152         }
0153         else if ((G4String)cmd[2] != "false") {
0154           G4ExceptionDescription errmsg;
0155           errmsg << "Invalid input. Third value must be a boolean"
0156                  << "written as true or false." << G4endl;
0157           G4Exception("DNAGeometryMessenger::SetNewValue", "DNAGeometry001", FatalException,
0158                       errmsg);
0159         }
0160       }
0161 
0162       if (utility::Path_exists(path)) {
0163         fpDNAGeometry->AddVoxelFile(vname, path, twist);
0164       }
0165       else {
0166         G4ExceptionDescription errmsg;
0167         errmsg << path << " has wrong path or does not exist :";
0168         G4Exception("DNAGeometryMessenger::SetNewValue", "DNAGeometry002", FatalException, errmsg);
0169       }
0170     }
0171   }
0172   else if (command == fpVoxelPlacementsFile.get()) {
0173     G4String path = newValue;
0174     if (utility::Path_exists(path)) {
0175       fpDNAGeometry->SetFractalFilename(path);
0176     }
0177     else {
0178       G4ExceptionDescription errmsg;
0179       errmsg << "The fractal file path has wrong path or does not exist : " << path << G4endl;
0180       G4Exception("DNAGeometryMessenger::SetNewValue", "DNAGeometry003", FatalException, errmsg);
0181     }
0182   }
0183   else if (command == fpVoxelSideLength.get()) {
0184     fpDNAGeometry->SetVoxelSideLength(G4UIcmdWith3VectorAndUnit::GetNew3VectorValue(newValue));
0185   }
0186   else if (command == fpFractalScaling.get()) {
0187     fpDNAGeometry->SetFractalScaling(G4UIcmdWith3VectorAndUnit::GetNew3VectorValue(newValue));
0188   }
0189   else if (command == fpAnglesAsPi.get()) {
0190     fpDNAGeometry->SetFractalAnglesAsPi(G4UIcmdWithABool::GetNewBoolValue(newValue));
0191   }
0192   else if (command == fpCustomMoleculeSizes.get()) {
0193     fpDNAGeometry->EnableCustomMoleculeSizes(G4UIcmdWithABool::GetNewBoolValue(newValue));
0194   }
0195   else if (command == fpAddMoleculeSize.get()) {
0196     std::vector<G4String> cmd = utility::Split(newValue, ' ');
0197     // input validation
0198     if (cmd.size() == 4) {
0199       G4cout << "Invalid input. Command takes 4  parameters only" << G4endl
0200              << "separated by a space, e.g. phosphate 1 2 3" << G4endl;
0201     }
0202     else {
0203       try {
0204         G4String name = cmd[0];
0205         G4int x_size = std::stod(cmd[1]);
0206         G4int y_size = std::stod(cmd[2]);
0207         G4int z_size = std::stod(cmd[3]);
0208         G4ThreeVector molecule_size = G4ThreeVector(x_size, y_size, z_size);
0209         fpDNAGeometry->AddChangeMoleculeSize(name, molecule_size);
0210       }
0211       catch (const std::invalid_argument& ia) {
0212         G4cerr << "Invalid argument to convert to double: " << ia.what() << G4endl;
0213       }
0214     }
0215   }
0216   // control
0217   else if (command == fpCheckDNAOverlaps.get()) {
0218     fpDNAGeometry->SetOverlaps(G4UIcmdWithABool::GetNewBoolValue(newValue));
0219   }
0220   else if (command == fpVerbosity.get()) {
0221     fpDNAGeometry->SetVerbosity(G4UIcmdWithAnInteger::GetNewIntValue(newValue));
0222   }
0223   else if (command == fpSmartless.get()) {
0224     fpDNAGeometry->SetSmartless(G4UIcmdWithAnInteger::GetNewIntValue(newValue));
0225   }
0226   else if (command == fpIntRangeDirect.get()) {
0227     fpDNAGeometry->SetDirectInteractionRange(
0228       G4UIcmdWithADoubleAndUnit::GetNewDoubleValue(newValue));
0229   }
0230   else if (command == fpRadicalKillDistance.get()) {
0231     fpDNAGeometry->SetRadicalKillDistance(G4UIcmdWithADoubleAndUnit::GetNewDoubleValue(newValue));
0232   }
0233   else if (command == fpUseHistoneScav.get()) {
0234     fpDNAGeometry->SetHistoneScav(G4UIcmdWithABool::GetNewBoolValue(newValue));
0235   }
0236   else if (command == fpDrawCellVolumes.get()) {
0237     fpDNAGeometry->SetDrawCellVolumes(G4UIcmdWithABool::GetNewBoolValue(newValue));
0238   }
0239   else if (command == fpChromosomeTest.get()) {
0240     fpDNAGeometry->ChromosomeTest();
0241   }
0242   else if (command == fpBasePairTest.get()) {
0243     fpDNAGeometry->BasePairIndexTest();
0244   }
0245   else if (command == fpUniqueIDTest.get()) {
0246     fpDNAGeometry->UniqueIDTest();
0247   }
0248 }
0249 
0250 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......