Flow123d  release_2.1.2-337-g6b7a56b
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 using namespace arma;
28 
29 
30 
31 
32 template<unsigned int dim, unsigned int spacedim>
34 {
35 }
36 
37 template<unsigned int dim, unsigned int spacedim>
39 {
41 
42  // Initialize the gradients of the canonical basis in the
43  // barycentric coordinates.
44  // In the case of P1 mapping the shape functions are linear,
45  // hence the gradients are constant and can be precomputed.
46  if ((flags & update_jacobians) |
47  (flags & update_volume_elements) |
48  (flags & update_JxW_values) |
49  (flags & update_side_JxW_values) |
50  (flags & update_inverse_jacobians))
51  {
52  grad.zeros();
53  for (unsigned int i=0; i<dim; i++)
54  {
55  grad(0,i) = -1;
56  grad(i+1,i) = 1;
57  }
58  }
59 
60  // barycentric coordinates of quadrature points
61  if (flags & update_quadrature_points)
62  {
63  vec::fixed<dim+1> basis;
64  data->bar_coords.resize(q.size());
65  for (unsigned int i=0; i<q.size(); i++)
66  {
67  basis[0] = 1;
68  for (unsigned int j=0; j<dim; j++)
69  {
70  basis[0] -= q.point(i)[j];
71  basis[j+1] = q.point(i)[j];
72  }
73  data->bar_coords[i] = basis;
74  }
75  }
76 
77 
78 
79  return data;
80 }
81 
82 template<unsigned int dim, unsigned int spacedim>
84 {
85  UpdateFlags f = flags;
86 
87  if (flags & update_normal_vectors)
89 
90  if ((flags & update_volume_elements) |
91  (flags & update_JxW_values) |
92  (flags & update_inverse_jacobians))
93  f |= update_jacobians;
94 
95  return f;
96 }
97 
98 
99 template<unsigned int dim, unsigned int spacedim>
101  const Quadrature<dim> &q,
102  MappingInternalData &data,
104 {
105  mat::fixed<spacedim,dim+1> coords;
106  mat::fixed<spacedim,dim> jac;
107 
108  if ((fv_data.update_flags & update_jacobians) |
110  (fv_data.update_flags & update_JxW_values) |
113  {
114  coords.zeros();
115  for (unsigned int n=0; n<dim+1; n++)
116  for (unsigned int c=0; c<spacedim; c++)
117  coords(c,n) = cell->node[n]->point()[c];
118  }
119 
120  // calculation of Jacobian dependent data
121  if ((fv_data.update_flags & update_jacobians) |
123  (fv_data.update_flags & update_JxW_values) |
125  {
126  jac = coords*grad;
127 
128  // update Jacobians
129  if (fv_data.update_flags & update_jacobians)
130  for (unsigned int i=0; i<q.size(); i++)
131  fv_data.jacobians[i] = jac;
132 
133  // calculation of determinant dependent data
135  {
136  double det = fabs(determinant(jac));
137 
138  // update determinants
139  if (fv_data.update_flags & update_volume_elements)
140  for (unsigned int i=0; i<q.size(); i++)
141  fv_data.determinants[i] = det;
142 
143  // update JxW values
144  if (fv_data.update_flags & update_JxW_values)
145  for (unsigned int i=0; i<q.size(); i++)
146  fv_data.JxW_values[i] = det*q.weight(i);
147  }
148 
149  // update inverse Jacobians
151  {
152  mat::fixed<dim,spacedim> ijac;
153  if (dim==spacedim)
154  {
155  ijac = inv(jac);
156  }
157  else
158  {
159  ijac = pinv(jac);
160  }
161  for (unsigned int i=0; i<q.size(); i++)
162  fv_data.inverse_jacobians[i] = ijac;
163  }
164  }
165 
166  // quadrature points in the actual cell coordinate system
168  {
169  vec::fixed<dim+1> basis;
170  for (unsigned int i=0; i<q.size(); i++)
171  fv_data.points[i] = coords*data.bar_coords[i];
172  }
173 }
174 
175 template<unsigned int dim, unsigned int spacedim>
177  unsigned int sid,
178  const Quadrature<dim> &q,
179  MappingInternalData &data,
181 {
182  mat::fixed<spacedim,dim+1> coords;
183 
184  if ((fv_data.update_flags & update_jacobians) |
187  (fv_data.update_flags & update_normal_vectors) |
189  {
190  coords.zeros();
191  for (unsigned int n=0; n<dim+1; n++)
192  for (unsigned int c=0; c<spacedim; c++)
193  coords(c,n) = cell->node[n]->point()[c];
194  }
195 
196  // calculation of cell Jacobians and dependent data
197  if ((fv_data.update_flags & update_jacobians) |
201  {
202  mat::fixed<spacedim,dim> jac = coords*grad;
203 
204  // update cell Jacobians
205  if (fv_data.update_flags & update_jacobians)
206  for (unsigned int i=0; i<q.size(); i++)
207  fv_data.jacobians[i] = jac;
208 
209  // update determinants of Jacobians
210  if (fv_data.update_flags & update_volume_elements)
211  {
212  double det = fabs(determinant(jac));
213  for (unsigned int i=0; i<q.size(); i++)
214  fv_data.determinants[i] = det;
215  }
216 
217  // inverse Jacobians
219  {
220  mat::fixed<dim,spacedim> ijac;
221  if (dim==spacedim)
222  {
223  ijac = inv(jac);
224  }
225  else
226  {
227  ijac = pinv(jac);
228  }
229  for (unsigned int i=0; i<q.size(); i++)
230  fv_data.inverse_jacobians[i] = ijac;
231 
232  // calculation of normal vectors to the side
233  if ((fv_data.update_flags & update_normal_vectors))
234  {
235  vec::fixed<spacedim> n_cell;
236  n_cell = trans(ijac)*RefElement<dim>::normal_vector(sid);
237  n_cell = n_cell/norm(n_cell,2);
238  for (unsigned int i=0; i<q.size(); i++)
239  fv_data.normal_vectors[i] = n_cell;
240  }
241  }
242  }
243 
244  // Quadrature points in the actual cell coordinate system.
245  // The points location can vary from side to side.
247  {
248  for (unsigned int i=0; i<q.size(); i++)
249  fv_data.points[i] = coords*data.bar_coords[i];
250  }
251 
252  if (fv_data.update_flags & update_side_JxW_values)
253  {
254  double side_det;
255  if (dim <= 1)
256  {
257  side_det = 1;
258  }
259  else
260  {
261  mat::fixed<spacedim,dim> side_coords;
262  mat::fixed<spacedim, MatrixSizes<dim>::dim_minus_one > side_jac; // some compilers complain for case dim==0
263 
264  // calculation of side Jacobian
265  side_coords.zeros();
266  for (unsigned int n=0; n<dim; n++)
267  for (unsigned int c=0; c<spacedim; c++)
268  side_coords(c,n) = cell->side(sid)->node(n)->point()[c];
269  side_jac = side_coords * grad.submat(0,0,dim-1,dim-2);
270 
271  // calculation of JxW
272  side_det = fabs(determinant(side_jac));
273  }
274  for (unsigned int i=0; i<q.size(); i++)
275  fv_data.JxW_values[i] = side_det*q.weight(i);
276  }
277 }
278 
279 template<>
281  unsigned int sid,
282  const Quadrature<0> &q,
283  MappingInternalData &data,
284  FEValuesData<0,3> &fv_data)
285 {}
286 
287 
288 
289 template class MappingP1<0,3>;
290 template class MappingP1<1,3>;
291 template class MappingP1<2,3>;
292 template class MappingP1<3,3>;
293 
294 
295 
296 
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.
MappingInternalData * initialize(const Quadrature< dim > &q, UpdateFlags flags)
Initializes the structures and computes static data.
Definition: mapping_p1.cc:38
Determinant of the Jacobian.
std::vector< arma::mat::fixed< dim, spacedim > > inverse_jacobians
Inverse Jacobians at the quadrature points.
Definition: fe_values.hh:82
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
Definition: fe_values.hh:121
UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
Definition: mapping_p1.cc:83
Class FEValues calculates finite element data on the actual cells such as shape function values...
std::vector< double > JxW_values
Transformed quadrature weights.
Definition: fe_values.hh:67
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:31
Volume element.
double determinant(const T &M)
Calculates determinant of a rectangular matrix.
std::vector< arma::vec::fixed< spacedim > > normal_vectors
Normal vectors to the element at the quadrature points lying on a side.
Definition: fe_values.hh:116
void fill_fe_values(const typename DOFHandlerBase::CellIterator &cell, const Quadrature< dim > &q, MappingInternalData &data, FEValuesData< dim, spacedim > &fv_data)
Calculates the mapping data on the actual cell.
Definition: mapping_p1.cc:100
Transformed quadrature points.
std::vector< arma::mat::fixed< spacedim, dim > > jacobians
Jacobians of the mapping at the quadrature points.
Definition: fe_values.hh:72
std::vector< arma::vec::fixed< spacedim > > points
Coordinates of quadrature points in the actual cell coordinate system.
Definition: fe_values.hh:87
double weight(const unsigned int i) const
Returns the ith weight.
Definition: quadrature.hh:147
void fill_fe_side_values(const typename DOFHandlerBase::CellIterator &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:176
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:33
Affine mapping between reference and actual cell.
Definition: mapping_p1.hh:53
Normal vectors.
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:46
const unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:125
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
Definition: fe_values.hh:77
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:130
Mapping data that can be precomputed on the actual cell.
Definition: mapping.hh:99
Transformed quadrature weights.