Geant4-11
Public Types | Public Member Functions | Private Attributes
G4DNAMesh Class Reference

#include <G4DNAMesh.hh>

Public Types

using Index = G4Voxel::Index
 
using Key = unsigned int
 
using VoxelMap = std::map< Key, G4Voxel * >
 

Public Member Functions

VoxelMap::iterator begin ()
 
VoxelMap::const_iterator begin () const
 
VoxelMap::iterator end ()
 
VoxelMap::const_iterator end () const
 
std::vector< IndexFindNeighboringVoxels (const Index &index) const
 
std::vector< IndexFindVoxelNeighbors (const Index &index) const
 
 G4DNAMesh (const G4DNABoundingBox &, G4int)
 
const G4DNABoundingBoxGetBoundingBox () const
 
G4DNABoundingBox GetBoundingBox (const Index &index)
 
Index GetIndex (const G4ThreeVector &position) const
 
Index GetIndex (const Index &index, int) const
 
Index GetIndex (Key key) const
 
Key GetKey (const G4ThreeVector &pos) const
 
Key GetKey (const Index &index) const
 
G4int GetNumberOfType (G4Voxel::MolType type) const
 
G4double GetResolution () const
 
G4VoxelGetVoxel (const Index &index)
 
G4VoxelGetVoxel (Key key)
 
G4Voxel::MapListGetVoxelMapList (const Index &index)
 
G4Voxel::MapListGetVoxelMapList (Key key)
 
void PrintMesh ()
 
void PrintVoxel (const Index &index)
 
void Reset ()
 
void SetVoxelMapList (const Key &key, G4Voxel::MapList &&mapList)
 
size_t size ()
 
 ~G4DNAMesh ()
 

Private Attributes

VoxelMap fMesh
 
const G4DNABoundingBoxfpBoundingMesh
 
G4double fResolution
 

Detailed Description

Definition at line 112 of file G4DNAMesh.hh.

Member Typedef Documentation

◆ Index

Definition at line 115 of file G4DNAMesh.hh.

◆ Key

using G4DNAMesh::Key = unsigned int

Definition at line 116 of file G4DNAMesh.hh.

◆ VoxelMap

using G4DNAMesh::VoxelMap = std::map<Key, G4Voxel*>

Definition at line 117 of file G4DNAMesh.hh.

Constructor & Destructor Documentation

◆ G4DNAMesh()

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

Definition at line 30 of file G4DNAMesh.cc.

31 : fpBoundingMesh(&boundingBox)
32 , fResolution((2 * boundingBox.halfSideLengthInY() / pixel))
33{}
G4double halfSideLengthInY() const
const G4DNABoundingBox * fpBoundingMesh
Definition: G4DNAMesh.hh:153
G4double fResolution
Definition: G4DNAMesh.hh:154

◆ ~G4DNAMesh()

G4DNAMesh::~G4DNAMesh ( )

Definition at line 35 of file G4DNAMesh.cc.

35{ Reset(); }
void Reset()
Definition: G4DNAMesh.cc:198

References Reset().

Member Function Documentation

◆ begin() [1/2]

VoxelMap::iterator G4DNAMesh::begin ( )
inline

Definition at line 138 of file G4DNAMesh.hh.

138{ return fMesh.begin(); }
VoxelMap fMesh
Definition: G4DNAMesh.hh:152

References fMesh.

Referenced by G4DNAGillespieDirectMethod::Initialize().

◆ begin() [2/2]

VoxelMap::const_iterator G4DNAMesh::begin ( ) const
inline

Definition at line 140 of file G4DNAMesh.hh.

140{ return fMesh.begin(); }

References fMesh.

◆ end() [1/2]

VoxelMap::iterator G4DNAMesh::end ( )
inline

Definition at line 137 of file G4DNAMesh.hh.

137{ return fMesh.end(); }

References fMesh.

Referenced by G4DNAGillespieDirectMethod::Initialize().

