Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAMesh Class Reference

#include <G4DNAMesh.hh>

+ Inheritance diagram for G4DNAMesh:

Public Types

using Box = G4DNABoundingBox
 
using MolType = const G4MolecularConfiguration *
 
using Data = std::map< MolType, size_t >
 
using Voxel = std::tuple< Index, Box, Data >
 
using IndexMap = std::unordered_map< Index, G4int, G4VDNAMesh::hashFunc >
 
using VoxelVector = std::vector< Voxel >
 

Public Member Functions

 G4DNAMesh (const G4DNABoundingBox &, G4int)
 
 ~G4DNAMesh () override
 
Index GetIndex (const G4ThreeVector &position) const
 
VoxelGetVoxel (const Index &index)
 
size_t size ()
 
Index ConvertIndex (const Index &index, const G4int &) const
 
std::vector< IndexFindNeighboringVoxels (const Index &index) const
 
void Reset ()
 
DataGetVoxelMapList (const Index &index)
 
auto end ()
 
auto begin ()
 
VoxelVector::const_iterator const_end () const
 
VoxelVector::const_iterator const_begin () const
 
void PrintMesh ()
 
void PrintVoxel (const Index &index)
 
const G4DNABoundingBoxGetBoundingBox () const
 
G4DNABoundingBox GetBoundingBox (const Index &index)
 
G4int GetNumberOfType (MolType type) const
 
void InitializeVoxel (const Index &key, Data &&mapList)
 
G4double GetResolution () const
 
- Public Member Functions inherited from G4VDNAMesh
 G4VDNAMesh ()=default
 
virtual ~G4VDNAMesh ()=default
 

Detailed Description

Definition at line 41 of file G4DNAMesh.hh.

Member Typedef Documentation

◆ Box

Definition at line 44 of file G4DNAMesh.hh.

◆ Data

using G4DNAMesh::Data = std::map<MolType, size_t>

Definition at line 46 of file G4DNAMesh.hh.

◆ IndexMap

using G4DNAMesh::IndexMap = std::unordered_map<Index, G4int, G4VDNAMesh::hashFunc>

Definition at line 48 of file G4DNAMesh.hh.

◆ MolType

Definition at line 45 of file G4DNAMesh.hh.

◆ Voxel

using G4DNAMesh::Voxel = std::tuple<Index, Box, Data>

Definition at line 47 of file G4DNAMesh.hh.

◆ VoxelVector

using G4DNAMesh::VoxelVector = std::vector<Voxel>

Definition at line 49 of file G4DNAMesh.hh.

Constructor & Destructor Documentation

◆ G4DNAMesh()

G4DNAMesh::G4DNAMesh ( const G4DNABoundingBox boundingBox,
G4int  pixel 
)

Definition at line 55 of file G4DNAMesh.cc.

56 : fpBoundingMesh(&boundingBox)
57 , fResolution((2 * boundingBox.halfSideLengthInY() / pixel))
58{}
G4double halfSideLengthInY() const

◆ ~G4DNAMesh()

G4DNAMesh::~G4DNAMesh ( )
override

Definition at line 60 of file G4DNAMesh.cc.

60{ Reset(); }
void Reset()
Definition: G4DNAMesh.cc:134

Member Function Documentation

◆ begin()

auto G4DNAMesh::begin ( )
inline

Definition at line 60 of file G4DNAMesh.hh.

60{ return fVoxelVector.begin(); }

Referenced by G4DNAGillespieDirectMethod::Initialize().

◆ const_begin()

VoxelVector::const_iterator G4DNAMesh::const_begin ( ) const
inline

Definition at line 62 of file G4DNAMesh.hh.

63 {
64 return fVoxelVector.begin();
65 }

◆ const_end()

VoxelVector::const_iterator G4DNAMesh::const_end ( ) const
inline

Definition at line 61 of file G4DNAMesh.hh.

61{ return fVoxelVector.end(); }

◆ ConvertIndex()

G4VDNAMesh::Index G4DNAMesh::ConvertIndex ( const Index index,
const G4int pixels 
) const

