Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAMesh.cc
Go to the documentation of this file.
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
30std::ostream& operator<<(std::ostream& stream, const G4VDNAMesh::Index& rhs)
31{
32 stream << "{" << rhs.x << ", " << rhs.y << ", " << rhs.z << "}";
33 return stream;
34}
35
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
54 : fpBoundingMesh(&boundingBox)
55 , fResolution((2 * boundingBox.halfSideLengthInY() / pixel))
56{}
57
59
61{
62 auto& pVoxel = GetVoxel(key);
63 return std::get<2>(pVoxel);
64}
65
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
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
100void 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
115void G4DNAMesh::InitializeVoxel(const Index& index, Data&& mapList)
116{
117 auto& pVoxel = GetVoxel(index);
118 std::get<2>(pVoxel) = std::move(mapList);
119}
120
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
133{
134 fIndexMap.clear();
135 fVoxelVector.clear();
136}
137
139{
140 return *fpBoundingMesh;
141}
142
143std::vector<G4DNAMesh::Index> // array is better ?
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
183G4double G4DNAMesh::GetResolution() const { return fResolution; }
184
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
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}
std::ostream & operator<<(std::ostream &stream, const G4VDNAMesh::Index &rhs)
Definition G4DNAMesh.cc:30
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double Getxlo() const
G4double Getyhi() const
G4double Getylo() const
G4double Getxhi() const
G4bool contains(const G4DNABoundingBox &other) const
G4double Getzlo() const
G4double Getzhi() const
void InitializeVoxel(const Index &key, Data &&mapList)
Definition G4DNAMesh.cc:115
std::tuple< Index, Box, Data > Voxel
Definition G4DNAMesh.hh:47
~G4DNAMesh() override
Definition G4DNAMesh.cc:58
void PrintVoxel(const Index &index)
Definition G4DNAMesh.cc:100
G4int GetNumberOfType(MolType type) const
Definition G4DNAMesh.cc:84
Index ConvertIndex(const Index &index, const G4int &) const
Definition G4DNAMesh.cc:213
void PrintMesh()
Definition G4DNAMesh.cc:66
Voxel & GetVoxel(const Index &index)
Definition G4DNAMesh.cc:36
G4double GetResolution() const
Definition G4DNAMesh.cc:183
void Reset()
Definition G4DNAMesh.cc:132
Data & GetVoxelMapList(const Index &index)
Definition G4DNAMesh.cc:60
std::map< MolType, size_t > Data
Definition G4DNAMesh.hh:46
const G4DNABoundingBox & GetBoundingBox() const
Definition G4DNAMesh.cc:138
G4DNAMesh(const G4DNABoundingBox &, G4int)
Definition G4DNAMesh.cc:53
Index GetIndex(const G4ThreeVector &position) const
Definition G4DNAMesh.cc:185
std::vector< Index > FindNeighboringVoxels(const Index &index) const
Definition G4DNAMesh.cc:144