◆ end() [2/2]

VoxelMap::const_iterator G4DNAMesh::end ( ) const
inline

Definition at line 139 of file G4DNAMesh.hh.

139{ return fMesh.end(); }

References fMesh.

◆ FindNeighboringVoxels()

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

Definition at line 276 of file G4DNAMesh.cc.

277{
278 std::vector<Index> neighbors;
279 // auto key = GetKey(index);
280 auto xMax = (int) (std::floor(
282 auto yMax = (int) (std::floor(
284 auto zMax = (int) (std::floor(
286
287 if(index.x - 1 >= 0)
288 {
289 neighbors.emplace_back(Index(index.x - 1, index.y, index.z));
290 }
291 if(index.y - 1 >= 0)
292 {
293 neighbors.emplace_back(Index(index.x, index.y - 1, index.z));
294 }
295 if(index.z - 1 >= 0)
296 {
297 neighbors.emplace_back(Index(index.x, index.y, index.z - 1));
298 }
299 if(index.x + 1 < xMax)
300 {
301 neighbors.emplace_back(Index(index.x + 1, index.y, index.z));
302 }
303 if(index.y + 1 < yMax)
304 {
305 neighbors.emplace_back(Index(index.x, index.y + 1, index.z));
306 }
307 if(index.z + 1 < zMax)
308 {
309 neighbors.emplace_back(Index(index.x, index.y, index.z + 1));
310 }
311
312#ifdef DEBUG
313 G4cout << "Neighbors of : " << index << G4endl;
314 for(const auto& it : neighbors)
315 {
316 G4cout << it << G4endl;
317 }
318#endif
319
320 if(neighbors.size() > 6)
321 {
322 G4ExceptionDescription exceptionDescription;
323 exceptionDescription << "neighbors.size() > 6";
324 G4Exception("G4DNAMesh::FindVoxelNeighbors", "G4DNAMesh002",
325 FatalErrorInArgument, exceptionDescription);
326 }
327 return neighbors;
328}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double Getxlo() const
G4double Getyhi() const
G4double Getylo() const
G4double Getxhi() const
G4double Getzlo() const
G4double Getzhi() const
G4Voxel::Index Index
Definition: G4DNAMesh.hh:115

References FatalErrorInArgument, fpBoundingMesh, fResolution, G4cout, G4endl, G4Exception(), G4DNABoundingBox::Getxhi(), G4DNABoundingBox::Getxlo(), G4DNABoundingBox::Getyhi(), G4DNABoundingBox::Getylo(), G4DNABoundingBox::Getzhi(), G4DNABoundingBox::Getzlo(), G4Voxel::Index::x, G4Voxel::Index::y, and G4Voxel::Index::z.

Referenced by G4DNAGillespieDirectMethod::DiffusiveJumping().

◆ FindVoxelNeighbors()

std::vector< G4Voxel::Index > G4DNAMesh::FindVoxelNeighbors ( const Index index) const

Definition at line 232 of file G4DNAMesh.cc.

233{
234 std::vector<Index> neighbors;
235
236 auto xMax = (int) (std::floor(
238 auto yMax = (int) (std::floor(
240 auto zMax = (int) (std::floor(
242
243 auto xmin = (index.x - 1) < 0 ? 0 : (index.x - 1);
244 auto ymin = (index.y - 1) < 0 ? 0 : (index.y - 1);
245 auto zmin = (index.z - 1) < 0 ? 0 : (index.z - 1);
246
247 auto xmax = (index.x + 1) > xMax ? xMax : (index.x + 1);
248 auto ymax = (index.y + 1) > yMax ? yMax : (index.y + 1);
249 auto zmax = (index.z + 1) > zMax ? zMax : (index.z + 1);
250 for(int ix = xmin; ix <= xmax; ix++)
251 {
252 for(int iy = ymin; iy <= ymax; iy++)
253 {
254 for(int iz = zmin; iz <= zmax; iz++)
255 {
256 auto key = iz * yMax * xMax + iy * xMax + ix;
257 if(GetIndex(key) != index)
258 { // deleting the middle element
259 neighbors.push_back(GetIndex(key));
260 }
261 }
262 }
263 }
264 if(neighbors.empty())
265 {
266 G4ExceptionDescription exceptionDescription;
267 exceptionDescription << "neighbors.empty()";
268 G4Exception("G4DNAMesh::FindVoxelNeighbors", "G4DNAMesh001",
269 FatalErrorInArgument, exceptionDescription);
270 }
271
272 return neighbors;
273}
Index GetIndex(Key key) const
Definition: G4DNAMesh.cc:164

References FatalErrorInArgument, fpBoundingMesh, fResolution, G4Exception(), GetIndex(), G4DNABoundingBox::Getxhi(), G4DNABoundingBox::Getxlo(), G4DNABoundingBox::Getyhi(), G4DNABoundingBox::Getylo(), G4DNABoundingBox::Getzhi(), G4DNABoundingBox::Getzlo(), G4Voxel::Index::x, G4Voxel::Index::y, and G4Voxel::Index::z.

◆ GetBoundingBox() [1/2]

const G4DNABoundingBox & G4DNAMesh::GetBoundingBox ( ) const

◆ GetBoundingBox() [2/2]

G4DNABoundingBox G4DNAMesh::GetBoundingBox ( const Index index)

Definition at line 152 of file G4DNAMesh.cc.

153{
154 auto xlo = fpBoundingMesh->Getxlo() + index.x * fResolution;
155 auto ylo = fpBoundingMesh->Getylo() + index.y * fResolution;
156 auto zlo = fpBoundingMesh->Getzlo() + index.z * fResolution;
157
158 auto xhi = fpBoundingMesh->Getxlo() + (index.x + 1) * fResolution;
159 auto yhi = fpBoundingMesh->Getylo() + (index.y + 1) * fResolution;
160 auto zhi = fpBoundingMesh->Getzlo() + (index.z + 1) * fResolution;
161 return G4DNABoundingBox({ xhi, xlo, yhi, ylo, zhi, zlo });
162}

References fpBoundingMesh, fResolution, G4DNABoundingBox::Getxlo(), G4DNABoundingBox::Getylo(), G4DNABoundingBox::Getzlo(), G4Voxel::Index::x, G4Voxel::Index::y, and G4Voxel::Index::z.

◆ GetIndex() [1/3]

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

Definition at line 129 of file G4DNAMesh.cc.

130{
131 int dx = std::floor((position.x() - fpBoundingMesh->Getxlo()) / fResolution);
132 int dy = std::floor((position.y() - fpBoundingMesh->Getylo()) / fResolution);
133 int dz = std::floor((position.z() - fpBoundingMesh->Getzlo()) / fResolution);
134 assert(dx >= 0 && dy >= 0 && dz >= 0);
135 return G4Voxel::Index{ dx, dy, dz };
136}

References fpBoundingMesh, fResolution, G4DNABoundingBox::Getxlo(), G4DNABoundingBox::Getylo(), and G4DNABoundingBox::Getzlo().

◆ GetIndex() [2/3]

G4Voxel::Index G4DNAMesh::GetIndex ( const Index index,
int  pixels 
) const

Definition at line 137 of file G4DNAMesh.cc.

138{
139 int xmax =
140 std::floor((fpBoundingMesh->Getxhi() - fpBoundingMesh->Getxlo()) / fResolution);
141 int ymax =
142 std::floor((fpBoundingMesh->Getyhi() - fpBoundingMesh->Getylo()) / fResolution);
143 int zmax =
144 std::floor((fpBoundingMesh->Getzhi() - fpBoundingMesh->Getzlo()) / fResolution);
145 int dx = (int) (index.x * pixels / xmax);
146 int dy = (int) (index.y * pixels / ymax);
147 int dz = (int) (index.z * pixels / zmax);
148 assert(dx >= 0 && dy >= 0 && dz >= 0);
149 return Index{ dx, dy, dz };
150}

References fpBoundingMesh, fResolution, G4DNABoundingBox::Getxhi(), G4DNABoundingBox::Getxlo(), G4DNABoundingBox::Getyhi(), G4DNABoundingBox::Getylo(), G4DNABoundingBox::Getzhi(), G4DNABoundingBox::Getzlo(), G4Voxel::Index::x, G4Voxel::Index::y, and G4Voxel::Index::z.

◆ GetIndex() [3/3]

G4Voxel::Index G4DNAMesh::GetIndex ( Key  key) const

Definition at line 164 of file G4DNAMesh.cc.

165{
166 G4int xmax =
167 std::floor((fpBoundingMesh->Getxhi() - fpBoundingMesh->Getxlo()) / fResolution);
168 G4int ymax =
169 std::floor((fpBoundingMesh->Getyhi() - fpBoundingMesh->Getylo()) / fResolution);
170 G4int id = key;
171 G4int x_ = id % xmax;
172 id /= xmax;
173 G4int y_ = id % ymax;
174 id /= ymax;
175 G4int z_ = id;
176
177 if(xmax != ymax)
178 {
179 G4cout << xmax << " " << ymax << " " << key << G4endl;
180 G4ExceptionDescription exceptionDescription;
181 exceptionDescription << "xmax != ymax";
182 G4Exception("G4DNAMesh::GetIndex", "G4DNAMesh006", FatalErrorInArgument,
183 exceptionDescription);
184 }
185
186 if(x_ < 0 || y_ < 0 || z_ < 0)
187 {
188 G4cout << xmax << " " << ymax << " " << key << G4endl;
189 G4cout << x_ << " " << y_ << " " << z_ << G4endl;
190 G4ExceptionDescription exceptionDescription;
191 exceptionDescription << "x_ < 0 || y_ < 0 || z_ < 0";
192 G4Exception("G4DNAMesh::GetIndex", "G4DNAMesh005", FatalErrorInArgument,
193 exceptionDescription);
194 }
195 return Index{ x_, y_, z_ };
196}
int G4int
Definition: G4Types.hh:85

References FatalErrorInArgument, fpBoundingMesh, fResolution, G4cout, G4endl, G4Exception(), G4DNABoundingBox::Getxhi(), G4DNABoundingBox::Getxlo(), G4DNABoundingBox::Getyhi(), and G4DNABoundingBox::Getylo().

Referenced by G4DNAGillespieDirectMethod::CreateEvent(), FindVoxelNeighbors(), G4DNAGillespieDirectMethod::Initialize(), and SetVoxelMapList().

◆ GetKey() [1/2]

G4DNAMesh::Key G4DNAMesh::GetKey ( const G4ThreeVector pos) const

Definition at line 332 of file G4DNAMesh.cc.

333{
335 {
336 G4ExceptionDescription exceptionDescription;
337 exceptionDescription << "the position: " << position
338 << " is not in the box";
339 G4Exception("G4DNAMesh::GetKey", "G4DNAMesh010", FatalErrorInArgument,
340 exceptionDescription);
341 }
342
343 auto dx = std::floor((position.x() - fpBoundingMesh->Getxlo()) / fResolution);
344 auto dy = std::floor((position.y() - fpBoundingMesh->Getylo()) / fResolution);
345 auto dz = std::floor((position.z() - fpBoundingMesh->Getzlo()) / fResolution);
346 auto xmax =
347 std::floor((fpBoundingMesh->Getxhi() - fpBoundingMesh->Getxlo()) / fResolution);
348 auto ymax =
349 std::floor((fpBoundingMesh->Getyhi() - fpBoundingMesh->Getylo()) / fResolution);
350 return dz * ymax * xmax + dy * xmax + dx;
351}
G4bool contains(const G4DNABoundingBox &other) const

References G4DNABoundingBox::contains(), FatalErrorInArgument, fpBoundingMesh, fResolution, G4Exception(), G4DNABoundingBox::Getxhi(), G4DNABoundingBox::Getxlo(), G4DNABoundingBox::Getyhi(), G4DNABoundingBox::Getylo(), and G4DNABoundingBox::Getzlo().

Referenced by GetVoxel(), GetVoxelMapList(), and PrintVoxel().

◆ GetKey() [2/2]

G4DNAMesh::Key G4DNAMesh::GetKey ( const Index index) const

Definition at line 58 of file G4DNAMesh.cc.

59{
60 auto xmax = (unsigned int) (std::floor(
62 auto ymax = (unsigned int) (std::floor(
64 return index.z * ymax * xmax + index.y * xmax + index.x;
65}

References fpBoundingMesh, fResolution, G4DNABoundingBox::Getxhi(), G4DNABoundingBox::Getxlo(), G4DNABoundingBox::Getyhi(), G4DNABoundingBox::Getylo(), G4Voxel::Index::x, G4Voxel::Index::y, and G4Voxel::Index::z.

◆ GetNumberOfType()

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

Definition at line 78 of file G4DNAMesh.cc.

79{
80 G4int output = 0;
81 auto iter = fMesh.begin();
82 for(; iter != fMesh.end(); iter++)
83 {
84 auto node = dynamic_cast<G4Voxel*>(iter->second);
85 if(node == nullptr)
86 {
87 continue;
88 }
89 auto it = node->GetMapList().find(type);
90 if(it != node->GetMapList().end())
91 {
92 output += it->second;
93 }
94 }
95 return output;
96}
MapList & GetMapList()
Definition: G4DNAMesh.hh:89

References fMesh, and G4Voxel::GetMapList().

Referenced by G4DNAUpdateSystemModel::UpdateSystem().

◆ GetResolution()

G4double G4DNAMesh::GetResolution ( ) const

Definition at line 330 of file G4DNAMesh.cc.

330{ return fResolution; }

References fResolution.

◆ GetVoxel() [1/2]

G4Voxel * G4DNAMesh::GetVoxel ( const Index index)

Definition at line 226 of file G4DNAMesh.cc.

227{
228 return GetVoxel(GetKey(index));
229}
Key GetKey(const G4ThreeVector &pos) const
Definition: G4DNAMesh.cc:332
G4Voxel * GetVoxel(Key key)
Definition: G4DNAMesh.cc:216

References GetKey(), and GetVoxel().

◆ GetVoxel() [2/2]

G4Voxel * G4DNAMesh::GetVoxel ( Key  key)

Definition at line 216 of file G4DNAMesh.cc.

217{
218 auto it = fMesh.find(key);
219 if(it != fMesh.end())
220 {
221 return it->second;
222 }
223 return nullptr;
224}

References fMesh.

Referenced by GetVoxel().

◆ GetVoxelMapList() [1/2]

G4Voxel::MapList & G4DNAMesh::GetVoxelMapList ( const Index index)

Definition at line 52 of file G4DNAMesh.cc.

53{
54 auto key = GetKey(index);
55 return GetVoxelMapList(key);
56}
G4Voxel::MapList & GetVoxelMapList(Key key)
Definition: G4DNAMesh.cc:37

References GetKey(), and GetVoxelMapList().

◆ GetVoxelMapList() [2/2]

G4Voxel::MapList & G4DNAMesh::GetVoxelMapList ( Key  key)

Definition at line 37 of file G4DNAMesh.cc.

38{
39 auto iter = fMesh.find(key);
40 if(iter == fMesh.end())
41 {
42 G4Voxel::MapList maplist;
43 SetVoxelMapList(key, std::move(maplist));
44 return GetVoxelMapList(key);
45 }
46 else
47 {
48 return iter->second->GetMapList();
49 }
50}
void SetVoxelMapList(const Key &key, G4Voxel::MapList &&mapList)
Definition: G4DNAMesh.cc:113
std::map< MolType, size_t > MapList
Definition: G4DNAMesh.hh:77

References fMesh, GetVoxelMapList(), and SetVoxelMapList().

Referenced by G4DNAGillespieDirectMethod::ComputeNumberInNode(), G4DNAUpdateSystemModel::CreateMolecule(), GetVoxelMapList(), G4DNAUpdateSystemModel::JumpIn(), G4DNAUpdateSystemModel::JumpTo(), G4DNAUpdateSystemModel::KillMolecule(), PrintVoxel(), and G4DNAGillespieDirectMethod::PropensityFunction().

◆ PrintMesh()

void G4DNAMesh::PrintMesh ( )

Definition at line 67 of file G4DNAMesh.cc.

68{
69 G4cout << "*********PrintMesh::Size : " << fMesh.size() << G4endl;
70 auto iter = fMesh.begin();
71 for(; iter != fMesh.end(); iter++)
72 {
73 auto index = iter->second->GetIndex();
74 PrintVoxel(index);
75 }
76 G4cout << G4endl;
77}
void PrintVoxel(const Index &index)
Definition: G4DNAMesh.cc:98

References fMesh, G4cout, G4endl, and PrintVoxel().

◆ PrintVoxel()

void G4DNAMesh::PrintVoxel ( const Index index)

Definition at line 98 of file G4DNAMesh.cc.

99{
100 G4cout << "*********PrintVoxel::";
101 G4cout << "key: " << GetKey(index) << " index : " << index
102 << " number of type : " << this->GetVoxelMapList(index).size()
103 << G4endl;
104
105 for(const auto& it : this->GetVoxelMapList(index))
106 {
107 G4cout << "_____________" << it.first->GetName() << " : " << it.second
108 << G4endl;
109 }
110 G4cout << G4endl;
111}

References G4cout, G4endl, GetKey(), and GetVoxelMapList().

Referenced by G4DNAGillespieDirectMethod::Initialize(), and PrintMesh().

◆ Reset()

void G4DNAMesh::Reset ( )

Definition at line 198 of file G4DNAMesh.cc.

199{
200 if(fMesh.empty())
201 {
202 return;
203 }
204 for(auto iter : fMesh) // should use smart ptr
205 {
206 delete iter.second;
207 }
208 fMesh.clear();
209}

References fMesh.

Referenced by ~G4DNAMesh().

◆ SetVoxelMapList()

void G4DNAMesh::SetVoxelMapList ( const Key key,
G4Voxel::MapList &&  mapList 
)

Definition at line 113 of file G4DNAMesh.cc.

114{
115 auto index = GetIndex(key);
116 auto pVoxel = fMesh[key];
117 if(nullptr == pVoxel)
118 {
119 pVoxel = new G4Voxel(std::move(mapList), index, GetBoundingBox(index));
120 fMesh[key] = pVoxel;
121 }
122 else
123 {
124 assert(pVoxel->GetMapList().empty()); // check if map list is empty
125 pVoxel->SetMapList(std::move(mapList));
126 }
127}
const G4DNABoundingBox & GetBoundingBox() const
Definition: G4DNAMesh.cc:211

References fMesh, GetBoundingBox(), and GetIndex().

Referenced by GetVoxelMapList().

◆ size()

size_t G4DNAMesh::size ( )
inline

Definition at line 127 of file G4DNAMesh.hh.

127{ return fMesh.size(); }

References fMesh.

Field Documentation

◆ fMesh

VoxelMap G4DNAMesh::fMesh
private

◆ fpBoundingMesh

const G4DNABoundingBox* G4DNAMesh::fpBoundingMesh
private

◆ fResolution

G4double G4DNAMesh::fResolution
private

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