Definition at line 215 of file G4DNAMesh.cc.

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}
@ 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
int G4int
Definition: G4Types.hh:85
G4double Getxlo() const
G4double Getyhi() const
G4double Getylo() const
G4double Getxhi() const
G4double Getzlo() const
G4double Getzhi() const

◆ end()

auto G4DNAMesh::end ( )
inline

Definition at line 59 of file G4DNAMesh.hh.

59{ return fVoxelVector.end(); }

Referenced by G4DNAGillespieDirectMethod::Initialize().

◆ FindNeighboringVoxels()

std::vector< G4DNAMesh::Index > G4DNAMesh::FindNeighboringVoxels ( const Index index) const

Definition at line 146 of file G4DNAMesh.cc.

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}

◆ GetBoundingBox() [1/2]

const G4DNABoundingBox & G4DNAMesh::GetBoundingBox ( ) const

Definition at line 140 of file G4DNAMesh.cc.

141{
142 return *fpBoundingMesh;
143}

Referenced by GetVoxel().

◆ GetBoundingBox() [2/2]

G4DNABoundingBox G4DNAMesh::GetBoundingBox ( const Index index)

Definition at line 123 of file G4DNAMesh.cc.

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}

◆ GetIndex()

G4DNAMesh::Index G4DNAMesh::GetIndex ( const G4ThreeVector position) const

Definition at line 187 of file G4DNAMesh.cc.

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}
G4bool contains(const G4DNABoundingBox &other) const

◆ GetNumberOfType()

G4int G4DNAMesh::GetNumberOfType ( G4DNAMesh::MolType  type) const

Definition at line 86 of file G4DNAMesh.cc.

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}

◆ GetResolution()

G4double G4DNAMesh::GetResolution ( ) const

Definition at line 185 of file G4DNAMesh.cc.

185{ return fResolution; }

◆ GetVoxel()

G4DNAMesh::Voxel & G4DNAMesh::GetVoxel ( const Index index)

Definition at line 36 of file G4DNAMesh.cc.

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}
std::tuple< Index, Box, Data > Voxel
Definition: G4DNAMesh.hh:47
std::map< MolType, size_t > Data
Definition: G4DNAMesh.hh:46
const G4DNABoundingBox & GetBoundingBox() const
Definition: G4DNAMesh.cc:140

Referenced by G4DNAGillespieDirectMethod::CreateEvent(), GetVoxelMapList(), and InitializeVoxel().

◆ GetVoxelMapList()

G4DNAMesh::Data & G4DNAMesh::GetVoxelMapList ( const Index index)

Definition at line 62 of file G4DNAMesh.cc.

63{
64 auto& pVoxel = GetVoxel(key);
65 return std::get<2>(pVoxel);
66}
Voxel & GetVoxel(const Index &index)
Definition: G4DNAMesh.cc:36

Referenced by PrintVoxel().

◆ InitializeVoxel()

void G4DNAMesh::InitializeVoxel ( const Index key,
Data &&  mapList 
)

Definition at line 117 of file G4DNAMesh.cc.

118{
119 auto& pVoxel = GetVoxel(index);
120 std::get<2>(pVoxel) = std::move(mapList);
121}

◆ PrintMesh()

void G4DNAMesh::PrintMesh ( )

Definition at line 68 of file G4DNAMesh.cc.

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}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ PrintVoxel()

void G4DNAMesh::PrintVoxel ( const Index index)

Definition at line 102 of file G4DNAMesh.cc.

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}
Data & GetVoxelMapList(const Index &index)
Definition: G4DNAMesh.cc:62

Referenced by G4DNAGillespieDirectMethod::Initialize().

◆ Reset()

void G4DNAMesh::Reset ( )

Definition at line 134 of file G4DNAMesh.cc.

135{
136 fIndexMap.clear();
137 fVoxelVector.clear();
138}

Referenced by ~G4DNAMesh().

◆ size()

size_t G4DNAMesh::size ( )
inline

Definition at line 54 of file G4DNAMesh.hh.

54{ return fVoxelVector.size(); };

The documentation for this class was generated from the following files: