Geant4-11
G4DNAMesh.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27#ifndef G4DNAMesh_hh
28#define G4DNAMesh_hh 1
29
30#include "globals.hh"
31#include "G4DNABoundingBox.hh"
32#include <vector>
33#include <array>
34#include <cmath>
35#include <map>
36#include <memory>
37#include <sstream>
39#include "G4Track.hh"
40
42{
43 public:
44 struct Index
45 {
47 : x(0)
48 , y(0)
49 , z(0)
50 {}
51 Index(int x0, int y0, int z0)
52 : x(x0)
53 , y(y0)
54 , z(z0)
55 {}
56
57 friend std::ostream& operator<<(std::ostream& stream, const Index& rhs)
58 {
59 stream << "(" << rhs.x << ", " << rhs.y << ", " << rhs.z << ")";
60 return stream;
61 }
62 G4bool operator==(const Index& rhs) const
63 {
64 return x == rhs.x && y == rhs.y && z == rhs.z;
65 }
66 G4bool operator!=(const Index& rhs) const
67 {
68 return x != rhs.x || y != rhs.y || z != rhs.z;
69 }
70 Index operator+(const Index& rhs) const
71 {
72 return { x + rhs.x, y + rhs.y, z + rhs.z };
73 }
74 int x, y, z;
75 };
77 using MapList = std::map<MolType, size_t>;
78
79 G4Voxel(MapList&& list, Index& index, G4DNABoundingBox&& box)
80 : fMapList(std::move(list))
81 , fIndex(index)
82 , fBox(std::move(box))
83 {}
84
85 ~G4Voxel() = default;
86
87 const Index& GetIndex() const { return fIndex; }
88
90
91 void SetMapList(MapList&& mapList) { fMapList = std::move(mapList); }
92
94 {
95 auto xlo = fBox.Getxlo();
96 auto ylo = fBox.Getylo();
97 auto zlo = fBox.Getzlo();
98
99 auto xhi = fBox.Getxhi();
100 auto yhi = fBox.Getyhi();
101 auto zhi = fBox.Getzhi();
102
103 return (xhi - xlo) * (yhi - ylo) * (zhi - zlo);
104 }
105
106 private:
110};
111
113{
114 public:
116 using Key = unsigned int;
117 using VoxelMap = std::map<Key, G4Voxel*>;
119 ~G4DNAMesh();
120
121 Key GetKey(const G4ThreeVector& pos) const;
122 Index GetIndex(Key key) const;
123 Index GetIndex(const G4ThreeVector& position) const;
124
125 G4Voxel* GetVoxel(Key key);
126 [[maybe_unused]] G4Voxel* GetVoxel(const Index& index);
127 size_t size() { return fMesh.size(); }
128 Index GetIndex(const Index& index, int) const;
129 std::vector<Index> FindVoxelNeighbors(const Index& index) const;
130 std::vector<Index> FindNeighboringVoxels(const Index& index) const;
131 void Reset();
132
135
136 Key GetKey(const Index& index) const;
137 VoxelMap::iterator end() { return fMesh.end(); }
138 VoxelMap::iterator begin() { return fMesh.begin(); }
139 VoxelMap::const_iterator end() const { return fMesh.end(); }
140 VoxelMap::const_iterator begin() const { return fMesh.begin(); }
141
142 void PrintMesh();
143 void PrintVoxel(const Index& index);
144 const G4DNABoundingBox& GetBoundingBox() const;
147 void SetVoxelMapList(const Key& key, G4Voxel::MapList&& mapList);
148
149 G4double GetResolution() const;
150
151 private:
155};
156#endif
static const G4double pos
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4double Getxlo() const
G4double Getyhi() const
G4double Getylo() const
G4double Getxhi() const
G4double Getzlo() const
G4double Getzhi() const
void PrintVoxel(const Index &index)
Definition: G4DNAMesh.cc:98
void SetVoxelMapList(const Key &key, G4Voxel::MapList &&mapList)
Definition: G4DNAMesh.cc:113
~G4DNAMesh()
Definition: G4DNAMesh.cc:35
size_t size()
Definition: G4DNAMesh.hh:127
unsigned int Key
Definition: G4DNAMesh.hh:116
void PrintMesh()
Definition: G4DNAMesh.cc:67
VoxelMap::const_iterator end() const
Definition: G4DNAMesh.hh:139
VoxelMap::const_iterator begin() const
Definition: G4DNAMesh.hh:140
VoxelMap::iterator end()
Definition: G4DNAMesh.hh:137
G4double GetResolution() const
Definition: G4DNAMesh.cc:330
Key GetKey(const G4ThreeVector &pos) const
Definition: G4DNAMesh.cc:332
std::vector< Index > FindNeighboringVoxels(const Index &index) const
Definition: G4DNAMesh.cc:276
void Reset()
Definition: G4DNAMesh.cc:198
G4Voxel * GetVoxel(Key key)
Definition: G4DNAMesh.cc:216
std::vector< Index > FindVoxelNeighbors(const Index &index) const
Definition: G4DNAMesh.cc:232
G4Voxel::MapList & GetVoxelMapList(Key key)
Definition: G4DNAMesh.cc:37
const G4DNABoundingBox * fpBoundingMesh
Definition: G4DNAMesh.hh:153
std::map< Key, G4Voxel * > VoxelMap
Definition: G4DNAMesh.hh:117
VoxelMap fMesh
Definition: G4DNAMesh.hh:152
const G4DNABoundingBox & GetBoundingBox() const
Definition: G4DNAMesh.cc:211
VoxelMap::iterator begin()
Definition: G4DNAMesh.hh:138
G4double fResolution
Definition: G4DNAMesh.hh:154
G4int GetNumberOfType(G4Voxel::MolType type) const
Definition: G4DNAMesh.cc:78
Index GetIndex(Key key) const
Definition: G4DNAMesh.cc:164
G4DNAMesh(const G4DNABoundingBox &, G4int)
Definition: G4DNAMesh.cc:30
G4Voxel::Index Index
Definition: G4DNAMesh.hh:115
G4double GetVolume() const
Definition: G4DNAMesh.hh:93
G4DNABoundingBox fBox
Definition: G4DNAMesh.hh:109
MapList & GetMapList()
Definition: G4DNAMesh.hh:89
MapList fMapList
Definition: G4DNAMesh.hh:107
const Index & GetIndex() const
Definition: G4DNAMesh.hh:87
void SetMapList(MapList &&mapList)
Definition: G4DNAMesh.hh:91
std::map< MolType, size_t > MapList
Definition: G4DNAMesh.hh:77
G4Voxel(MapList &&list, Index &index, G4DNABoundingBox &&box)
Definition: G4DNAMesh.hh:79
Index fIndex
Definition: G4DNAMesh.hh:108
~G4Voxel()=default
G4bool operator==(const Index &rhs) const
Definition: G4DNAMesh.hh:62
G4bool operator!=(const Index &rhs) const
Definition: G4DNAMesh.hh:66
Index operator+(const Index &rhs) const
Definition: G4DNAMesh.hh:70
friend std::ostream & operator<<(std::ostream &stream, const Index &rhs)
Definition: G4DNAMesh.hh:57
Index(int x0, int y0, int z0)
Definition: G4DNAMesh.hh:51