Geant4 11.1.1
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 else
49 {
50 auto index = fIndexMap[key];
51 return fVoxelVector[index];
52 }
53}
54
56 : fpBoundingMesh(&boundingBox)
57 , fResolution((2 * boundingBox.halfSideLengthInY() / pixel))
58{}
59
61
63{
64 auto& pVoxel = GetVoxel(key);
65 return std::get<2>(pVoxel);
66}
67
69{
70 G4cout << "*********PrintMesh::Size : " << fVoxelVector.size() << G4endl;
71 for(const auto& iter : fVoxelVector)
72 {
73 auto data = std::get<2>(iter);
74 G4cout << "Index : " << std::get<0>(iter)
75 << " number of type : " << std::get<2>(iter).size() << G4endl;
76 for(const auto& it : data)
77 {
78 G4cout << "_____________" << it.first->GetName() << " : " << it.second
79 << G4endl;
80 }
81 G4cout << G4endl;
82 }
83 G4cout << G4endl;
84}
85
87{
88 G4int output = 0;
89
90 for(const auto& iter : fVoxelVector)
91 {
92 auto data = std::get<2>(iter);
93 auto it = data.find(type);
94 if(it != data.end())
95 {
96 output += it->second;
97 }
98 }
99 return output;
100}
101
102void G4DNAMesh::PrintVoxel(const Index& index)
103{
104 G4cout << "*********PrintVoxel::";
105 G4cout << " index : " << index
106 << " number of type : " << this->GetVoxelMapList(index).size()
107 << G4endl;
108
109 for(const auto& it : this->GetVoxelMapList(index))
110 {
111 G4cout << "_____________" << it.first->GetName() << " : " << it.second
112 << G4endl;
113 }
114 G4cout << G4endl;
115}
116
117void G4DNAMesh::InitializeVoxel(const Index& index, Data&& mapList)
118{
119 auto& pVoxel = GetVoxel(index);
120 std::get<2>(pVoxel) = std::move(mapList);
121}
122
124{
125 auto xlo = fpBoundingMesh->Getxlo() + index.x * fResolution;
126 auto ylo = fpBoundingMesh->Getylo() + index.y * fResolution;
127 auto zlo = fpBoundingMesh->Getzlo() + index.z * fResolution;
128 auto xhi = fpBoundingMesh->Getxlo() + (index.x + 1) * fResolution;
129 auto yhi = fpBoundingMesh->Getylo() + (index.y + 1) * fResolution;
130 auto zhi = fpBoundingMesh->Getzlo() + (index.z + 1) * fResolution;
131 return G4DNABoundingBox({ xhi, xlo, yhi, ylo, zhi, zlo });
132}
133
135{
136 fIndexMap.clear();
137 fVoxelVector.clear();
138}
139
141{
142 return *fpBoundingMesh;
143}
144
145std::vector<G4DNAMesh::Index> // array is better ?
147{
148 std::vector<Index> neighbors;
149 neighbors.reserve(6);
150 auto xMax = (G4int) (std::floor(
151 (fpBoundingMesh->Getxhi() - fpBoundingMesh->Getxlo()) / fResolution));
152 auto yMax = (G4int) (std::floor(
153 (fpBoundingMesh->Getyhi() - fpBoundingMesh->Getylo()) / fResolution));
154 auto zMax = (G4int) (std::floor(
155 (fpBoundingMesh->Getzhi() - fpBoundingMesh->Getzlo()) / fResolution));
156
157 if(index.x - 1 >= 0)
158 {
159 neighbors.push_back(Index(index.x - 1, index.y, index.z));
160 }
161 if(index.y - 1 >= 0)
162 {
163 neighbors.push_back(Index(index.x, index.y - 1, index.z));
164 }
165 if(index.z - 1 >= 0)
166 {
167 neighbors.push_back(Index(index.x, index.y, index.z - 1));
168 }
169 if(index.x + 1 < xMax)
170 {
171 neighbors.push_back(Index(index.x + 1, index.y, index.z));
172 }
173 if(index.y + 1 < yMax)
174 {
175 neighbors.push_back(Index(index.x, index.y + 1, index.z));
176 }
177 if(index.z + 1 < zMax)
178 {
179 neighbors.push_back(Index(index.x, index.y, index.z + 1));
180 }
181
182 return neighbors;
183}
184
185G4double G4DNAMesh::GetResolution() const { return fResolution; }
186
188{
189 if(!fpBoundingMesh->contains(position))
190 {
191 G4ExceptionDescription exceptionDescription;
192 exceptionDescription << "the position: " << position
193 << " is not in the box : " << *fpBoundingMesh;
194 G4Exception("G4DNAMesh::GetKey", "G4DNAMesh010", FatalErrorInArgument,
195 exceptionDescription);
196 }
197
198 G4int dx =
199 std::floor((position.x() - fpBoundingMesh->Getxlo()) / fResolution);
200 G4int dy =
201 std::floor((position.y() - fpBoundingMesh->Getylo()) / fResolution);
202 G4int dz =
203 std::floor((position.z() - fpBoundingMesh->Getzlo()) / fResolution);
204 if(dx < 0 || dy < 0 || dz < 0)
205 {
206 G4ExceptionDescription exceptionDescription;
207 exceptionDescription << "the old index: " << position
208 << " to new index : " << Index(dx, dx, dx);
209 G4Exception("G4DNAMesh::CheckIndex", "G4DNAMesh015", FatalErrorInArgument,
210 exceptionDescription);
211 }
212 return Index{ dx, dy, dz };
213}
214
216 const G4int& pixels) const
217{
218 G4int xmax = std::floor(
219 (fpBoundingMesh->Getxhi() - fpBoundingMesh->Getxlo()) / fResolution);
220 G4int ymax = std::floor(
221 (fpBoundingMesh->Getyhi() - fpBoundingMesh->Getylo()) / fResolution);
222 G4int zmax = std::floor(
223 (fpBoundingMesh->Getzhi() - fpBoundingMesh->Getzlo()) / fResolution);
224 auto dx = (G4int) (index.x * pixels / xmax);
225 auto dy = (G4int) (index.y * pixels / ymax);
226 auto dz = (G4int) (index.z * pixels / zmax);
227 if(dx < 0 || dy < 0 || dz < 0)
228 {
229 G4ExceptionDescription exceptionDescription;
230 exceptionDescription << "the old index: " << index
231 << " to new index : " << Index(dx, dx, dx);
232 G4Exception("G4DNAMesh::CheckIndex", "G4DNAMesh013", FatalErrorInArgument,
233 exceptionDescription);
234 }
235 return Index{ dx, dy, dz };
236}
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)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
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:117
std::tuple< Index, Box, Data > Voxel
Definition: G4DNAMesh.hh:47
~G4DNAMesh() override
Definition: G4DNAMesh.cc:60
void PrintVoxel(const Index &index)
Definition: G4DNAMesh.cc:102
G4int GetNumberOfType(MolType type) const
Definition: G4DNAMesh.cc:86
Index ConvertIndex(const Index &index, const G4int &) const
Definition: G4DNAMesh.cc:215
void PrintMesh()
Definition: G4DNAMesh.cc:68
Voxel & GetVoxel(const Index &index)
Definition: G4DNAMesh.cc:36
G4double GetResolution() const
Definition: G4DNAMesh.cc:185
void Reset()
Definition: G4DNAMesh.cc:134
Data & GetVoxelMapList(const Index &index)
Definition: G4DNAMesh.cc:62
std::map< MolType, size_t > Data
Definition: G4DNAMesh.hh:46
const G4DNABoundingBox & GetBoundingBox() const
Definition: G4DNAMesh.cc:140
G4DNAMesh(const G4DNABoundingBox &, G4int)
Definition: G4DNAMesh.cc:55
Index GetIndex(const G4ThreeVector &position) const
Definition: G4DNAMesh.cc:187
std::vector< Index > FindNeighboringVoxels(const Index &index) const
Definition: G4DNAMesh.cc:146