DOLFIN
DOLFIN C++ interface
ParallelRefinement.h
1 // Copyright (C) 2012-2014 Chris Richardson
2 //
3 // This file is part of DOLFIN.
4 //
5 // DOLFIN is free software: you can redistribute it and/or modify
6 // it under the terms of the GNU Lesser General Public License as published by
7 // the Free Software Foundation, either version 3 of the License, or
8 // (at your option) any later version.
9 //
10 // DOLFIN is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 // GNU Lesser General Public License for more details.
14 //
15 // You should have received a copy of the GNU Lesser General Public License
16 // along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
17 //
18 //
19 // First Added: 2013-01-02
20 
21 #ifndef __PARALLEL_REFINEMENT_H
22 #define __PARALLEL_REFINEMENT_H
23 
24 #include <unordered_map>
25 #include <vector>
26 
27 namespace dolfin
28 {
29 
30  // Forward declarations
31  class Mesh;
32  template<typename T> class MeshFunction;
33 
35 
41 
43  {
44  public:
45 
48 
51 
53  const Mesh& mesh() const
54  { return _mesh; }
55 
58  bool is_marked(std::size_t edge_index) const;
59 
63  void mark(std::size_t edge_index);
64 
66  void mark_all();
67 
71  void mark(const MeshFunction<bool>& refinement_marker);
72 
75  void mark(const MeshEntity& cell);
76 
80  std::vector<std::size_t> marked_edge_list(const MeshEntity& cell) const;
81 
84 
88  void create_new_vertices();
89 
92  std::shared_ptr<const std::map<std::size_t, std::size_t> > edge_to_new_vertex() const;
93 
96  void new_cell(const Cell& cell);
97 
103  void new_cell(std::size_t i0, std::size_t i1, std::size_t i2,
104  std::size_t i3);
105 
110  void new_cell(std::size_t i0, std::size_t i1, std::size_t i2);
111 
114  void new_cells(const std::vector<std::size_t>& idx);
115 
119  void partition(Mesh& new_mesh, bool redistribute) const;
120 
123  void build_local(Mesh& new_mesh) const;
124 
125  private:
126 
127  // Mesh reference
128  const Mesh& _mesh;
129 
130  // Shared edges between processes. In R^2, vector size is 1
131  std::unordered_map<unsigned int, std::vector<std::pair<unsigned int,
132  unsigned int> > > shared_edges;
133 
134  // Mapping from old local edge index to new global vertex, needed
135  // to create new topology
136  std::shared_ptr<std::map<std::size_t, std::size_t> > local_edge_to_new_vertex;
137 
138  // New storage for all coordinates when creating new vertices
139  std::vector<double> new_vertex_coordinates;
140 
141  // New storage for all cells when creating new topology
142  std::vector<std::size_t> new_cell_topology;
143 
144  // Management of marked edges
145  std::vector<bool> marked_edges;
146 
147  // Temporary storage for edges that have been recently marked
148  std::vector<std::vector<std::size_t> > marked_for_update;
149  };
150 
151 }
152 
153 #endif
std::vector< std::size_t > marked_edge_list(const MeshEntity &cell) const
Definition: ParallelRefinement.cpp:116
std::shared_ptr< const std::map< std::size_t, std::size_t > > edge_to_new_vertex() const
Definition: ParallelRefinement.cpp:90
Definition: adapt.h:41
~ParallelRefinement()
Destructor.
Definition: ParallelRefinement.cpp:53
Definition: adapt.h:29
Data structure and methods for refining meshes in parallel.
Definition: ParallelRefinement.h:42
void update_logical_edgefunction()
Transfer marked edges between processes.
Definition: ParallelRefinement.cpp:130
void partition(Mesh &new_mesh, bool redistribute) const
Definition: ParallelRefinement.cpp:290
const Mesh & mesh() const
Original mesh associated with this refinement.
Definition: ParallelRefinement.h:53
void mark(std::size_t edge_index)
Definition: ParallelRefinement.cpp:64
A Cell is a MeshEntity of topological codimension 0.
Definition: Cell.h:42
ParallelRefinement(const Mesh &mesh)
Constructor.
Definition: ParallelRefinement.cpp:44
void new_cells(const std::vector< std::size_t > &idx)
Definition: ParallelRefinement.cpp:363
bool is_marked(std::size_t edge_index) const
Definition: ParallelRefinement.cpp:58
void create_new_vertices()
Definition: ParallelRefinement.cpp:148
void build_local(Mesh &new_mesh) const
Definition: ParallelRefinement.cpp:245
Definition: MeshEntity.h:42
void mark_all()
Mark all edges in mesh.
Definition: ParallelRefinement.cpp:84
void new_cell(const Cell &cell)
Definition: ParallelRefinement.cpp:340
Definition: Mesh.h:82