/usr/include/deal.II/fe/fe_poly.h is in libdeal.ii-dev 8.1.0-6ubuntu1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 | // ---------------------------------------------------------------------
// $Id: fe_poly.h 30036 2013-07-18 16:55:32Z maier $
//
// Copyright (C) 2004 - 2013 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the deal.II distribution.
//
// ---------------------------------------------------------------------
#ifndef __deal2__fe_poly_h
#define __deal2__fe_poly_h
#include <deal.II/fe/fe.h>
DEAL_II_NAMESPACE_OPEN
/*!@addtogroup febase */
/*@{*/
/**
* This class gives a unified framework for the implementation of
* FiniteElement classes based on polynomial spaces like the
* TensorProductPolynomials or PolynomialSpace classes.
*
* Every class conforming to the following interface can be used as
* template parameter POLY.
*
* @code
* static const unsigned int dimension;
*
* double compute_value (const unsigned int i,
* const Point<dim> &p) const;
*
* Tensor<1,dim> compute_grad (const unsigned int i,
* const Point<dim> &p) const;
*
* Tensor<2,dim> compute_grad_grad (const unsigned int i,
* const Point<dim> &p) const;
* @endcode
* Example classes are TensorProductPolynomials, PolynomialSpace or
* PolynomialsP.
*
* This class is not a fully implemented FiniteElement class. Instead
* there are several pure virtual functions declared in the
* FiniteElement and FiniteElement classes which cannot
* implemented by this class but are left for implementation in
* derived classes.
*
* Furthermore, this class assumes that shape functions of the
* FiniteElement under consideration do <em>not</em> depend on the
* actual shape of the cells in real space, i.e. update_once()
* includes <tt>update_values</tt>. For FiniteElements whose shape
* functions depend on the cells in real space, the update_once() and
* update_each() functions must be overloaded.
*
* @todo Since nearly all functions for spacedim != dim are
* specialized, this class needs cleaning up.
*
* @author Ralf Hartmann 2004, Guido Kanschat, 2009
**/
template <class POLY, int dim=POLY::dimension, int spacedim=dim>
class FE_Poly : public FiniteElement<dim,spacedim>
{
public:
/**
* Constructor.
*/
FE_Poly (const POLY &poly_space,
const FiniteElementData<dim> &fe_data,
const std::vector<bool> &restriction_is_additive_flags,
const std::vector<ComponentMask> &nonzero_components);
/**
* Return the polynomial degree of this finite element, i.e. the value
* passed to the constructor.
*/
unsigned int get_degree () const;
/**
* Return the numbering of the underlying polynomial space compared to
* lexicographic ordering of the basis functions. Returns
* POLY::get_numbering().
*/
std::vector<unsigned int> get_poly_space_numbering() const;
/**
* Return the inverse numbering of the underlying polynomial space. Returns
* POLY::get_numbering_inverse().
*/
std::vector<unsigned int> get_poly_space_numbering_inverse() const;
/**
* Return the value of the <tt>i</tt>th shape function at the point
* <tt>p</tt>. See the FiniteElement base class for more information about
* the semantics of this function.
*/
virtual double shape_value (const unsigned int i,
const Point<dim> &p) const;
/**
* Return the value of the <tt>component</tt>th vector component of the
* <tt>i</tt>th shape function at the point <tt>p</tt>. See the
* FiniteElement base class for more information about the semantics of this
* function.
*
* Since this element is scalar, the returned value is the same as if the
* function without the <tt>_component</tt> suffix were called, provided
* that the specified component is zero.
*/
virtual double shape_value_component (const unsigned int i,
const Point<dim> &p,
const unsigned int component) const;
/**
* Return the gradient of the <tt>i</tt>th shape function at the point
* <tt>p</tt>. See the FiniteElement base class for more information about
* the semantics of this function.
*/
virtual Tensor<1,dim> shape_grad (const unsigned int i,
const Point<dim> &p) const;
/**
* Return the gradient of the <tt>component</tt>th vector component of the
* <tt>i</tt>th shape function at the point <tt>p</tt>. See the
* FiniteElement base class for more information about the semantics of this
* function.
*
* Since this element is scalar, the returned value is the same as if the
* function without the <tt>_component</tt> suffix were called, provided
* that the specified component is zero.
*/
virtual Tensor<1,dim> shape_grad_component (const unsigned int i,
const Point<dim> &p,
const unsigned int component) const;
/**
* Return the tensor of second derivatives of the <tt>i</tt>th shape
* function at point <tt>p</tt> on the unit cell. See the FiniteElement base
* class for more information about the semantics of this function.
*/
virtual Tensor<2,dim> shape_grad_grad (const unsigned int i,
const Point<dim> &p) const;
/**
* Return the second derivative of the <tt>component</tt>th vector component
* of the <tt>i</tt>th shape function at the point <tt>p</tt>. See the
* FiniteElement base class for more information about the semantics of this
* function.
*
* Since this element is scalar, the returned value is the same as if the
* function without the <tt>_component</tt> suffix were called, provided
* that the specified component is zero.
*/
virtual Tensor<2,dim> shape_grad_grad_component (const unsigned int i,
const Point<dim> &p,
const unsigned int component) const;
protected:
virtual
typename Mapping<dim,spacedim>::InternalDataBase *
get_data (const UpdateFlags,
const Mapping<dim,spacedim> &mapping,
const Quadrature<dim> &quadrature) const ;
virtual void
fill_fe_values (const Mapping<dim,spacedim> &mapping,
const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const Quadrature<dim> &quadrature,
typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
typename Mapping<dim,spacedim>::InternalDataBase &fe_internal,
FEValuesData<dim,spacedim> &data,
CellSimilarity::Similarity &cell_similarity) const;
virtual void
fill_fe_face_values (const Mapping<dim,spacedim> &mapping,
const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int face_no,
const Quadrature<dim-1> &quadrature,
typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
typename Mapping<dim,spacedim>::InternalDataBase &fe_internal,
FEValuesData<dim,spacedim> &data) const ;
virtual void
fill_fe_subface_values (const Mapping<dim,spacedim> &mapping,
const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int face_no,
const unsigned int sub_no,
const Quadrature<dim-1> &quadrature,
typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
typename Mapping<dim,spacedim>::InternalDataBase &fe_internal,
FEValuesData<dim,spacedim> &data) const ;
/**
* Determine the values that need to be computed on the unit cell to be able
* to compute all values required by <tt>flags</tt>.
*
* For the purpuse of this function, refer to the documentation in
* FiniteElement.
*
* This class assumes that shape functions of this FiniteElement do
* <em>not</em> depend on the actual shape of the cells in real
* space. Therefore, the effect in this element is as follows: if
* <tt>update_values</tt> is set in <tt>flags</tt>, copy it to the
* result. All other flags of the result are cleared, since everything else
* must be computed for each cell.
*/
virtual UpdateFlags update_once (const UpdateFlags flags) const;
/**
* Determine the values that need to be computed on every cell to be able to
* compute all values required by <tt>flags</tt>.
*
* For the purpuse of this function, refer to the documentation in
* FiniteElement.
*
* This class assumes that shape functions of this FiniteElement do
* <em>not</em> depend on the actual shape of the cells in real space.
*
* The effect in this element is as follows:
* <ul>
*
* <li> if <tt>update_gradients</tt> is set, the result will contain
* <tt>update_gradients</tt> and <tt>update_covariant_transformation</tt>.
* The latter is required to transform the gradient on the unit cell to the
* real cell. Remark, that the action required by
* <tt>update_covariant_transformation</tt> is actually performed by the
* Mapping object used in conjunction with this finite element.
*
* <li> if <tt>update_hessians</tt> is set, the result will contain
* <tt>update_hessians</tt> and <tt>update_covariant_transformation</tt>.
* The rationale is the same as above and no higher derivatives of the
* transformation are required, since we use difference quotients for the
* actual computation.
*
* </ul>
*/
virtual UpdateFlags update_each (const UpdateFlags flags) const;
/**
* Fields of cell-independent data.
*
* For information about the general purpose of this class, see the
* documentation of the base class.
*/
class InternalData : public FiniteElement<dim,spacedim>::InternalDataBase
{
public:
/**
* Array with shape function values in quadrature points. There is one row
* for each shape function, containing values for each quadrature point.
*
* In this array, we store the values of the shape function in the
* quadrature points on the unit cell. Since these values do not change
* under transformation to the real cell, we only need to copy them over
* when visiting a concrete cell.
*/
std::vector<std::vector<double> > shape_values;
/**
* Array with shape function gradients in quadrature points. There is one
* row for each shape function, containing values for each quadrature
* point.
*
* We store the gradients in the quadrature points on the unit cell. We
* then only have to apply the transformation (which is a matrix-vector
* multiplication) when visiting an actual cell.
*/
std::vector<std::vector<Tensor<1,dim> > > shape_gradients;
};
/**
* The polynomial space. Its type is given by the template parameter POLY.
*/
POLY poly_space;
};
/*@}*/
DEAL_II_NAMESPACE_CLOSE
#endif
|