Flow123d  jenkins-Flow123d-linux-release-multijob-282
finite_element.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 Abstract class for description of finite elements.
27  * @author Jan Stebel
28  */
29 
30 
31 #include "system/system.hh"
32 #include "quadrature/quadrature.hh"
33 #include "fem/dofhandler.hh"
34 #include "fem/finite_element.hh"
35 #include "fem/fe_values.hh"
36 
37 
38 
39 using namespace std;
40 
41 
42 
43 
44 template<unsigned int dim, unsigned int spacedim>
46 {
47  init();
48 }
49 
50 template<unsigned int dim, unsigned int spacedim>
52 {
53  number_of_dofs = 0;
54  is_scalar_fe = true;
55  for (unsigned int i = 0; i <= dim; i++)
56  {
57  number_of_single_dofs[i] = 0;
58  number_of_pairs[i] = 0;
59  number_of_triples[i] = 0;
60  number_of_sextuples[i] = 0;
61  }
62 }
63 
64 template<unsigned int dim, unsigned int spacedim> inline
65 const unsigned int FiniteElement<dim,spacedim>::n_dofs() const
66 {
67  return number_of_dofs;
68 }
69 
70 template<unsigned int dim, unsigned int spacedim> inline
72  unsigned int object_dim, DofMultiplicity multiplicity)
73 {
74  ASSERT(object_dim >= 0 && object_dim <= dim,
75  "Object type number is out of range.");
76  switch (multiplicity)
77  {
78  case DOF_SINGLE:
79  return number_of_single_dofs[object_dim];
80  case DOF_PAIR:
81  return number_of_pairs[object_dim];
82  case DOF_TRIPLE:
83  return number_of_triples[object_dim];
84  case DOF_SEXTUPLE:
85  return number_of_sextuples[object_dim];
86  }
87 
88  return 0;
89 }
90 
91 template<unsigned int dim, unsigned int spacedim> inline
93 {
94  ASSERT_EQUAL(get_generalized_support_points().size(), number_of_dofs);
95 
96  arma::mat M(number_of_dofs, number_of_dofs);
97 
98  for (unsigned int i = 0; i < number_of_dofs; i++)
99  for (unsigned int j = 0; j < number_of_dofs; j++) {
100  M(j, i) = basis_value(j, get_generalized_support_points()[i]);
101 
102  }
103  node_matrix = arma::inv(M);
104 }
105 
106 template<unsigned int dim, unsigned int spacedim>
108 {
109  FEInternalData *data = new FEInternalData;
110 
111  if (flags & update_values)
112  {
113  arma::vec values(number_of_dofs);
114  data->basis_values.resize(q.size());
115  for (unsigned int i=0; i<q.size(); i++)
116  {
117  for (unsigned int j=0; j<number_of_dofs; j++)
118  values[j] = basis_value(j, q.point(i));
119  data->basis_values[i] = node_matrix * values;
120  }
121  }
122 
123  if (flags & update_gradients)
124  {
125  arma::mat grads(number_of_dofs, dim);
126  data->basis_grads.resize(q.size());
127  for (unsigned int i=0; i<q.size(); i++)
128  {
129  for (unsigned int j=0; j<number_of_dofs; j++)
130  grads.row(j) = arma::trans(basis_grad(j, q.point(i)));
131  data->basis_grads[i] = node_matrix * grads;
132  }
133  }
134 
135  return data;
136 }
137 
138 template<unsigned int dim, unsigned int spacedim> inline
140 {
141  UpdateFlags f = flags;
142 
143  if (flags & update_gradients)
145 
146  return f;
147 }
148 
149 template<unsigned int dim, unsigned int spacedim> inline
151  const Quadrature<dim> &q,
152  FEInternalData &data,
154 {
155  // shape values
156  if (fv_data.update_flags & update_values)
157  {
158  for (unsigned int i = 0; i < q.size(); i++)
159  fv_data.shape_values[i] = data.basis_values[i];
160  }
161 
162  // shape gradients
163  if (fv_data.update_flags & update_gradients)
164  {
165  for (unsigned int i = 0; i < q.size(); i++)
166  {
167  fv_data.shape_gradients[i] = data.basis_grads[i] * fv_data.inverse_jacobians[i];
168  }
169  }
170 }
171 
172 template<unsigned int dim, unsigned int spacedim>
174 {
175  if (generalized_support_points.size() > 0)
176  {
177  return generalized_support_points;
178  }
179  else
180  {
181  return unit_support_points;
182  }
183 }
184 
185 
186 template<unsigned int dim, unsigned int spacedim>
188 {}
189 
190 
191 template class FiniteElement<0,3>;
192 template class FiniteElement<1,3>;
193 template class FiniteElement<2,3>;
194 template class FiniteElement<3,3>;
195 
196 
std::vector< arma::vec > shape_values
Shape functions evaluated at the quadrature points.
Definition: fe_values.hh:103
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
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
const std::vector< arma::vec::fixed< dim > > & get_generalized_support_points()
Returns either the generalized support points (if they are defined) or the unit support points...
void init()
Clears all internal structures.
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
Class FEValues calculates finite element data on the actual cells such as shape function values...
std::vector< arma::mat > basis_grads
Precomputed gradients of basis functions at the quadrature points.
Base class for quadrature rules on simplices in arbitrary dimensions.
Definition: fe_values.hh:42
virtual UpdateFlags update_each(UpdateFlags flags)
Decides which additional quantities have to be computed for each cell.
virtual ~FiniteElement()
Destructor.
#define ASSERT(...)
Definition: global_defs.h:121
#define ASSERT_EQUAL(a, b)
Definition: global_defs.h:136
virtual void compute_node_matrix()
Initializes the node_matrix for computing the coefficients of the raw basis functions from values at ...
std::vector< arma::mat > shape_gradients
Gradients of shape functions evaluated at the quadrature points.
Definition: fe_values.hh:110
Basic definitions of numerical quadrature rules.
Shape function gradients.
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Definition: fe_values.hh:57
const unsigned int n_object_dofs(unsigned int object_dim, DofMultiplicity multiplicity)
Returns the number of single dofs/dof pairs/triples/sextuples that lie on a single geometric entity o...
const unsigned int size() const
Returns number of quadrature points.
Definition: quadrature.hh:136
std::vector< arma::vec > basis_values
Precomputed values of basis functions at the quadrature points.
DofMultiplicity
Multiplicity of finite element dofs.
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
Definition: quadrature.hh:141
FiniteElement()
Constructor.
const unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
virtual FEInternalData * initialize(const Quadrature< dim > &q, UpdateFlags flags)
Calculates the data on the reference cell.
Structure for storing the precomputed finite element data.
Abstract class for description of finite elements.
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...
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Definition: dofhandler.hh:40