WaveBEM: Unsteady Nonlinear Potential Flow Solver for Ship-Wave Interaction.
octree_block.h
Go to the documentation of this file.
1 #ifndef octree_block_h
2 #define octree_block_h
3 
10 #include <deal.II/base/point.h>
11 
14 #include <deal.II/lac/matrix_lib.h>
15 #include <deal.II/lac/vector.h>
19 
20 #include <deal.II/grid/tria.h>
24 #include <deal.II/grid/grid_in.h>
25 #include <deal.II/grid/grid_out.h>
27 
30 #include <deal.II/dofs/dof_tools.h>
32 #include <deal.II/fe/fe_q.h>
33 #include <deal.II/fe/fe_values.h>
34 #include <deal.II/fe/fe_system.h>
36 #include <deal.II/fe/mapping_q1.h>
37 
41 
42 // And here are a few C++ standard header
43 // files that we will need:
44 #include <cmath>
45 #include <iostream>
46 #include <fstream>
47 #include <string>
48 #include <set>
49 #include <map>
50 
51 
52 using namespace dealii;
53 
54 template <int dim>
56 {
57 
58 public :
59 
60  typedef typename DoFHandler<dim-1,dim>::active_cell_iterator cell_it;
61 
62 
63 private :
64 
65  unsigned int level;
66 
67  unsigned int parentId;
68 
69  unsigned int numChildren;
70 
71  unsigned int childrenId[8];
72 
73  std::vector <std::set <unsigned int> > nearNeigh;
74 
75  std::vector <std::set <unsigned int> > intList;
76 
77  std::vector <std::set <unsigned int> > nonIntList;
78 
80 
81  double delta;
82 
83  std::vector <unsigned int> nodesId;
84 
85  std::map <cell_it, std::vector<unsigned int> > quadPointsId;
86 
87 
88 public:
89 
90  OctreeBlock();
91 
92  OctreeBlock(unsigned int level, unsigned int parent, Point<dim> pMin, double delta);
93 
94  OctreeBlock(const OctreeBlock<dim> &other);
95 
96  ~OctreeBlock();
97 
98  void CopyContent(const OctreeBlock *other);
99 
100  void AddNode(unsigned int nodeId);
101 
102  void AddQuadPoint(cell_it elemPointer, unsigned int quadPointId);
103 
104  std::vector <unsigned int> GetBlockNodeList() const;
105 
106  void DelNodeList();
107 
108  std::map <cell_it, std::vector<unsigned int> > GetBlockQuadPointsList() const;
109 
110  void DelQuadPointsList();
111 
112  unsigned int GetBlockNodesNum() const;
113 
114  unsigned int GetBlockChildrenNum() const;
115 
116  unsigned int GetParentId() const;
117 
118  void AddChild(unsigned int childId);
119 
120  unsigned int GetChildId(unsigned int idInList) const;
121 
122  Point<dim> GetPMin() const;
123 
124  double GetDelta() const;
125 
126  void AddNearNeigh(unsigned int sublevel, const unsigned int nnBlockId);
127 
128  unsigned int NumNearNeigh(unsigned int sublevel) const;
129 
130  unsigned int NumNearNeighLevels() const;
131 
132  std::set <unsigned int> GetNearNeighs(unsigned int sublevel) const;
133 
134  void AddBlockToIntList(unsigned int sublevel, const unsigned int intListBlockId);
135 
136  unsigned int NumIntList(unsigned int sublevel) const;
137 
138  unsigned int NumIntListLevels() const;
139 
140  std::set <unsigned int> GetIntList(unsigned int sublevel) const;
141 
142  std::vector<std::set <unsigned int> > GetIntList() const;
143 
144  void AddBlockToNonIntList(unsigned int sublevel, const unsigned int intListBlockId);
145 
146  unsigned int NumNonIntList(unsigned int sublevel) const;
147 
148  unsigned int NumNonIntListLevels() const;
149 
150  std::set <unsigned int> GetNonIntList(unsigned int sublevel) const;
151 
152  void SetNearNeighSize(unsigned int sublevels);
153 
154  void SetIntListSize(unsigned int sublevels);
155 
156  void SetNonIntListSize(unsigned int sublevels);
157 
158  unsigned int GetNearNeighSize() const;
159 
160  unsigned int GetIntListSize() const;
161 
162  unsigned int GetNonIntListSize() const;
163 
164 };
165 
166 
167 #endif
std::vector< std::set< unsigned int > > intList
Definition: octree_block.h:75
std::vector< std::set< unsigned int > > nonIntList
Definition: octree_block.h:77
std::vector< unsigned int > nodesId
Definition: octree_block.h:83
double delta
Definition: octree_block.h:81
std::vector< std::set< unsigned int > > nearNeigh
Definition: octree_block.h:73
unsigned int level
Definition: octree_block.h:65
unsigned int numChildren
Definition: octree_block.h:69
DoFHandler< dim-1, dim >::active_cell_iterator cell_it
Definition: octree_block.h:60
unsigned int parentId
Definition: octree_block.h:67
Point< dim > pMin
Definition: octree_block.h:79
std::map< cell_it, std::vector< unsigned int > > quadPointsId
Definition: octree_block.h:85