Geant4 Cross Reference |
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 // 27 /// \file DNAHit.cc 28 /// \brief Hit class for a hit interacting with a DNA molecule 29 30 #include "DNAHit.hh" 31 32 #include <utility> 33 34 G4ThreadLocal G4Allocator<DNAHit>* MolecularDNAHitAllocator = nullptr; 35 36 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 37 38 DNAHit::DNAHit(const molecule& mol, G4int placement_idx, // ORG 39 G4int chain, G4int strand, // ORG 40 int64_t bp, const G4ThreeVector& pos, 41 const G4ThreeVector& localpos, // dousatsu 42 const G4double& energy, const G4double& d, const G4String& chromo, 43 const G4MolecularConfiguration* radical) 44 : fMoleculeEnum(mol), 45 fPlacementIdx(placement_idx), 46 fChainIdx(chain), 47 fStrandIdx(strand), 48 fBasePairIdx(bp), 49 fPosition(pos), 50 fLocalPosition(localpos), 51 fEnergy(energy), 52 fDistance(d), 53 fChromosome(chromo), 54 fRadical(radical) 55 { 56 // Computed quantities 57 if (fStrandIdx == 0) { 58 if ((fMoleculeEnum == SUGAR) || (fMoleculeEnum == PHOSPHATE)) { 59 fStrand1Rad = fRadical; 60 fStrand1Energy = fEnergy; 61 } 62 else if ((fMoleculeEnum == CYTOSINE) || (fMoleculeEnum == GUANINE) || (fMoleculeEnum == ADENINE) 63 || (fMoleculeEnum == THYMINE)) 64 { 65 fBase1Rad = fRadical; 66 fBP1Energy = fEnergy; 67 } 68 else { 69 G4Exception("DNAHit", "ERR_UNKNOWN_MOLECULE", JustWarning, 70 "Chemical Reaction with unknown molecule"); 71 } 72 } 73 else if (fStrandIdx == 1) { 74 if ((fMoleculeEnum == SUGAR) || (fMoleculeEnum == PHOSPHATE)) { 75 fStrand2Rad = fRadical; 76 fStrand2Energy = fEnergy; 77 } 78 else if ((fMoleculeEnum == CYTOSINE) || (fMoleculeEnum == GUANINE) || (fMoleculeEnum == ADENINE) 79 || (fMoleculeEnum == THYMINE)) 80 { 81 fBase2Rad = fRadical; 82 fBP2Energy = fEnergy; 83 } 84 else { 85 G4Exception("DNAHit", "ERR_UNKNOWN_MOLECULE", JustWarning, "Hit with unknown molecule"); 86 } 87 } 88 } 89 90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 91 92 // This method is used to combine the computed quantities of two hits. 93 // It keeps the settable parameters of the current hit, but the computed 94 // parameters come from both, letting the hit reflect what happened on each 95 // base pair. 96 void DNAHit::AddHit(const DNAHit& right) 97 { 98 this->fStrand1Energy += right.GetStrand1Energy(); 99 this->fStrand2Energy += right.GetStrand2Energy(); 100 this->fBP1Energy += right.GetBP1Energy(); 101 this->fBP2Energy += right.GetBP2Energy(); 102 if (right.GetStrand1Rad() != nullptr) { 103 this->fStrand1Rad = right.GetStrand1Rad(); 104 } 105 if (right.GetBase1Rad() != nullptr) { 106 this->fBase1Rad = right.GetBase1Rad(); 107 } 108 if (right.GetStrand2Rad() != nullptr) { 109 this->fStrand2Rad = right.GetStrand2Rad(); 110 } 111 if (right.GetBase2Rad() != nullptr) { 112 this->fBase2Rad = right.GetBase2Rad(); 113 } 114 } 115 116 DNAHit::DNAHit(const DNAHit& right) 117 { 118 this->SetPlacementIdx(right.GetPlacementIdx()); 119 this->SetMolecule(right.GetMolecule()); 120 this->SetChainIdx(right.GetChainIdx()); 121 this->SetStrandIdx(right.GetStrandIdx()); 122 this->SetBasePairIdx(right.GetBasePairIdx()); 123 this->SetPosition(right.GetPosition()); 124 this->SetLocalPosition(right.GetLocalPosition()); 125 this->SetEnergy(right.GetEnergy()); 126 this->SetDistance(right.GetDistance()); 127 this->SetChromosome(right.GetChromosome()); 128 this->SetRadical(right.GetRadical()); 129 130 // Computed Quantities, no setters. 131 this->fStrand1Energy = right.GetStrand1Energy(); 132 this->fStrand2Energy = right.GetStrand2Energy(); 133 this->fBP1Energy = right.GetBP1Energy(); 134 this->fBP2Energy = right.GetBP2Energy(); 135 this->fStrand1Rad = right.GetStrand1Rad(); 136 this->fBase1Rad = right.GetBase1Rad(); 137 this->fStrand2Rad = right.GetStrand2Rad(); 138 this->fBase2Rad = right.GetBase2Rad(); 139 } 140 141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 142 143 const DNAHit& DNAHit::operator=(const DNAHit& right) 144 { 145 this->SetMolecule(right.GetMolecule()); 146 this->SetPlacementIdx(right.GetPlacementIdx()); 147 this->SetChainIdx(right.GetChainIdx()); 148 this->SetStrandIdx(right.GetStrandIdx()); 149 this->SetBasePairIdx(right.GetBasePairIdx()); 150 this->SetPosition(right.GetPosition()); 151 this->SetLocalPosition(right.GetLocalPosition()); 152 this->SetEnergy(right.GetEnergy()); 153 this->SetDistance(right.GetDistance()); 154 this->SetChromosome(right.GetChromosome()); 155 this->SetRadical(right.GetRadical()); 156 157 this->fStrand1Energy = right.GetStrand1Energy(); 158 this->fStrand2Energy = right.GetStrand2Energy(); 159 this->fBP1Energy = right.GetBP1Energy(); 160 this->fBP2Energy = right.GetBP2Energy(); 161 this->fStrand1Rad = right.GetStrand1Rad(); 162 this->fBase1Rad = right.GetBase1Rad(); 163 this->fStrand2Rad = right.GetStrand2Rad(); 164 this->fBase2Rad = right.GetBase2Rad(); 165 return *this; 166 } 167 168 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 169 170 G4int DNAHit::operator==(const DNAHit& right) const 171 { 172 return (this == &right) ? 1 : 0; 173 } 174 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 175