Flow123d  release_2.2.0-914-gf1a3a4f
elements.h
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 elements.h
15  * @brief
16 */
17 
18 #ifndef ELEMENTS_H
19 #define ELEMENTS_H
20 
21 #include "mesh/nodes.hh"
22 #include "mesh/region.hh"
23 #include "mesh/bounding_box.hh"
24 #include "mesh/ref_element.hh"
25 
26 template <int spacedim>
27 class ElementAccessor;
28 
29 class Mesh;
30 class Side;
31 class SideIter;
32 class Neighbour;
33 
34 
35 
36 //=============================================================================
37 // STRUCTURE OF THE ELEMENT OF THE MESH
38 //=============================================================================
39 class Element
40 {
41 public:
42  Element();
43  Element(unsigned int dim, Mesh *mesh_in, RegionIdx reg);
44  void init(unsigned int dim, Mesh *mesh_in, RegionIdx reg);
45  ~Element();
46 
47 
48  inline unsigned int dim() const;
49  inline unsigned int index() const;
50  unsigned int n_sides() const; // Number of sides
51  unsigned int n_nodes() const; // Number of nodes
52 
53  ///Gets ElementAccessor of this element
55 
56  /// Computes the measure of the element.
57  double measure() const;
58 
59  /** Computes the Jacobian of the element.
60  * J = det ( 1 1 1 1 )
61  * x1 x2 x3 x4
62  * y1 y2 y3 y4
63  * z1 z2 z3 z4
64  */
65  double tetrahedron_jacobian() const;
66 
67  /// Computes the barycenter.
68  arma::vec3 centre() const;
69  /**
70 * Quality of the element based on the smooth and scale-invariant quality measures proposed in:
71 * J. R. Schewchuk: What is a Good Linear Element?
72 *
73 * We scale the measure so that is gives value 1 for regular elements. Line 1d elements
74 * have always quality 1.
75 */
76  double quality_measure_smooth();
77 
78  unsigned int n_sides_by_dim(unsigned int side_dim);
79  inline SideIter side(const unsigned int loc_index);
80  inline const SideIter side(const unsigned int loc_index) const;
81  Region region() const;
82  inline RegionIdx region_idx() const
83  { return region_idx_; }
84 
85  unsigned int id() const;
86 
87  int pid; // Id # of mesh partition
88 
89  // Type specific data
90  Node** node; // Element's nodes
91 
92 
93  unsigned int *edge_idx_; // Edges on sides
94  unsigned int *boundary_idx_; // Possible boundaries on sides (REMOVE) all bcd assembly should be done through iterating over boundaries
95  // ?? deal.ii has this not only boundary iterators
96  /**
97  * Indices of permutations of nodes on sides.
98  * It determines, in which order to take the nodes of the side so as to obtain
99  * the same order as on the reference side (side 0 on the particular edge).
100  *
101  * Permutations are defined in RefElement::side_permutations.
102  */
103  unsigned int *permutation_idx_;
104 
105  /**
106  * Computes bounding box of element (OBSOLETE) ??
107  */
109 
110  /// Return precomputed bounding box.
112 
113  /**
114  * Return bounding box of the element.
115  * Simpler code, but need to check performance penelty.
116  */
118  return BoundingBox(this->vertex_list());
119  }
120 
121  /**
122  * Return list of element vertices.
123  */
125  vector<arma::vec3> vertices(this->n_nodes());
126  for(unsigned int i=0; i<n_nodes(); i++) vertices[i]=node[i]->point();
127  return vertices;
128  }
129 
130  unsigned int get_proc() const;
131 
132 
133  unsigned int n_neighs_vb; // # of neighbours, V-B type (comp.)
134  // only ngh from this element to higher dimension edge
135  Neighbour **neigh_vb; // List og neighbours, V-B type (comp.)
136 
137 
138  Mesh *mesh_; // should be removed as soon as the element is also an Accessor
139 
140 
141 protected:
142  // Data readed from mesh file
144  unsigned int dim_;
145 
146  friend class Mesh;
147 
148  template<int spacedim, class Value>
149  friend class Field;
150 
151 };
152 
153 
154 
155 
156 #define FOR_ELEMENT_NODES(i,j) for((j)=0;(j)<(i)->n_nodes();(j)++)
157 #define FOR_ELEMENT_SIDES(i,j) for(unsigned int j=0; j < (i)->n_sides(); j++)
158 #define FOR_ELM_NEIGHS_VB(i,j) for((j)=0;(j)<(i)->n_neighs_vb;(j)++)
159 
160 
161 #endif
162 //-----------------------------------------------------------------------------
163 // vim: set cindent:
Definition: sides.h:31
Bounding box in 3d ambient space.
Definition: bounding_box.hh:45
double measure() const
Computes the measure of the element.
Definition: elements.cc:90
unsigned int n_nodes() const
unsigned int * boundary_idx_
Definition: elements.h:94
Definition: nodes.hh:32
Nodes of a mesh.
int pid
Definition: elements.h:87
unsigned int * permutation_idx_
Definition: elements.h:103
RegionIdx region_idx() const
Definition: elements.h:82
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:62
unsigned int get_proc() const
Definition: elements.cc:221
unsigned int index() const
Definition: mesh.h:99
double tetrahedron_jacobian() const
Definition: elements.cc:116
Node ** node
Definition: elements.h:90
unsigned int id() const
Definition: elements.cc:170
Element()
Definition: elements.cc:35
unsigned int dim() const
BoundingBox & get_bounding_box_fast(BoundingBox &bounding_box) const
Return precomputed bounding box.
Definition: elements.cc:214
unsigned int * edge_idx_
Definition: elements.h:93
~Element()
Definition: elements.cc:81
Neighbour ** neigh_vb
Definition: elements.h:135
Region region() const
Definition: elements.cc:165
unsigned int n_sides() const
vector< arma::vec3 > vertex_list()
Definition: elements.h:124
void get_bounding_box(BoundingBox &bounding_box) const
Definition: elements.cc:206
SideIter side(const unsigned int loc_index)
unsigned int dim_
Definition: elements.h:144
void init(unsigned int dim, Mesh *mesh_in, RegionIdx reg)
Definition: elements.cc:61
Mesh * mesh_
Definition: elements.h:138
RegionIdx region_idx_
Definition: elements.h:143
arma::vec3 centre() const
Computes the barycenter.
Definition: elements.cc:130
unsigned int n_neighs_vb
Definition: elements.h:133
double quality_measure_smooth()
Definition: elements.cc:174
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
ElementAccessor< 3 > element_accessor() const
Gets ElementAccessor of this element.
Definition: elements.cc:158
unsigned int n_sides_by_dim(unsigned int side_dim)
Definition: elements.cc:147
BoundingBox bounding_box()
Definition: elements.h:117