Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 <unordered_map>
36#include <memory>
37#include <sstream>
39#include "G4VDNAMesh.hh"
40
41class G4DNAMesh : public G4VDNAMesh
42{
43 public:
46 using Data = std::map<MolType, size_t>;
47 using Voxel = std::tuple<Index, Box, Data>;
48 using IndexMap = std::unordered_map<Index, G4int, G4VDNAMesh::hashFunc>;
49 using VoxelVector = std::vector<Voxel>;
51 ~G4DNAMesh() override;
53 Voxel& GetVoxel(const Index& index); // GetorCreateVoxel
54 size_t size() { return fVoxelVector.size(); };
55 Index ConvertIndex(const Index& index, const G4int&) const;
56 std::vector<Index> FindNeighboringVoxels(const Index& index) const;
57 void Reset();
58 Data& GetVoxelMapList(const Index& index);
59 auto end() { return fVoxelVector.end(); }
60 auto begin() { return fVoxelVector.begin(); }
61 VoxelVector::const_iterator const_end() const { return fVoxelVector.end(); }
62 VoxelVector::const_iterator const_begin() const
63 {
64 return fVoxelVector.begin();
65 }
66 void PrintMesh();
67 void PrintVoxel(const Index& index);
68 const G4DNABoundingBox& GetBoundingBox() const;
69 G4DNABoundingBox GetBoundingBox(const Index& index);
70 G4int GetNumberOfType(MolType type) const;
71 void InitializeVoxel(const Index& key, Data&& mapList);
72 G4double GetResolution() const;
73
74 private:
75 IndexMap fIndexMap;
76 VoxelVector fVoxelVector;
77 const G4DNABoundingBox* fpBoundingMesh;
78 G4double fResolution;
79};
80#endif
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void InitializeVoxel(const Index &key, Data &&mapList)
Definition G4DNAMesh.cc:115
std::tuple< Index, Box, Data > Voxel
Definition G4DNAMesh.hh:47
VoxelVector::const_iterator const_end() const
Definition G4DNAMesh.hh:61
~G4DNAMesh() override
Definition G4DNAMesh.cc:58
void PrintVoxel(const Index &index)
Definition G4DNAMesh.cc:100
std::vector< Voxel > VoxelVector
Definition G4DNAMesh.hh:49
auto begin()
Definition G4DNAMesh.hh:60
G4int GetNumberOfType(MolType type) const
Definition G4DNAMesh.cc:84
Index ConvertIndex(const Index &index, const G4int &) const
Definition G4DNAMesh.cc:213
std::unordered_map< Index, G4int, G4VDNAMesh::hashFunc > IndexMap
Definition G4DNAMesh.hh:48
size_t size()
Definition G4DNAMesh.hh:54
void PrintMesh()
Definition G4DNAMesh.cc:66
Voxel & GetVoxel(const Index &index)
Definition G4DNAMesh.cc:36
G4double GetResolution() const
Definition G4DNAMesh.cc:183
VoxelVector::const_iterator const_begin() const
Definition G4DNAMesh.hh:62
const G4MolecularConfiguration * MolType
Definition G4DNAMesh.hh:45
void Reset()
Definition G4DNAMesh.cc:132
auto end()
Definition G4DNAMesh.hh:59
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