Flow123d  JS_before_hm-1-ge159291
accessors.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 accessors.hh
15  * @brief
16  */
17 
18 #ifndef ACCESSORS_HH_
19 #define ACCESSORS_HH_
20 
21 #include "mesh/bounding_box.hh"
22 #include "mesh/region.hh"
23 #include "mesh/elements.h"
24 #include "mesh/mesh.h"
25 #include "mesh/nodes.hh"
26 #include "mesh/node_accessor.hh"
27 #include "la/distribution.hh"
28 #include <vector>
29 #include <armadillo>
30 
31 
32 /**
33  * Due to circular dependence of return parameters in mesh accessors,
34  * and the intention to have all methods inlined,
35  * we gathered the accessors into a single file.
36  * This way, it is possible to implement it, see the simple test below.
37  *
38  * Do not include accessors_impl.hh file anywhere but at the end of this file!
39  *
40  * previous include loops:
41  * side -> boundary -> elemen accessor
42  * element accessor <-> side
43  * edge <-> side
44  * boundary -> edge
45  */
46 
47 // Compilable test of class loop dependence of return parameters.
48 // class A;
49 // class B;
50 
51 // class A{
52 // public:
53 // A()
54 // {};
55 
56 // B create_b();
57 // };
58 
59 // class B{
60 // public:
61 // B()
62 // {};
63 
64 // A create_a();
65 // };
66 
67 
68 // B A::create_b()
69 // { return B(); }
70 
71 // A B::create_a()
72 // { return A(); }
73 
74 
75 class Side;
76 class SiderIter;
77 class Edge;
78 class Boundary;
79 
80 /**
81  * Element accessor templated just by dimension of the embedding space, used by Fields.
82  * This should allow algorithms over elements where dimension of particular element is runtime parameter.
83  *
84  * This class suites as interface of Fields to the mesh elements, in particular this accessor knows directly
85  * the region, and also can be used as an accessor that works on the whole region if used by Fields that do not depend on
86  * particular elements as FieldConstant, FiledFormula, and FieldPython.
87  *
88  * TODO:
89  * - make this kind of accessor subclass of FieldCommonBase or at least move it into src/fields
90  * since it has functionality particular for Fields
91  *
92  * Ideas:
93  * need function to calculate intersection (object) of two ElementAccessors, but this definitely should be templated by
94  * dimension of the ref. element (or rather shape of ref. element), here we can have case dispatch
95  *
96  */
97 template <int spacedim>
98 class ElementAccessor {
99 public:
100  /// Default invalid accessor.
101  ElementAccessor();
102 
103  /// Regional accessor.
104  ElementAccessor(const Mesh *mesh, RegionIdx r_idx);
105 
106  /// Element accessor.
107  ElementAccessor(const Mesh *mesh, unsigned int idx);
108 
109  /// Incremental function of the Element iterator.
110  void inc();
111 
112  /// Return list of element vertices.
114 
115  /// Computes the measure of the element.
116  double measure() const;
117 
118  /** Computes the Jacobian of the element.
119  * J = det ( 1 1 1 1 )
120  * x1 x2 x3 x4
121  * y1 y2 y3 y4
122  * z1 z2 z3 z4
123  */
124  double tetrahedron_jacobian() const;
125 
126  /// Computes the barycenter.
127  arma::vec::fixed<spacedim> centre() const;
128 
129  /**
130  * Quality of the element based on the smooth and scale-invariant quality measures proposed in:
131  * J. R. Schewchuk: What is a Good Linear Element?
132  *
133  * We scale the measure so that is gives value 1 for regular elements. Line 1d elements
134  * have always quality 1.
135  */
136  double quality_measure_smooth(SideIter side) const;
137 
138  SideIter side(const unsigned int loc_index);
139 
140  const SideIter side(const unsigned int loc_index) const;
141 
142 
143 
144  bool is_regional() const {
145  return dim_ == undefined_dim_;
146  }
147 
148  bool is_elemental() const {
149  return ( is_valid() && ! is_regional() );
150  }
151 
152  bool is_valid() const {
153  return mesh_ != NULL;
154  }
155 
156  unsigned int dim() const
157  { return dim_; }
158 
159  const Element * element() const {
160  return &(mesh_->element_vec_[element_idx_]);
161  }
162 
163 
164  Region region() const
165  { return Region( r_idx_, mesh_->region_db()); }
166 
168  { return r_idx_; }
169 
170  /// We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular mesh
171  //unsigned int region_id() const {
172  // return region().id();
173  //}
174 
175  bool is_boundary() const {
176  return boundary_;
177  }
178 
179  /// Return local idx of element in boundary / bulk part of element vector
180  unsigned int idx() const {
181  if (boundary_) return ( element_idx_ - mesh_->bulk_size_ );
182  else return element_idx_;
183  }
184 
185  /// Return global idx of element in full element vector
186  unsigned int mesh_idx() const {
187  return element_idx_;
188  }
189 
190  unsigned int index() const {
191  return (unsigned int)mesh_->find_elem_id(element_idx_);
192  }
193 
194  unsigned int proc() const {
196  }
197 
198  const Node * node(unsigned int ni) const {
199  return &(mesh_->node_vec_[element()->node_idx(ni)]);
200  }
201 
202  NodeAccessor<3> node_accessor(unsigned int ni) const {
203  return mesh_->node_accessor( element()->node_idx(ni) );
204  }
205 
206  /**
207  * Return bounding box of the element.
208  * Simpler code, but need to check performance penelty.
209  */
211  return BoundingBox(this->vertex_list());
212  }
213 
215  return (element_idx_ == other.element_idx_);
216  }
217 
218  /**
219  * -> dereference operator
220  *
221  * Allow simplify calling of element() method. Example:
222  @code
223  ElementAccessor<3> elm_ac(mesh, index);
224  arma::vec centre;
225  centre = elm_ac.element()->node_idx(0); // full format of access to element
226  centre = elm_ac->node_idx(0); // short format with dereference operator
227  @endcode
228  */
229  const Element * operator ->() const {
230  return &(mesh_->element_vec_[element_idx_]);
231  }
232 
233 
234 private:
235  /**
236  * When dim_ == undefined_dim_ ; the value of element_idx_ is invalid.
237  * Is used for ElementAccessors for whole region
238  */
239  static const unsigned int undefined_dim_ = 100;
240 
241  /// Dimension of reference element.
242  unsigned int dim_;
243 
244  /// Pointer to the mesh owning the element.
245  const Mesh *mesh_;
246  /// True if the element is boundary
247  bool boundary_;
248 
249  /// Index into Mesh::element_vec_ array.
250  unsigned int element_idx_;
251 
252  /// Region index.
254 };
255 
256 
257 
258 
259 
260 //=============================================================================
261 // Edge class
262 //=============================================================================
263 class Edge
264 {
265 public:
266  /// Default invalid edge accessor constructor.
267  Edge();
268 
269  /// Valid edge accessor constructor.
270  Edge(const Mesh *mesh, unsigned int edge_idx);
271 
272  /// Gets side iterator of the @p i -th side.
273  SideIter side(const unsigned int i) const;
274 
275 
276 
277  bool is_valid() const
278  { return mesh_ != nullptr; }
279 
280  /// Returns edge global index.
281  unsigned int idx() const {
282  ASSERT_DBG(is_valid());
283  return edge_idx_;
284  }
285 
286  /// Incremental function of the Edge iterator.
287  void inc() {
288  ASSERT_DBG(is_valid()).error("Do not call inc() for invalid accessor!");
289  edge_idx_++;
290  }
291 
292  /// Comparison operator of the iterator.
293  bool operator==(const Edge& other) const{
294  return (edge_idx_ == other.edge_idx_);
295  }
296 
297  /// Returns number of sides aligned with the edge.
298  unsigned int n_sides() const
299  { return edge_data()->n_sides;}
300 
301 private:
302  /// Pointer to the mesh owning the node.
303  const Mesh *mesh_;
304  /// Index into Mesh::edges vector.
305  unsigned int edge_idx_;
306 
307  /// Getter for edge data from mesh.
308  const EdgeData* edge_data() const;
309 };
310 
311 
312 
313 
314 
315 //=============================================================================
316 // Boundary class
317 //=============================================================================
318 class Boundary
319 {
320 public:
321  Boundary();
322  Boundary(BoundaryData* boundary_data);
323 
324  Edge edge();
325  ElementAccessor<3> element_accessor();
326  Region region();
327  Element * element();
328 
329  bool is_valid() const {
330  return boundary_data_ != nullptr;
331  }
332 
333  Mesh* mesh() {
334  ASSERT_DBG(is_valid());
335  return boundary_data_->mesh_;
336  }
337 
339  ASSERT_DBG(is_valid());
340  return boundary_data_->edge_idx_;
341  }
342 
344  ASSERT_DBG(is_valid());
345  return boundary_data_->bc_ele_idx_;
346  }
347 
348 private:
350 };
351 
352 
353 
354 
355 
356 //=============================================================================
357 // Side class
358 //=============================================================================
359 class Side {
360 public:
361  /// Default invalid side accessor constructor.
362  Side();
363 
364  /// Valid edge accessor constructor.
365  Side(const Mesh * mesh, unsigned int elem_idx, unsigned int set_lnum);
366 
367  double measure() const; ///< Calculate metrics of the side
368  arma::vec3 centre() const; ///< Centre of side
369  arma::vec3 normal() const; ///< Vector of (generalized) normal
370  double diameter() const; ///< Calculate the side diameter.
371 
372  /// Returns dimension of the side, that is dimension of the element minus one.
373  unsigned int dim() const;
374 
375  /// Returns true for all sides either on boundary or connected to vb neigboring.
376  bool is_external() const;
377 
378  /// Returns true for side on the boundary.
379  bool is_boundary() const;
380 
381  /// Returns node for given local index @p i on the side.
382  NodeAccessor<3> node(unsigned int i) const;
383 
384  /// Returns iterator to the element of the side.
385  ElementAccessor<3> element() const;
386 
387  /// Returns global index of the edge connected to the side.
388  unsigned int edge_idx() const;
389 
390  /// Returns pointer to the edge connected to the side.
391  Edge edge() const;
392 
393  /**
394  * Returns boundary condition prescribed on the side.
395  * Fails on assert if side if not on boundary and no BC is prescribed.
396  */
397  Boundary cond() const;
398 
399  /// Returns global index of the prescribed boundary condition.
400  unsigned int cond_idx() const;
401 
402 
403 
404  /// Returns number of nodes of the side.
405  unsigned int n_nodes() const
406  { return dim()+1; }
407 
408  /// Returns pointer to the mesh.
409  const Mesh * mesh() const
410  { return this->mesh_; }
411 
412  /// Returns local index of the side on the element.
413  unsigned int side_idx() const
414  { return side_idx_; }
415 
416  /// Returns index of element in Mesh::element_vec_.
417  unsigned int elem_idx() const
418  { return elem_idx_; }
419 
420  /// Returns true if the side has assigned element.
421  bool is_valid() const
422  { return mesh_!= nullptr; }
423 
424  /// Iterate over local sides of the element.
425  void inc() {
426  ASSERT_DBG(is_valid()).error("Do not call inc() for invalid accessor!");
427  side_idx_++;
428  }
429 
430  /// This is necessary by current DofHandler, should change this
431  //void *make_ptr() const;
432 private:
433 
434  arma::vec3 normal_point() const;
435  arma::vec3 normal_line() const;
436  arma::vec3 normal_triangle() const;
437 
438  // Topology of the mesh
439 
440  const Mesh * mesh_; ///< Pointer to Mesh to which belonged
441  unsigned int elem_idx_; ///< Index of element in Mesh::element_vec_
442  unsigned int side_idx_; ///< Local # of side in element (to remove it, we heve to remove calc_side_rhs)
443 
444 };
445 
446 
447 /*
448  * Iterator to a side.
449  */
450 class SideIter {
451 public:
453  {}
454 
456  : side_(side)
457  {}
458 
459  bool operator==(const SideIter &other) {
460  return (side_.mesh() == other.side_.mesh() ) && ( side_.elem_idx() == other.side_.elem_idx() )
461  && ( side_.side_idx() == other.side_.side_idx() );
462  }
463 
464 
465  bool operator!=(const SideIter &other) {
466  return !( *this == other);
467  }
468 
469  /// * dereference operator
470  const Side & operator *() const
471  { return side_; }
472 
473  /// -> dereference operator
474  const Side * operator ->() const
475  { return &side_; }
476 
477  /// prefix increment iterate only on local element
478  SideIter &operator ++ () {
479  side_.inc();
480  return (*this);
481  }
482 
483 private:
485 };
486 
487 
488 
489 #include "mesh/accessors_impl.hh"
490 
491 #endif /* ACCESSORS_HH_ */
unsigned int get_proc(unsigned int idx) const
get processor of the given index
unsigned int n_sides() const
Returns number of sides aligned with the edge.
Definition: accessors.hh:298
const Element * element() const
Definition: accessors.hh:159
Bounding box in 3d ambient space.
Definition: bounding_box.hh:54
vector< Element > element_vec_
Definition: mesh.h:512
void inc()
Incremental function of the Edge iterator.
Definition: accessors.hh:287
LongIdx * get_row_4_el() const
Definition: mesh.h:165
unsigned int uint
NodeAccessor< 3 > node_accessor(unsigned int ni) const
Definition: accessors.hh:202
unsigned int dim_
Dimension of reference element.
Definition: accessors.hh:242
Definition: nodes.hh:31
Nodes of a mesh.
unsigned int side_idx() const
Returns local index of the side on the element.
Definition: accessors.hh:413
unsigned int idx() const
Returns edge global index.
Definition: accessors.hh:281
uint edge_idx()
Definition: accessors.hh:338
bool is_boundary() const
We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular...
Definition: accessors.hh:175
unsigned int mesh_idx() const
Return global idx of element in full element vector.
Definition: accessors.hh:186
NodeAccessor< 3 > node_accessor(unsigned int idx) const
Create and return NodeAccessor to node of given idx.
Definition: mesh.cc:745
const Mesh * mesh_
Pointer to the mesh owning the node.
Definition: accessors.hh:303
Definition: mesh.h:76
bool is_valid() const
Definition: accessors.hh:152
const Mesh * mesh_
Pointer to the mesh owning the element.
Definition: accessors.hh:245
arma::Mat< Type >::template fixed< resRows, resCols > operator*(const Mat< Type, resRows, commonDimension > &a, const Mat< Type, commonDimension, resCols > &b)
Definition: armor.hh:191
uint bc_ele_idx()
Definition: accessors.hh:343
SideIter side(const unsigned int loc_index)
BoundaryData * boundary_data_
Definition: accessors.hh:349
const RegionDB & region_db() const
Definition: mesh.h:141
unsigned int element_idx_
Index into Mesh::element_vec_ array.
Definition: accessors.hh:250
unsigned int proc() const
Definition: accessors.hh:194
Side side_
Definition: accessors.hh:484
Mesh * mesh()
Definition: accessors.hh:333
bool operator==(const ElementAccessor< spacedim > &other)
Definition: accessors.hh:214
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
unsigned int elem_idx_
Index of element in Mesh::element_vec_.
Definition: accessors.hh:441
bool operator==(const Edge &other) const
Comparison operator of the iterator.
Definition: accessors.hh:293
static const unsigned int undefined_dim_
Definition: accessors.hh:239
double quality_measure_smooth(SideIter side) const
bool is_elemental() const
Definition: accessors.hh:148
const Mesh * mesh() const
Returns pointer to the mesh.
Definition: accessors.hh:409
int find_elem_id(unsigned int pos) const
Return element id (in GMSH file) of element of given position in element vector.
Definition: mesh.h:359
unsigned int node_idx(unsigned int ni) const
Return index (in Mesh::node_vec) of ni-th node.
Definition: elements.h:74
bool is_valid() const
Definition: accessors.hh:329
unsigned int index() const
Definition: accessors.hh:190
bool is_valid() const
Returns true if the side has assigned element.
Definition: accessors.hh:421
bool operator!=(const SideIter &other)
Definition: accessors.hh:465
void inc()
Incremental function of the Element iterator.
BoundingBox bounding_box() const
Definition: accessors.hh:210
vector< Node > node_vec_
Definition: mesh.h:529
double measure() const
Computes the measure of the element.
bool boundary_
True if the element is boundary.
Definition: accessors.hh:247
Distribution * get_el_ds() const
Definition: mesh.h:162
ElementAccessor()
Default invalid accessor.
Region region() const
Definition: accessors.hh:164
const Element * operator->() const
Definition: accessors.hh:229
unsigned int side_idx_
Local # of side in element (to remove it, we heve to remove calc_side_rhs)
Definition: accessors.hh:442
RegionIdx r_idx_
Region index.
Definition: accessors.hh:253
Support classes for parallel programing.
unsigned int edge_idx_
Index into Mesh::edges vector.
Definition: accessors.hh:305
vector< arma::vec3 > vertex_list() const
Return list of element vertices.
const Mesh * mesh_
Pointer to Mesh to which belonged.
Definition: accessors.hh:440
RegionIdx region_idx() const
Definition: accessors.hh:167
#define ASSERT_DBG(expr)
double tetrahedron_jacobian() const
unsigned int elem_idx() const
Returns index of element in Mesh::element_vec_.
Definition: accessors.hh:417
void inc()
Iterate over local sides of the element.
Definition: accessors.hh:425
bool operator==(const SideIter &other)
Definition: accessors.hh:459
unsigned int dim() const
Definition: accessors.hh:156
bool is_valid() const
Definition: accessors.hh:277
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:180
unsigned int bulk_size_
Count of bulk elements.
Definition: mesh.h:518
unsigned int n_nodes() const
Returns number of nodes of the side.
Definition: accessors.hh:405
SideIter(const Side &side)
Definition: accessors.hh:455
const Node * node(unsigned int ni) const
Definition: accessors.hh:198
bool is_regional() const
Definition: accessors.hh:144