Flow123d  release_3.0.0-680-gbed5aba
mapping_p1.cc
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 mapping_p1.cc
15  * @brief Class MappingP1 implements the affine transformation of
16  * the unit cell onto the actual cell.
17  * @author Jan Stebel
18  */
19 
20 #include "fem/mapping_p1.hh"
21 #include "quadrature/quadrature.hh"
22 #include "fem/fe_values.hh"
23 
24 
25 
26 using namespace std;
27 
28 
29 
30 template<unsigned int dim, unsigned int spacedim>
32 {
33 }
34 
35 template<unsigned int dim, unsigned int spacedim>
37 {
39 
40  // Initialize the gradients of the canonical basis in the
41  // barycentric coordinates.
42  // In the case of P1 mapping the shape functions are linear,
43  // hence the gradients are constant and can be precomputed.
44  if ((flags & update_jacobians) |
45  (flags & update_volume_elements) |
46  (flags & update_JxW_values) |
47  (flags & update_side_JxW_values) |
48  (flags & update_inverse_jacobians))
49  {
50  grad.zeros();
51  for (unsigned int i=0; i<dim; i++)
52  {
53  grad(0,i) = -1;
54  grad(i+1,i) = 1;
55  }
56  }
57 
58  // barycentric coordinates of quadrature points
59  if (flags & update_quadrature_points)
60  {
61  data->bar_coords.resize(q.size());
62  for (unsigned int i=0; i<q.size(); i++)
64  }
65 
66 
67 
68  return data;
69 }
70 
71 template<unsigned int dim, unsigned int spacedim>
73 {
74  UpdateFlags f = flags;
75 
76  if (flags & update_normal_vectors)
78 
79  if ((flags & update_volume_elements) |
80  (flags & update_JxW_values) |
81  (flags & update_inverse_jacobians))
82  f |= update_jacobians;
83 
84  return f;
85 }
86 
87 
88 template<unsigned int dim, unsigned int spacedim>
90  const Quadrature<dim> &q,
91  MappingInternalData &data,
93 {
94  ElementMap coords;
95  arma::mat::fixed<spacedim,dim> jac;
96 
97  if ((fv_data.update_flags & update_jacobians) |
99  (fv_data.update_flags & update_JxW_values) |
102  {
103  coords = element_map(cell);
104  }
105 
106  // calculation of Jacobian dependent data
107  if ((fv_data.update_flags & update_jacobians) |
109  (fv_data.update_flags & update_JxW_values) |
111  {
112  jac = coords*grad;
113 
114  // update Jacobians
115  if (fv_data.update_flags & update_jacobians)
116  for (unsigned int i=0; i<q.size(); i++)
117  fv_data.jacobians[i] = jac;
118 
119  // calculation of determinant dependent data
121  {
122  double det = fabs(determinant(jac));
123 
124  // update determinants
125  if (fv_data.update_flags & update_volume_elements)
126  for (unsigned int i=0; i<q.size(); i++)
127  fv_data.determinants[i] = det;
128 
129  // update JxW values
130  if (fv_data.update_flags & update_JxW_values)
131  for (unsigned int i=0; i<q.size(); i++)
132  fv_data.JxW_values[i] = det*q.weight(i);
133  }
134 
135  // update inverse Jacobians
137  {
138  arma::mat::fixed<dim,spacedim> ijac;
139  if (dim==spacedim)
140  {
141  ijac = inv(jac);
142  }
143  else
144  {
145  ijac = pinv(jac);
146  }
147  for (unsigned int i=0; i<q.size(); i++)
148  fv_data.inverse_jacobians[i] = ijac;
149  }
150  }
151 
152  // quadrature points in the actual cell coordinate system
154  {
155  BaryPoint basis;
156  for (unsigned int i=0; i<q.size(); i++)
157  fv_data.points[i] = coords*data.bar_coords[i];
158  }
159 }
160 
161 template<unsigned int dim, unsigned int spacedim>
163  unsigned int sid,
164  const Quadrature<dim> &q,
165  MappingInternalData &data,
167 {
168  ElementMap coords;
169 
170  if ((fv_data.update_flags & update_jacobians) |
173  (fv_data.update_flags & update_normal_vectors) |
175  {
176  coords = element_map(cell);
177  }
178 
179  // calculation of cell Jacobians and dependent data
180  if ((fv_data.update_flags & update_jacobians) |
184  {
185  arma::mat::fixed<spacedim,dim> jac = coords*grad;
186 
187  // update cell Jacobians
188  if (fv_data.update_flags & update_jacobians)
189  for (unsigned int i=0; i<q.size(); i++)
190  fv_data.jacobians[i] = jac;
191 
192  // update determinants of Jacobians
193  if (fv_data.update_flags & update_volume_elements)
194  {
195  double det = fabs(determinant(jac));
196  for (unsigned int i=0; i<q.size(); i++)
197  fv_data.determinants[i] = det;
198  }
199 
200  // inverse Jacobians
202  {
203  arma::mat::fixed<dim,spacedim> ijac;
204  if (dim==spacedim)
205  {
206  ijac = inv(jac);
207  }
208  else
209  {
210  ijac = pinv(jac);
211  }
212  ASSERT_LE_DBG(q.size(), fv_data.inverse_jacobians.size());
213  for (unsigned int i=0; i<q.size(); i++)
214  fv_data.inverse_jacobians[i] = ijac;
215 
216  // calculation of normal vectors to the side
217  if ((fv_data.update_flags & update_normal_vectors))
218  {
219  arma::vec::fixed<spacedim> n_cell;
220  n_cell = trans(ijac)*RefElement<dim>::normal_vector(sid);
221  n_cell = n_cell/norm(n_cell,2);
222  for (unsigned int i=0; i<q.size(); i++)
223  fv_data.normal_vectors[i] = n_cell;
224  }
225  }
226  }
227 
228  // Quadrature points in the actual cell coordinate system.
229  // The points location can vary from side to side.
231  {
232  for (unsigned int i=0; i<q.size(); i++)
233  fv_data.points[i] = coords*data.bar_coords[i];
234  }
235 
236  if (fv_data.update_flags & update_side_JxW_values)
237  {
238  double side_det;
239  if (dim <= 1)
240  {
241  side_det = 1;
242  }
243  else
244  {
245  arma::mat::fixed<spacedim,dim> side_coords;
246  arma::mat::fixed<spacedim, MatrixSizes<dim>::dim_minus_one > side_jac; // some compilers complain for case dim==0
247 
248  // calculation of side Jacobian
249  side_coords.zeros();
250  for (unsigned int n=0; n<dim; n++)
251  for (unsigned int c=0; c<spacedim; c++)
252  side_coords(c,n) = cell.side(sid)->node(n)->point()[c];
253  side_jac = side_coords * grad.submat(0,0,dim-1,dim-2);
254 
255  // calculation of JxW
256  side_det = fabs(determinant(side_jac));
257  }
258  for (unsigned int i=0; i<q.size(); i++)
259  fv_data.JxW_values[i] = side_det*q.weight(i);
260  }
261 }
262 
263 
264 template<unsigned int dim, unsigned int spacedim>
266 {
267  ElementMap coords;
268  for (unsigned int i=0; i<dim+1; i++)
269  coords.col(i) = elm.node(i)->point();
270  return coords;
271 }
272 
273 
274 template<unsigned int dim, unsigned int spacedim>
276 {
277  arma::mat::fixed<3, dim> A = map.cols(1,dim);
278  for(unsigned int i=0; i < dim; i++ ) {
279  A.col(i) -= map.col(0);
280  }
281 
282  arma::mat::fixed<dim, dim> AtA = A.t()*A;
283  arma::vec::fixed<dim> Atb = A.t()*(point - map.col(0));
284  arma::vec::fixed<dim+1> bary_coord;
285  bary_coord.rows(1, dim) = arma::solve(AtA, Atb);
286  bary_coord( 0 ) = 1.0 - arma::sum( bary_coord.rows(1,dim) );
287  return bary_coord;
288 }
289 
290 template<unsigned int dim, unsigned int spacedim>
292 {
293  return map * point;
294 }
295 
296 template<unsigned int dim, unsigned int spacedim>
298  return RefElement<dim>::clip(barycentric);
299 }
300 
301 template <unsigned int dim, unsigned int spacedim>
303 {
304  arma::vec projection = this->project_real_to_unit(point, this->element_map(elm));
305  return (projection.min() >= -BoundingBox::epsilon);
306 }
307 
308 
309 
310 template class MappingP1<1,3>;
311 template class MappingP1<2,3>;
312 template class MappingP1<3,3>;
313 
314 
315 
316 
arma::mat::fixed< spacedim, dim+1 > ElementMap
Definition: mapping_p1.hh:71
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:67
static LocalPoint normal_vector(unsigned int sid)
Transformed quadrature weight for cell sides.
RealPoint project_unit_to_real(const BaryPoint &point, const ElementMap &map) const
Definition: mapping_p1.cc:291
MappingInternalData * initialize(const Quadrature< dim > &q, UpdateFlags flags)
Initializes the structures and computes static data.
Definition: mapping_p1.cc:36
BaryPoint project_real_to_unit(const RealPoint &point, const ElementMap &map) const
Definition: mapping_p1.cc:275
Determinant of the Jacobian.
std::vector< arma::mat::fixed< dim, spacedim > > inverse_jacobians
Inverse Jacobians at the quadrature points.
Definition: fe_values.hh:122
void fill_fe_side_values(const ElementAccessor< 3 > &cell, unsigned int sid, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data)
Calculates the mapping data on a side of a cell.
Definition: mapping_p1.cc:162
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
Definition: fe_values.hh:161
UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
Definition: mapping_p1.cc:72
Class FEValues calculates finite element data on the actual cells such as shape function values...
static const double epsilon
stabilization parameter
Definition: bounding_box.hh:66
arma::vec::fixed< spacedim > RealPoint
Definition: mapping_p1.hh:70
SideIter side(const unsigned int loc_index)
Definition: accessors.hh:137
std::vector< double > JxW_values
Transformed quadrature weights.
Definition: fe_values.hh:107
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:35
Volume element.
double determinant(const T &M)
Calculates determinant of a rectangular matrix.
std::vector< arma::vec::fixed< spacedim > > normal_vectors
Shape functions (for vectorial finite elements) evaluated at quadrature points.
Definition: fe_values.hh:156
Transformed quadrature points.
std::vector< arma::mat::fixed< spacedim, dim > > jacobians
Jacobians of the mapping at the quadrature points.
Definition: fe_values.hh:112
std::vector< arma::vec::fixed< spacedim > > points
Coordinates of quadrature points in the actual cell coordinate system.
Definition: fe_values.hh:127
double weight(const unsigned int i) const
Returns the ith weight.
Definition: quadrature.hh:161
arma::vec::fixed< dim+1 > BaryPoint
Definition: mapping_p1.hh:69
std::vector< arma::vec > bar_coords
Auxiliary array of barycentric coordinates of quadrature points.
Definition: mapping.hh:105
Basic definitions of numerical quadrature rules.
MappingP1()
Constructor.
Definition: mapping_p1.cc:31
NodeAccessor< 3 > node(unsigned int i) const
Definition: side_impl.hh:47
Normal vectors.
#define ASSERT_LE_DBG(a, b)
Definition of comparative assert macro (Less or Equal) only for debug mode.
Definition: asserts.hh:307
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:86
const unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:139
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
Definition: fe_values.hh:117
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:144
static BaryPoint clip(const BaryPoint &barycentric)
Definition: ref_element.cc:431
void fill_fe_values(const ElementAccessor< 3 > &cell, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data)
Calculates the mapping data on the actual cell.
Definition: mapping_p1.cc:89
BaryPoint clip_to_element(BaryPoint &barycentric)
Definition: mapping_p1.cc:297
bool contains_point(arma::vec point, ElementAccessor< 3 > elm)
Test if element contains given point.
Definition: mapping_p1.cc:302
ElementMap element_map(ElementAccessor< 3 > elm) const
Definition: mapping_p1.cc:265
Mapping data that can be precomputed on the actual cell.
Definition: mapping.hh:99
Transformed quadrature weights.