Flow123d  release_3.0.0-688-g6b683cf
dofhandler.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file dofhandler.hh
15  * @brief Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between local and global dofs.
16  * @author Jan Stebel
17  */
18 
19 #ifndef DOFHANDLER_HH_
20 #define DOFHANDLER_HH_
21 
22 #include <vector> // for vector
23 #include <unordered_map> // for unordered_map
24 #include "mesh/side_impl.hh"
25 #include "mesh/mesh.h"
26 #include "mesh/accessors.hh"
27 #include "mesh/long_idx.hh" // for LongIdx
28 #include "mesh/range_wrapper.hh"
30 #include "fem/discrete_space.hh" // for DiscreteSpace
31 #include "petscvec.h" // for Vec
32 
33 template<unsigned int dim> class FiniteElement;
34 class DHCellAccessor;
35 class Mesh;
36 class Distribution;
37 class Dof;
38 
39 
40 /**
41  * Class DOFHandlerBase provides an abstract interface for various dof handlers:
42  * - basic handler for a given spatial dimension
43  * - multi-dimensional handler for all spatial dimensions (1D-3D)
44  * - handler for a specific region / part of mesh
45  */
47 public:
48 
49  /**
50  * @brief Constructor.
51  * @param _mesh The mesh.
52  */
54  : n_global_dofs_(0), lsize_(0), loffset_(0), max_elem_dofs_(0), mesh_(&_mesh), dof_ds_(0) {}
55 
56  /**
57  * @brief Getter for the number of all mesh dofs required by the given
58  * finite element.
59  */
60  const unsigned int n_global_dofs() const { return n_global_dofs_; }
61 
62  /**
63  * @brief Returns the number of dofs on the current process.
64  */
65  const unsigned int lsize() const { return lsize_; }
66 
67  /**
68  * @brief Returns the offset of the local part of dofs.
69  */
70  const unsigned int loffset() const { return loffset_; }
71 
72  /**
73  * @brief Returns max. number of dofs on one element.
74  */
75  const unsigned int max_elem_dofs() const { return max_elem_dofs_; }
76 
77  std::shared_ptr<Distribution> distr() const { return dof_ds_; }
78 
79  /**
80  * @brief Returns the mesh.
81  */
82  Mesh *mesh() const { return mesh_; }
83 
84  /**
85  * @brief Fill vector of the global indices of dofs associated to the @p cell.
86  *
87  * @param cell The cell.
88  * @param indices Vector of dof indices on the cell.
89  */
90  virtual unsigned int get_dof_indices(const ElementAccessor<3> &cell, std::vector<LongIdx> &indices) const = 0;
91 
92  /**
93  * @brief Fill vector of the indices of dofs associated to the @p cell on the local process.
94  *
95  * @param cell The cell.
96  * @param indices Vector of dof indices on the cell.
97  */
98  virtual unsigned int get_loc_dof_indices(const ElementAccessor<3> &cell, std::vector<LongIdx> &indices) const =0;
99 
100  /**
101  * @brief Compute hash value of DOF handler.
102  */
103  virtual std::size_t hash() const =0;
104 
105  /// Destructor.
106  virtual ~DOFHandlerBase();
107 
108 protected:
109 
110  /**
111  * @brief Number of global dofs assigned by the handler.
112  */
113  unsigned int n_global_dofs_;
114 
115  /**
116  * @brief Number of dofs associated to local process.
117  */
118  unsigned int lsize_;
119 
120  /**
121  * @brief Index of the first dof on the local process.
122  */
123  unsigned int loffset_;
124 
125  /// Max. number of dofs per element.
126  unsigned int max_elem_dofs_;
127 
128  /**
129  * @brief Pointer to the mesh to which the dof handler is associated.
130  */
132 
133  /**
134  * @brief Distribution of dofs associated to local process.
135  */
136  std::shared_ptr<Distribution> dof_ds_;
137 
138 };
139 
140 
141 
142 
143 /**
144  * @brief Provides the numbering of the finite element degrees of freedom
145  * on the computational mesh.
146  *
147  * Class DOFHandlerMultiDim distributes the degrees of freedom (dof) for
148  * a particular triplet of 1d, 2d and 3d finite elements on the computational mesh
149  * and provides mappings between local and global dofs.
150  * The template parameter @p dim denotes the spatial dimension of
151  * the reference finite element.
152  *
153  * Currently the functionality is restricted to finite elements with internal and nodal dofs,
154  * i.e. the neighboring elements can share only dofs on nodes.
155  */
157 public:
158 
159  /**
160  * @brief Constructor.
161  * @param _mesh The mesh.
162  * @param make_elem_part Allow switch off make_element_partitioning, necessary for boundary DOF handler.
163  */
164  DOFHandlerMultiDim(Mesh &_mesh, bool make_elem_part = true);
165 
166 
167  /**
168  * @brief Distributes degrees of freedom on the mesh needed
169  * for the given discrete space.
170  *
171  * By default, the dof handler is parallel, meaning that each
172  * processor has access to dofs on the local elements and on one
173  * layer of ghost elements (owned by neighbouring elements).
174  * This can be changed by setting @p sequential to true.
175  *
176  * @param ds The discrete space consisting of finite elements for each mesh element.
177  * @param sequential If true then each processor will have information about all dofs.
178  */
179  void distribute_dofs(std::shared_ptr<DiscreteSpace> ds);
180 
181  /** @brief Returns sequential version of the current dof handler.
182  *
183  * Collective on all processors.
184  */
185  std::shared_ptr<DOFHandlerMultiDim> sequential();
186 
187  /**
188  * @brief Returns scatter context from parallel to sequential vectors.
189  *
190  * For sequential dof handler it returns null pointer.
191  * Collective on all processors.
192  */
193  std::shared_ptr<VecScatter> sequential_scatter();
194 
195  /**
196  * @brief Returns the global indices of dofs associated to the @p cell.
197  *
198  * @param cell The cell.
199  * @param indices Array of dof indices on the cell.
200  */
201  unsigned int get_dof_indices(const ElementAccessor<3> &cell,
202  std::vector<LongIdx> &indices) const override;
203 
204  /**
205  * @brief Returns the indices of dofs associated to the @p cell on the local process.
206  *
207  * @param cell The cell.
208  * @param indices Array of dof indices on the cell.
209  */
210  unsigned int get_loc_dof_indices(const ElementAccessor<3> &cell,
211  std::vector<LongIdx> &indices) const override;
212 
213  /**
214  * @brief Returns the global index of local element.
215  *
216  * @param loc_el Local index of element.
217  */
218  inline int el_index(int loc_el) const { return el_4_loc[loc_el]; }
219 
220  /**
221  * @brief Returns the global index of local edge.
222  *
223  * @param loc_edg Local index of edge.
224  */
225  inline LongIdx edge_index(int loc_edg) const { return edg_4_loc[loc_edg]; }
226 
227  /**
228  * @brief Returns the global index of local neighbour.
229  *
230  * @param loc_nb Local index of neighbour.
231  */
232  inline LongIdx nb_index(int loc_nb) const { return nb_4_loc[loc_nb]; }
233 
234  /**
235  * @brief Return number of dofs on given cell.
236  *
237  * @param cell Cell accessor.
238  */
239  unsigned int n_dofs(ElementAccessor<3> cell) const;
240 
241  /**
242  * @brief Returns number of local edges.
243  */
244  inline unsigned int n_loc_edges() const { return edg_4_loc.size(); }
245 
246  /**
247  * @brief Returns number of local neighbours.
248  */
249  inline unsigned int n_loc_nb() const { return nb_4_loc.size(); }
250 
251  /**
252  * Returns true if element is on local process.
253  * @param index Global element index.
254  */
255  bool el_is_local(int index) const;
256 
257  /**
258  * @brief Returns finite element object for given space dimension.
259  *
260  * @param cell Cell accessor.
261  */
262  template<unsigned int dim>
263  FiniteElement<dim> *fe(const ElementAccessor<3> &cell) const { return ds_->fe<dim>(cell); }
264 
265  /**
266  * @brief Return dof on a given cell.
267  * @param cell Mesh cell.
268  * @param idof Number of dof on the cell.
269  */
270  const Dof &cell_dof(ElementAccessor<3> cell,
271  unsigned int idof) const;
272 
273  /// Output structure of dof handler.
274  void print() const;
275 
276  /**
277  * Implements @p DOFHandlerBase::hash.
278  */
279  std::size_t hash() const override;
280 
281  /// Returns range of DOF handler cells (only range of own without ghost cells)
282  Range<DHCellAccessor> own_range() const;
283 
284  /// Returns range over own and ghost cells of DOF handler
285  Range<DHCellAccessor> local_range() const;
286 
287  /// Returns range over ghosts DOF handler cells
288  Range<DHCellAccessor> ghost_range() const;
289 
290  /// Return DHCellAccessor appropriate to ElementAccessor of given idx
291  DHCellAccessor cell_accessor_from_element(unsigned int elm_idx) const;
292 
293  /// Destructor.
294  ~DOFHandlerMultiDim() override;
295 
296 
297 
298  friend class DHCellAccessor;
299 
300 private:
301 
302  /**
303  * @brief Prepare parallel distribution of elements, edges and neighbours.
304  */
305  void make_elem_partitioning();
306 
307  /**
308  * @brief Initialize vector of starting indices for elements.
309  */
310  void init_cell_starts();
311 
312  /**
313  * @brief Initialize auxiliary vector of starting indices of nodal dofs.
314  *
315  * @param node_dof_starts Vector of starting indices (output).
316  */
317  void init_node_dof_starts(std::vector<LongIdx> &node_dof_starts);
318 
319  /**
320  * @brief Initialize node_status.
321  *
322  * Set VALID_NODE for nodes owned by local elements and
323  * INVALID_NODE for nodes owned by ghost elements.
324  *
325  * @param node_status Vector of nodal status (output).
326  */
327  void init_node_status(std::vector<short int> &node_status);
328 
329  /**
330  * @brief Obtain dof numbers on ghost elements from other processor.
331  * @param proc Neighbouring processor.
332  * @param dofs Array where dofs are stored (output).
333  */
334  void receive_ghost_dofs(unsigned int proc,
335  std::vector<LongIdx> &dofs);
336 
337  /**
338  * @brief Send dof numbers to other processor.
339  * @param proc Neighbouring processor.
340  */
341  void send_ghost_dofs(unsigned int proc);
342 
343  /**
344  * @brief Update dofs on local elements from ghost element dofs.
345  *
346  * @param proc Neighbouring processor.
347  * @param update_cells Vector of global indices of elements which need to be updated
348  * from ghost elements.
349  * @param dofs Vector of dof indices on ghost elements from processor @p proc.
350  * @param node_dof_starts Vector of starting indices of nodal dofs.
351  * @param node_dofs Vector of nodal dof indices (output).
352  */
353  void update_local_dofs(unsigned int proc,
354  const std::vector<bool> &update_cells,
355  const std::vector<LongIdx> &dofs,
356  const std::vector<LongIdx> &node_dof_starts,
357  std::vector<LongIdx> &node_dofs
358  );
359 
360  /**
361  * @brief Communicate local dof indices to all processors and create new sequential dof handler.
362  *
363  * Collective on all processors.
364  */
365  void create_sequential();
366 
367 
368  /**
369  * Flags used during distribution of dofs to mark node and dof status.
370  */
371  static const int INVALID_NODE = 1;
372  static const int VALID_NODE = 2;
373  static const int ASSIGNED_NODE = 3;
374  static const int INVALID_DOF = -1;
375 
376 
377  /// Pointer to the discrete space for which the handler distributes dofs.
378  std::shared_ptr<DiscreteSpace> ds_;
379 
380  /// Indicator for parallel/sequential dof handler.
382 
383  /// Sequential dof handler associated to the current (parallel) one.
384  std::shared_ptr<DOFHandlerMultiDim> dh_seq_;
385 
386  /// Scatter context for parallel to sequential vectors.
387  std::shared_ptr<VecScatter> scatter_to_seq_;
388 
389  /**
390  * @brief Starting indices for element dofs.
391  *
392  * E.g. dof_indices[cell_starts[idx]] = dof number for first dof on the
393  * cell with index idx within the parallel structure. To use with element
394  * accessor use the following:
395  *
396  * ElementAccessor<3> cell;
397  * ...
398  * // i-th dof number on the cell
399  * dof_indices[cell_starts[row_4_el[cell.idx()]]+i] = ...
400  *
401  * For parallel dof handler, only local and ghost elements are stored,
402  * but the vector has size mesh_->n_elements()+1.
403  */
405 
406  /**
407  * @brief Dof numbers on local and ghost elements.
408  *
409  * Dofs are ordered accordingly with cell_starts and local dof order
410  * given by the finite element. See cell_starts for more description.
411  */
413 
414  /**
415  * @brief Maps local and ghost dof indices to global ones.
416  *
417  * First lsize_ entries correspond to dofs owned by local processor,
418  * the remaining entries are ghost dofs sorted by neighbouring processor id.
419  */
421 
422  /**
423  * @brief Maps global element index into local/ghost index (obsolete).
424  */
425  std::unordered_map<LongIdx,LongIdx> global_to_local_el_idx_;
426 
427  /// Global element index -> index according to partitioning
429 
430  /// Local element index -> global element index
432 
433  /// Distribution of elements
435 
436  /// Local edge index -> global edge index
438 
439  /// Local neighbour index -> global neighbour index
441 
442  /// Indices of local nodes in mesh tree.
444 
445  /// Indices of ghost cells (neighbouring with local elements).
447 
448  /// Processors of ghost elements.
449  set<unsigned int> ghost_proc;
450 
451  /// Arrays of ghost cells for each neighbouring processor.
453 
454 };
455 
456 
457 
458 
459 #endif /* DOFHANDLER_HH_ */
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: long_idx.hh:22
virtual ~DOFHandlerBase()
Destructor.
Definition: dofhandler.cc:34
LongIdx * el_4_loc
Local element index -> global element index.
Definition: dofhandler.hh:431
const unsigned int lsize() const
Returns the number of dofs on the current process.
Definition: dofhandler.hh:65
FMT_API void print(std::FILE *f, CStringRef format_str, ArgList args)
Definition: format.cc:489
Range helper class.
Definition: mesh.h:52
Distribution * el_ds_
Distribution of elements.
Definition: dofhandler.hh:434
Definition: mesh.h:80
Template Iter serves as general template for internal iterators.
std::vector< LongIdx > local_to_global_dof_idx_
Maps local and ghost dof indices to global ones.
Definition: dofhandler.hh:420
unsigned int n_loc_nb() const
Returns number of local neighbours.
Definition: dofhandler.hh:249
std::shared_ptr< Distribution > dof_ds_
Distribution of dofs associated to local process.
Definition: dofhandler.hh:136
int el_index(int loc_el) const
Returns the global index of local element.
Definition: dofhandler.hh:218
Mesh * mesh() const
Returns the mesh.
Definition: dofhandler.hh:82
std::shared_ptr< DOFHandlerMultiDim > dh_seq_
Sequential dof handler associated to the current (parallel) one.
Definition: dofhandler.hh:384
unsigned int lsize_
Number of dofs associated to local process.
Definition: dofhandler.hh:118
map< unsigned int, vector< LongIdx > > ghost_proc_el
Arrays of ghost cells for each neighbouring processor.
Definition: dofhandler.hh:452
virtual unsigned int get_dof_indices(const ElementAccessor< 3 > &cell, std::vector< LongIdx > &indices) const =0
Fill vector of the global indices of dofs associated to the cell.
std::vector< LongIdx > cell_starts
Starting indices for element dofs.
Definition: dofhandler.hh:404
const unsigned int loffset() const
Returns the offset of the local part of dofs.
Definition: dofhandler.hh:70
const unsigned int max_elem_dofs() const
Returns max. number of dofs on one element.
Definition: dofhandler.hh:75
Mesh * mesh_
Pointer to the mesh to which the dof handler is associated.
Definition: dofhandler.hh:131
std::vector< LongIdx > dof_indices
Dof numbers on local and ghost elements.
Definition: dofhandler.hh:412
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:156
vector< LongIdx > nb_4_loc
Local neighbour index -> global neighbour index.
Definition: dofhandler.hh:440
virtual unsigned int get_loc_dof_indices(const ElementAccessor< 3 > &cell, std::vector< LongIdx > &indices) const =0
Fill vector of the indices of dofs associated to the cell on the local process.
unsigned int max_elem_dofs_
Max. number of dofs per element.
Definition: dofhandler.hh:126
vector< LongIdx > ghost_4_loc
Indices of ghost cells (neighbouring with local elements).
Definition: dofhandler.hh:446
std::unordered_map< LongIdx, LongIdx > global_to_local_el_idx_
Maps global element index into local/ghost index (obsolete).
Definition: dofhandler.hh:425
Declaration of class which provides the finite element for every mesh cell.
unsigned int n_loc_edges() const
Returns number of local edges.
Definition: dofhandler.hh:244
virtual std::size_t hash() const =0
Compute hash value of DOF handler.
std::shared_ptr< Distribution > distr() const
Definition: dofhandler.hh:77
LongIdx edge_index(int loc_edg) const
Returns the global index of local edge.
Definition: dofhandler.hh:225
vector< LongIdx > edg_4_loc
Local edge index -> global edge index.
Definition: dofhandler.hh:437
unsigned int loffset_
Index of the first dof on the local process.
Definition: dofhandler.hh:123
set< unsigned int > ghost_proc
Processors of ghost elements.
Definition: dofhandler.hh:449
vector< LongIdx > node_4_loc
Indices of local nodes in mesh tree.
Definition: dofhandler.hh:443
std::shared_ptr< DiscreteSpace > ds_
Pointer to the discrete space for which the handler distributes dofs.
Definition: dofhandler.hh:378
unsigned int n_global_dofs_
Number of global dofs assigned by the handler.
Definition: dofhandler.hh:113
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
bool is_parallel_
Indicator for parallel/sequential dof handler.
Definition: dofhandler.hh:381
FiniteElement< dim > * fe(const ElementAccessor< 3 > &cell) const
Returns finite element object for given space dimension.
Definition: dofhandler.hh:263
DOFHandlerBase(Mesh &_mesh)
Constructor.
Definition: dofhandler.hh:53
std::shared_ptr< VecScatter > scatter_to_seq_
Scatter context for parallel to sequential vectors.
Definition: dofhandler.hh:387
LongIdx * row_4_el
Global element index -> index according to partitioning.
Definition: dofhandler.hh:428
LongIdx nb_index(int loc_nb) const
Returns the global index of local neighbour.
Definition: dofhandler.hh:232
Implementation of range helper class.
const unsigned int n_global_dofs() const
Getter for the number of all mesh dofs required by the given finite element.
Definition: dofhandler.hh:60