Geant4 Cross Reference |
1 // ******************************************************************** 2 // * License and Disclaimer * 3 // * * 4 // * The Geant4 software is copyright of the Copyright Holders of * 5 // * the Geant4 Collaboration. It is provided under the terms and * 6 // * conditions of the Geant4 Software License, included in the file * 7 // * LICENSE and available at http://cern.ch/geant4/license . These * 8 // * include a list of copyright holders. * 9 // * * 10 // * Neither the authors of this software system, nor their employing * 11 // * institutes,nor the agencies providing financial support for this * 12 // * work make any representation or warranty, express or implied, * 13 // * regarding this software system or assume any liability for its * 14 // * use. Please see the license in the file LICENSE and URL above * 15 // * for the full disclaimer and the limitation of liability. * 16 // * * 17 // * This code implementation is the result of the scientific and * 18 // * technical work of the GEANT4 collaboration. * 19 // * By using, copying, modifying or distributing the software (or * 20 // * any work based on the software) you agree to acknowledge its * 21 // * use in resulting scientific publications, and indicate your * 22 // * acceptance of all terms of the Geant4 Software license. * 23 // ******************************************************************** 24 // 25 #include "G4DNAMesh.hh" 26 #include <algorithm> 27 #include <ostream> 28 #include "G4ITTrackHolder.hh" 29 30 std::ostream& operator<<(std::ostream& stream, const G4VDNAMesh::Index& rhs) 31 { 32 stream << "{" << rhs.x << ", " << rhs.y << ", " << rhs.z << "}"; 33 return stream; 34 } 35 36 G4DNAMesh::Voxel& G4DNAMesh::GetVoxel(const Index& key) 37 { 38 auto iter = fIndexMap.find(key); 39 if(iter == fIndexMap.end()) 40 { 41 auto box = GetBoundingBox(key); 42 Data mapList; 43 G4DNAMesh::Voxel& voxel = 44 fVoxelVector.emplace_back(std::make_tuple(key, box, std::move(mapList))); 45 fIndexMap[key] = G4int(fVoxelVector.size() - 1); 46 return voxel; 47 } 48 49 auto index = fIndexMap[key]; 50 return fVoxelVector[index]; 51 } 52 53 G4DNAMesh::G4DNAMesh(const G4DNABoundingBox& boundingBox, G4int pixel) 54 : fpBoundingMesh(&boundingBox) 55 , fResolution((2 * boundingBox.halfSideLengthInY() / pixel)) 56 {} 57 58 G4DNAMesh::~G4DNAMesh() { Reset(); } 59 60 G4DNAMesh::Data& G4DNAMesh::GetVoxelMapList(const Index& key) 61 { 62 auto& pVoxel = GetVoxel(key); 63 return std::get<2>(pVoxel); 64 } 65 66 void G4DNAMesh::PrintMesh() 67 { 68 G4cout << "*********PrintMesh::Size : " << fVoxelVector.size() << G4endl; 69 for(const auto& iter : fVoxelVector) 70 { 71 auto data = std::get<2>(iter); 72 G4cout << "Index : " << std::get<0>(iter) 73 << " number of type : " << std::get<2>(iter).size() << G4endl; 74 for(const auto& it : data) 75 { 76 G4cout << "_____________" << it.first->GetName() << " : " << it.second 77 << G4endl; 78 } 79 G4cout << G4endl; 80 } 81 G4cout << G4endl; 82 } 83 84 G4int G4DNAMesh::GetNumberOfType(G4DNAMesh::MolType type) const 85 { 86 G4int output = 0; 87 88 for(const auto& iter : fVoxelVector) 89 { 90 auto data = std::get<2>(iter); 91 auto it = data.find(type); 92 if(it != data.end()) 93 { 94 output += it->second; 95 } 96 } 97 return output; 98 } 99 100 void G4DNAMesh::PrintVoxel(const Index& index) 101 { 102 G4cout << "*********PrintVoxel::"; 103 G4cout << " index : " << index 104 << " number of type : " << this->GetVoxelMapList(index).size() 105 << G4endl; 106 107 for(const auto& it : this->GetVoxelMapList(index)) 108 { 109 G4cout << "_____________" << it.first->GetName() << " : " << it.second 110 << G4endl; 111 } 112 G4cout << G4endl; 113 } 114 115 void G4DNAMesh::InitializeVoxel(const Index& index, Data&& mapList) 116 { 117 auto& pVoxel = GetVoxel(index); 118 std::get<2>(pVoxel) = std::move(mapList); 119 } 120 121 G4DNABoundingBox G4DNAMesh::GetBoundingBox(const Index& index) 122 { 123 auto xlo = fpBoundingMesh->Getxlo() + index.x * fResolution; 124 auto ylo = fpBoundingMesh->Getylo() + index.y * fResolution; 125 auto zlo = fpBoundingMesh->Getzlo() + index.z * fResolution; 126 auto xhi = fpBoundingMesh->Getxlo() + (index.x + 1) * fResolution; 127 auto yhi = fpBoundingMesh->Getylo() + (index.y + 1) * fResolution; 128 auto zhi = fpBoundingMesh->Getzlo() + (index.z + 1) * fResolution; 129 return G4DNABoundingBox({ xhi, xlo, yhi, ylo, zhi, zlo }); 130 } 131 132 void G4DNAMesh::Reset() 133 { 134 fIndexMap.clear(); 135 fVoxelVector.clear(); 136 } 137 138 const G4DNABoundingBox& G4DNAMesh::GetBoundingBox() const 139 { 140 return *fpBoundingMesh; 141 } 142 143 std::vector<G4DNAMesh::Index> // array is better ? 144 G4DNAMesh::FindNeighboringVoxels(const Index& index) const 145 { 146 std::vector<Index> neighbors; 147 neighbors.reserve(6); 148 auto xMax = (G4int) (std::floor( 149 (fpBoundingMesh->Getxhi() - fpBoundingMesh->Getxlo()) / fResolution)); 150 auto yMax = (G4int) (std::floor( 151 (fpBoundingMesh->Getyhi() - fpBoundingMesh->Getylo()) / fResolution)); 152 auto zMax = (G4int) (std::floor( 153 (fpBoundingMesh->Getzhi() - fpBoundingMesh->Getzlo()) / fResolution)); 154 155 if(index.x - 1 >= 0) 156 { 157 neighbors.emplace_back(index.x - 1, index.y, index.z); 158 } 159 if(index.y - 1 >= 0) 160 { 161 neighbors.emplace_back(index.x, index.y - 1, index.z); 162 } 163 if(index.z - 1 >= 0) 164 { 165 neighbors.emplace_back(index.x, index.y, index.z - 1); 166 } 167 if(index.x + 1 < xMax) 168 { 169 neighbors.emplace_back(index.x + 1, index.y, index.z); 170 } 171 if(index.y + 1 < yMax) 172 { 173 neighbors.emplace_back(index.x, index.y + 1, index.z); 174 } 175 if(index.z + 1 < zMax) 176 { 177 neighbors.emplace_back(index.x, index.y, index.z + 1); 178 } 179 180 return neighbors; 181 } 182 183 G4double G4DNAMesh::GetResolution() const { return fResolution; } 184 185 G4DNAMesh::Index G4DNAMesh::GetIndex(const G4ThreeVector& position) const 186 { 187 if(!fpBoundingMesh->contains(position)) 188 { 189 G4ExceptionDescription exceptionDescription; 190 exceptionDescription << "the position: " << position 191 << " is not in the box : " << *fpBoundingMesh; 192 G4Exception("G4DNAMesh::GetKey", "G4DNAMesh010", FatalErrorInArgument, 193 exceptionDescription); 194 } 195 196 G4int dx = 197 std::floor((position.x() - fpBoundingMesh->Getxlo()) / fResolution); 198 G4int dy = 199 std::floor((position.y() - fpBoundingMesh->Getylo()) / fResolution); 200 G4int dz = 201 std::floor((position.z() - fpBoundingMesh->Getzlo()) / fResolution); 202 if(dx < 0 || dy < 0 || dz < 0) 203 { 204 G4ExceptionDescription exceptionDescription; 205 exceptionDescription << "the old index: " << position 206 << " to new index : " << Index(dx, dx, dx); 207 G4Exception("G4DNAMesh::CheckIndex", "G4DNAMesh015", FatalErrorInArgument, 208 exceptionDescription); 209 } 210 return Index{ dx, dy, dz }; 211 } 212 213 G4VDNAMesh::Index G4DNAMesh::ConvertIndex(const Index& index, 214 const G4int& pixels) const 215 { 216 G4int xmax = std::floor( 217 (fpBoundingMesh->Getxhi() - fpBoundingMesh->Getxlo()) / fResolution); 218 G4int ymax = std::floor( 219 (fpBoundingMesh->Getyhi() - fpBoundingMesh->Getylo()) / fResolution); 220 G4int zmax = std::floor( 221 (fpBoundingMesh->Getzhi() - fpBoundingMesh->Getzlo()) / fResolution); 222 auto dx = (G4int) (index.x * pixels / xmax); 223 auto dy = (G4int) (index.y * pixels / ymax); 224 auto dz = (G4int) (index.z * pixels / zmax); 225 if(dx < 0 || dy < 0 || dz < 0) 226 { 227 G4ExceptionDescription exceptionDescription; 228 exceptionDescription << "the old index: " << index 229 << " to new index : " << Index(dx, dx, dx); 230 G4Exception("G4DNAMesh::CheckIndex", "G4DNAMesh013", FatalErrorInArgument, 231 exceptionDescription); 232 } 233 return Index{ dx, dy, dz }; 234 } 235