Flow123d  jenkins-Flow123d-linux-release-multijob-198
fe_rt.hh
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 Definitions of Raviart-Thomas finite elements.
27  * @author Jan Stebel
28  */
29 
30 #ifndef FE_RT_HH_
31 #define FE_RT_HH_
32 
33 #include "fem/finite_element.hh"
34 
35 
36 
37 /**
38  * @brief Raviart-Thomas element of order 0.
39  *
40  * The lowest order Raviart-Thomas finite element with linear basis functions
41  * and continuous normal components across element sides.
42  */
43 template <unsigned int dim, unsigned int spacedim>
44 class FE_RT0 : public FiniteElement<dim,spacedim>
45 {
55 
56 public:
57 
58  /// Number of raw basis functions.
59  static const unsigned int n_raw_functions = dim+1;
60 
61  /**
62  * @brief Constructor.
63  */
64  FE_RT0();
65 
66  /**
67  * @brief The scalar variant of basis_vector must be implemented but may not be used.
68  */
69  double basis_value(const unsigned int i, const arma::vec::fixed<dim> &p) const;
70 
71  /**
72  * @brief The scalar variant of basis_grad_vector must be implemented but may not be used.
73  */
74  arma::vec::fixed<dim> basis_grad(const unsigned int i, const arma::vec::fixed<dim> &p) const;
75 
76  /**
77  * @brief Returns the @p ith basis function evaluated at the point @p p.
78  * @param i Number of the basis function.
79  * @param p Point of evaluation.
80  */
81  arma::vec::fixed<dim> basis_vector(const unsigned int i, const arma::vec::fixed<dim> &p) const;
82 
83  /**
84  * @brief Returns the gradient of the @p ith basis function at the point @p p.
85  * @param i Number of the basis function.
86  * @param p Point of evaluation.
87  */
88  arma::mat::fixed<dim,dim> basis_grad_vector(const unsigned int i, const arma::vec::fixed<dim> &p) const;
89 
90  /**
91  * @brief Computes the conversion matrix from internal basis to shape functions.
92  *
93  * Initializes the @p node_matrix for computing the coefficients
94  * of the raw basis functions from values at support points.
95  * Since Raviart-Thomas element is not Lagrangean, the method has
96  * to be reimplemented.
97  */
98  void compute_node_matrix();
99 
100  /**
101  * @brief Calculates the data on the reference cell.
102  *
103  * @param q Quadrature.
104  * @param flags Flags that indicate what quantities should be calculated.
105  */
107 
108  /**
109  * @brief Decides which additional quantities have to be computed
110  * for each cell.
111  */
113 
114  /**
115  * @brief Computes the shape function values and gradients on the actual cell
116  * and fills the FEValues structure.
117  *
118  * @param q Quadrature.
119  * @param data The precomputed finite element data on the reference cell.
120  * @param fv_data The data to be computed.
121  */
122  virtual void fill_fe_values(const Quadrature<dim> &q,
123  FEInternalData &data,
124  FEValuesData<dim,spacedim> &fv_data);
125 
126  /**
127  * @brief Destructor.
128  */
129  virtual ~FE_RT0() {};
130 
131 };
132 
133 
134 
135 
136 
137 
138 
139 
140 
141 
142 template<unsigned int dim, unsigned int spacedim>
144 {
145  arma::vec::fixed<dim> sp;
146 
147  this->init();
148 
149  number_of_dofs = dim+1;
150  number_of_single_dofs[dim] = dim+1;
151 
152  sp.fill(1./max(1.,(double)dim));
153  generalized_support_points.push_back(sp);
154  for (unsigned int i=0; i<dim; i++)
155  {
156  sp.fill(1./max(1.,(double)dim));
157  sp[i] = 0;
158  generalized_support_points.push_back(sp);
159  }
160 
161  order = 1;
162 
163  is_scalar_fe = false;
164 
165  compute_node_matrix();
166 }
167 
168 template<unsigned int dim, unsigned int spacedim>
169 double FE_RT0<dim,spacedim>::basis_value(const unsigned int i, const arma::vec::fixed<dim> &p) const
170 {
171  ASSERT(false, "basis_value() may not be called for vectorial finite element.");
172 
173  return 0.0;
174 }
175 
176 template<unsigned int dim, unsigned int spacedim>
177 arma::vec::fixed<dim> FE_RT0<dim,spacedim>::basis_grad(const unsigned int i, const arma::vec::fixed<dim> &p) const
178 {
179  ASSERT(false, "basis_grad() may not be called for vectorial finite element.");
180  return arma::vec::fixed<dim>();
181 }
182 
183 template<unsigned int dim, unsigned int spacedim>
184 arma::vec::fixed<dim> FE_RT0<dim,spacedim>::basis_vector(const unsigned int i, const arma::vec::fixed<dim> &p) const
185 {
186  ASSERT(i<n_raw_functions, "Index of basis function is out of range.");
187 
188  arma::vec::fixed<dim> v(p);
189 
190  if (i > 0)
191  v[i-1] -= 1;
192 
193  return v;
194 }
195 
196 template<unsigned int dim, unsigned int spacedim>
197 arma::mat::fixed<dim,dim> FE_RT0<dim,spacedim>::basis_grad_vector(const unsigned int i, const arma::vec::fixed<dim> &p) const
198 {
199  ASSERT(i<n_raw_functions, "Index of basis function is out of range.");
200 
201  return arma::eye(dim,dim);
202 }
203 
204 template<unsigned int dim, unsigned int spacedim>
206 {
207  arma::mat::fixed<n_raw_functions,dim+1> F;
208  arma::vec::fixed<dim> r;
209 
210  /*
211  * Node matrix helps creating the shape functions $\{b_k\}$ from
212  * the raw basis $\{r_i\}$:
213  *
214  * $$ b_k = \sum_{i=0}^{dim*(dim+1)-1} N_{ki} r_i. $$
215  *
216  * The shape functions must obey the flux condition
217  *
218  * $$ b_k\cdot n_j |\Gamma_j| = \delta_{kj}, $$
219  *
220  * where $n_j$, $|\Gamma_j|$ is the unit outward normal vector and
221  * the area of the $j$-th side, respectively. Consequently,
222  * the node matrix $N$ is determined as the Moon-Penrose
223  * pseudoinverse of the flux matrix, i.e.:
224  *
225  * $$ NF = I,\quad F_{ij} = r_i\cdot n_j |\Gamma_j|. $$
226  *
227  */
228 
229  for (unsigned int i=0; i<n_raw_functions; i++)
230  {
231  /*
232  * For the 0-th side we have:
233  *
234  * $$ |\Gamma_0| = \frac{\sqrt{dim}}{dim-1}, $$
235  * $$ n_0 = \frac{\sqrt{dim}}{dim}(1,\ldots,1). $$
236  *
237  * If dim==1, we set |\Gamma_0| = 1.
238  *
239  */
240 
241  r = basis_vector(i,generalized_support_points[0]);
242  if (dim == 1)
243  {
244  F(i,0) = sum(r);
245  }
246  else
247  {
248  F(i,0) = sum(r)/(1.*dim-1);
249  }
250 
251  /*
252  * For the remaining sides:
253  *
254  * $$ |\Gamma_j| = \sqrt(dim-1), $$
255  * $$ (n_j)_l = -\delta_{jl}. $$
256  *
257  * If dim==1, we set |\Gamma_1| = 1.
258  *
259  */
260 
261  if (dim == 1)
262  {
263  r = basis_vector(i,generalized_support_points[1]);
264  F(i,1) = -r(0);
265  }
266  else
267  {
268  for (unsigned int j=1; j<dim+1; j++)
269  {
270  r = basis_vector(i,generalized_support_points[j]);
271  F(i,j) = -r(j-1)/(1.*dim-1);
272  }
273  }
274  }
275 
276  if (dim>0) node_matrix = inv(F);
277 
278 }
279 
280 template<unsigned int dim, unsigned int spacedim>
282 {
283  FEInternalData *data = new FEInternalData;
284 
285  if (flags & update_values)
286  {
287  arma::mat::fixed<n_raw_functions,dim> raw_values;
288  arma::mat::fixed<dim+1,dim> shape_values;
289  vector<arma::vec> values;
290 
291  data->basis_vectors.resize(q.size());
292  values.resize(dim+1);
293  for (unsigned int i=0; i<q.size(); i++)
294  {
295  for (unsigned int j=0; j<n_raw_functions; j++)
296  raw_values.row(j) = trans(basis_vector(j, q.point(i)));
297 
298  shape_values = node_matrix * raw_values;
299 
300  for (unsigned int j=0; j<dim+1; j++)
301  values[j] = trans(shape_values.row(j));
302 
303  data->basis_vectors[i] = values;
304  }
305  }
306 
307  if (flags & update_gradients)
308  {
309  arma::mat::fixed<dim,dim> grad;
310  arma::mat::fixed<dim,dim> shape_grads;
311  vector<arma::mat> grads;
312 
313  data->basis_grad_vectors.resize(q.size());
314  grads.resize(dim+1);
315  for (unsigned int i=0; i<q.size(); i++)
316  {
317  for (unsigned int k=0; k<dim+1; k++)
318  {
319  grad.zeros();
320  for (unsigned int l=0; l<n_raw_functions; l++)
321  grad += basis_grad_vector(l, q.point(i)) * node_matrix(k,l);
322  grads[k] = grad;
323  }
324 
325  data->basis_grad_vectors[i] = grads;
326  }
327  }
328 
329  return data;
330 }
331 
332 template<unsigned int dim, unsigned int spacedim> inline
334 {
335  UpdateFlags f = flags;
336 
337  if (flags & update_values)
339 
340  if (flags & update_gradients)
342 
343  return f;
344 }
345 
346 template<unsigned int dim, unsigned int spacedim> inline
348  const Quadrature<dim> &q,
349  FEInternalData &data,
351 {
352  // shape values
353  if (fv_data.update_flags & update_values)
354  {
355  vector<arma::vec::fixed<spacedim> > vectors;
356  vectors.resize(dim+1);
357  for (unsigned int i = 0; i < q.size(); i++)
358  {
359  for (unsigned int k=0; k<dim+1; k++)
360  vectors[k] = fv_data.jacobians[i]*data.basis_vectors[i][k]/fv_data.determinants[i];
361 
362  fv_data.shape_vectors[i] = vectors;
363  }
364  }
365 
366  // shape gradients
367  if (fv_data.update_flags & update_gradients)
368  {
369  vector<arma::mat::fixed<spacedim,spacedim> > grads;
370  grads.resize(dim+1);
371  for (unsigned int i = 0; i < q.size(); i++)
372  {
373  for (unsigned int k=0; k<dim+1; k++)
374  grads[k] = fv_data.jacobians[i]*data.basis_grad_vectors[i][k]*fv_data.inverse_jacobians[i]/fv_data.determinants[i];
375 
376  fv_data.shape_grad_vectors[i] = grads;
377  }
378  }
379 }
380 
381 
382 
383 
384 
385 
386 
387 
388 
389 
390 
391 
392 
393 #endif /* FE_RT_HH_ */
Shape function values.
Definition: update_flags.hh:98
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Definition: update_flags.hh:78
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
double basis_value(const unsigned int i, const arma::vec::fixed< dim > &p) const
The scalar variant of basis_vector must be implemented but may not be used.
Definition: fe_rt.hh:169
FE_RT0()
Constructor.
Definition: fe_rt.hh:143
arma::mat::fixed< dim, dim > basis_grad_vector(const unsigned int i, const arma::vec::fixed< dim > &p) const
Returns the gradient of the ith basis function at the point p.
Definition: fe_rt.hh:197
std::vector< std::vector< arma::mat::fixed< spacedim, spacedim > > > shape_grad_vectors
Gradients of shape functions (for vectorial finite elements).
Definition: fe_values.hh:121
FEInternalData * initialize(const Quadrature< dim > &q, UpdateFlags flags)
Calculates the data on the reference cell.
Definition: fe_rt.hh:281
std::vector< std::vector< arma::vec > > basis_vectors
Precomputed values of basis functions at the quadrature points.
virtual void fill_fe_values(const Quadrature< dim > &q, FEInternalData &data, FEValuesData< dim, spacedim > &fv_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure...
Definition: fe_rt.hh:347
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:42
Volume element.
std::vector< arma::mat::fixed< spacedim, dim > > jacobians
Jacobians of the mapping at the quadrature points.
Definition: fe_values.hh:83
#define ASSERT(...)
Definition: global_defs.h:121
Raviart-Thomas element of order 0.
Definition: fe_rt.hh:44
arma::vec::fixed< dim > basis_grad(const unsigned int i, const arma::vec::fixed< dim > &p) const
The scalar variant of basis_grad_vector must be implemented but may not be used.
Definition: fe_rt.hh:177
virtual ~FE_RT0()
Destructor.
Definition: fe_rt.hh:129
UpdateFlags update_each(UpdateFlags flags)
Decides which additional quantities have to be computed for each cell.
Definition: fe_rt.hh:333
void compute_node_matrix()
Computes the conversion matrix from internal basis to shape functions.
Definition: fe_rt.hh:205
Shape function gradients.
std::vector< std::vector< arma::mat > > basis_grad_vectors
Precomputed gradients of basis functions at the quadrature points.
arma::vec::fixed< dim > basis_vector(const unsigned int i, const arma::vec::fixed< dim > &p) const
Returns the ith basis function evaluated at the point p.
Definition: fe_rt.hh:184
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_vectors
Shape functions (for vectorial finite elements) evaluated at quadrature points.
Definition: fe_values.hh:116
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
Structure for storing the precomputed finite element data.
Abstract class for description of finite elements.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Definition: dofhandler.hh:40
static const unsigned int n_raw_functions
Number of raw basis functions.
Definition: fe_rt.hh:59