Flow123d  jenkins-Flow123d-windows32-release-multijob-51
mapping_p1.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief Class MappingP1 implements the affine transformation of
27  * the unit cell onto the actual cell.
28  * @author Jan Stebel
29  */
30 
31 
32 #include "fem/mapping_p1.hh"
33 #include "quadrature/quadrature.hh"
34 #include "fem/fe_values.hh"
35 
36 
37 
38 using namespace std;
39 using namespace arma;
40 
41 
42 
43 
44 template<unsigned int dim, unsigned int spacedim>
46 {
47 }
48 
49 template<unsigned int dim, unsigned int spacedim>
51 {
53 
54  // Initialize the gradients of the canonical basis in the
55  // barycentric coordinates.
56  // In the case of P1 mapping the shape functions are linear,
57  // hence the gradients are constant and can be precomputed.
58  if ((flags & update_jacobians) |
59  (flags & update_volume_elements) |
60  (flags & update_JxW_values) |
61  (flags & update_side_JxW_values) |
62  (flags & update_inverse_jacobians))
63  {
64  grad.zeros();
65  for (unsigned int i=0; i<dim; i++)
66  {
67  grad(0,i) = -1;
68  grad(i+1,i) = 1;
69  }
70  }
71 
72  // barycentric coordinates of quadrature points
73  if (flags & update_quadrature_points)
74  {
75  vec::fixed<dim+1> basis;
76  data->bar_coords.resize(q.size());
77  for (unsigned int i=0; i<q.size(); i++)
78  {
79  basis[0] = 1;
80  for (unsigned int j=0; j<dim; j++)
81  {
82  basis[0] -= q.point(i)[j];
83  basis[j+1] = q.point(i)[j];
84  }
85  data->bar_coords[i] = basis;
86  }
87  }
88 
89 
90 
91  return data;
92 }
93 
94 template<unsigned int dim, unsigned int spacedim>
96 {
97  UpdateFlags f = flags;
98 
99  if (flags & update_normal_vectors)
101 
102  if ((flags & update_volume_elements) |
103  (flags & update_JxW_values) |
104  (flags & update_inverse_jacobians))
105  f |= update_jacobians;
106 
107  return f;
108 }
109 
110 
111 template<unsigned int dim, unsigned int spacedim>
113  const Quadrature<dim> &q,
114  MappingInternalData &data,
116 {
117  mat::fixed<spacedim,dim+1> coords;
118  mat::fixed<spacedim,dim> jac;
119 
120  if ((fv_data.update_flags & update_jacobians) |
122  (fv_data.update_flags & update_JxW_values) |
125  {
126  coords.zeros();
127  for (unsigned int n=0; n<dim+1; n++)
128  for (unsigned int c=0; c<spacedim; c++)
129  coords(c,n) = cell->node[n]->point()[c];
130  }
131 
132  // calculation of Jacobian dependent data
133  if ((fv_data.update_flags & update_jacobians) |
135  (fv_data.update_flags & update_JxW_values) |
137  {
138  jac = coords*grad;
139 
140  // update Jacobians
141  if (fv_data.update_flags & update_jacobians)
142  for (unsigned int i=0; i<q.size(); i++)
143  fv_data.jacobians[i] = jac;
144 
145  // calculation of determinant dependent data
147  {
148  double det = fabs(determinant(jac));
149 
150  // update determinants
151  if (fv_data.update_flags & update_volume_elements)
152  for (unsigned int i=0; i<q.size(); i++)
153  fv_data.determinants[i] = det;
154 
155  // update JxW values
156  if (fv_data.update_flags & update_JxW_values)
157  for (unsigned int i=0; i<q.size(); i++)
158  fv_data.JxW_values[i] = det*q.weight(i);
159  }
160 
161  // update inverse Jacobians
163  {
164  mat::fixed<dim,spacedim> ijac;
165  if (dim==spacedim)
166  {
167  ijac = inv(jac);
168  }
169  else
170  {
171  ijac = pinv(jac);
172  }
173  for (unsigned int i=0; i<q.size(); i++)
174  fv_data.inverse_jacobians[i] = ijac;
175  }
176  }
177 
178  // quadrature points in the actual cell coordinate system
180  {
181  vec::fixed<dim+1> basis;
182  for (unsigned int i=0; i<q.size(); i++)
183  fv_data.points[i] = coords*data.bar_coords[i];
184  }
185 }
186 
187 template<unsigned int dim, unsigned int spacedim>
189  unsigned int sid,
190  const Quadrature<dim> &q,
191  MappingInternalData &data,
193 {
194  mat::fixed<spacedim,dim+1> coords;
195 
196  if ((fv_data.update_flags & update_jacobians) |
199  (fv_data.update_flags & update_normal_vectors) |
201  {
202  coords.zeros();
203  for (unsigned int n=0; n<dim+1; n++)
204  for (unsigned int c=0; c<spacedim; c++)
205  coords(c,n) = cell->node[n]->point()[c];
206  }
207 
208  // calculation of cell Jacobians and dependent data
209  if ((fv_data.update_flags & update_jacobians) |
213  {
214  mat::fixed<spacedim,dim> jac = coords*grad;
215 
216  // update cell Jacobians
217  if (fv_data.update_flags & update_jacobians)
218  for (unsigned int i=0; i<q.size(); i++)
219  fv_data.jacobians[i] = jac;
220 
221  // update determinants of Jacobians
222  if (fv_data.update_flags & update_volume_elements)
223  {
224  double det = fabs(determinant(jac));
225  for (unsigned int i=0; i<q.size(); i++)
226  fv_data.determinants[i] = det;
227  }
228 
229  // inverse Jacobians
231  {
232  mat::fixed<dim,spacedim> ijac;
233  if (dim==spacedim)
234  {
235  ijac = inv(jac);
236  }
237  else
238  {
239  ijac = pinv(jac);
240  }
241  for (unsigned int i=0; i<q.size(); i++)
242  fv_data.inverse_jacobians[i] = ijac;
243 
244  // calculation of normal vectors to the side
245  if ((fv_data.update_flags & update_normal_vectors))
246  {
247  vec::fixed<spacedim> n_cell;
248  n_cell = trans(ijac)*RefElement<dim>::normal_vector(sid);
249  n_cell = n_cell/norm(n_cell,2);
250  for (unsigned int i=0; i<q.size(); i++)
251  fv_data.normal_vectors[i] = n_cell;
252  }
253  }
254  }
255 
256  // Quadrature points in the actual cell coordinate system.
257  // The points location can vary from side to side.
259  {
260  for (unsigned int i=0; i<q.size(); i++)
261  fv_data.points[i] = coords*data.bar_coords[i];
262  }
263 
264  if (fv_data.update_flags & update_side_JxW_values)
265  {
266  double side_det;
267  if (dim <= 1)
268  {
269  side_det = 1;
270  }
271  else
272  {
273  mat::fixed<spacedim,dim> side_coords;
274  mat::fixed<spacedim, MatrixSizes<dim>::dim_minus_one > side_jac; // some compilers complain for case dim==0
275 
276  // calculation of side Jacobian
277  side_coords.zeros();
278  for (unsigned int n=0; n<dim; n++)
279  for (unsigned int c=0; c<spacedim; c++)
280  side_coords(c,n) = cell->side(sid)->node(n)->point()[c];
281  side_jac = side_coords * grad.submat(0,0,dim-1,dim-2);
282 
283  // calculation of JxW
284  side_det = fabs(determinant(side_jac));
285  }
286  for (unsigned int i=0; i<q.size(); i++)
287  fv_data.JxW_values[i] = side_det*q.weight(i);
288  }
289 }
290 
291 template<>
293  unsigned int sid,
294  const Quadrature<0> &q,
295  MappingInternalData &data,
296  FEValuesData<0,3> &fv_data)
297 {}
298 
299 
300 
301 template class MappingP1<0,3>;
302 template class MappingP1<1,3>;
303 template class MappingP1<2,3>;
304 template class MappingP1<3,3>;
305 
306 
307 
308 
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:78
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:50
Determinant of the Jacobian.
std::vector< arma::mat::fixed< dim, spacedim > > inverse_jacobians
Inverse Jacobians at the quadrature points.
Definition: fe_values.hh:93
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
Definition: fe_values.hh:132
UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
Definition: mapping_p1.cc:95
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:78
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:42
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:127
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:112
Transformed quadrature points.
std::vector< arma::mat::fixed< spacedim, dim > > jacobians
Jacobians of the mapping at the quadrature points.
Definition: fe_values.hh:83
std::vector< arma::vec::fixed< spacedim > > points
Coordinates of quadrature points in the actual cell coordinate system.
Definition: fe_values.hh:98
double weight(const unsigned int i) const
Returns the ith weight.
Definition: quadrature.hh:158
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:188
std::vector< arma::vec > bar_coords
Auxiliary array of barycentric coordinates of quadrature points.
Definition: mapping.hh:116
Basic definitions of numerical quadrature rules.
MappingP1()
Constructor.
Definition: mapping_p1.cc:45
Affine mapping between reference and actual cell.
Definition: mapping_p1.hh:64
Normal vectors.
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:57
const unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:136
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
Definition: fe_values.hh:88
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:141
Mapping data that can be precomputed on the actual cell.
Definition: mapping.hh:110
static arma::vec::fixed< dim > normal_vector(unsigned int sid)
Definition: ref_element.cc:109
Transformed quadrature weights.