/usr/include/deal.II/fe/fe_poly.h is in libdeal.ii-dev 8.4.2-2+b1.
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 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 | // ---------------------------------------------------------------------
//
// Copyright (C) 2004 - 2016 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 dealii__fe_poly_h
#define dealii__fe_poly_h
#include <deal.II/fe/fe.h>
#include <deal.II/base/quadrature.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 PolynomialType.
*
* @code
* static const unsigned int dimension;
*
* void compute (const Point<dim> &unit_point,
* std::vector<double> &values,
* std::vector<Tensor<1,dim> > &grads,
* std::vector<Tensor<2,dim> > &grad_grads,
* std::vector<Tensor<3,dim> > &third_derivatives,
* std::vector<Tensor<4,dim> > &fourth_derivatives) const;
*
* double compute_value (const unsigned int i,
* const Point<dim> &p) const;
*
* template <int order>
* Tensor<order,dim> compute_derivative (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 be implemented by this class but are
* left for implementation in derived classes.
*
* @todo Since nearly all functions for spacedim != dim are specialized, this
* class needs cleaning up.
*
* @author Ralf Hartmann 2004, Guido Kanschat, 2009
*/
template <class PolynomialType, int dim=PolynomialType::dimension, int spacedim=dim>
class FE_Poly : public FiniteElement<dim,spacedim>
{
public:
/**
* Constructor.
*/
FE_Poly (const PolynomialType &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;
// for documentation, see the FiniteElement base class
virtual
UpdateFlags
requires_update_flags (const UpdateFlags update_flags) const;
/**
* Return the numbering of the underlying polynomial space compared to
* lexicographic ordering of the basis functions. Returns
* PolynomialType::get_numbering().
*/
std::vector<unsigned int> get_poly_space_numbering() const;
/**
* Return the inverse numbering of the underlying polynomial space. Returns
* PolynomialType::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;
/**
* Return the tensor of third 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<3,dim> shape_3rd_derivative (const unsigned int i,
const Point<dim> &p) const;
/**
* Return the third 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<3,dim> shape_3rd_derivative_component (const unsigned int i,
const Point<dim> &p,
const unsigned int component) const;
/**
* Return the tensor of fourth 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<4,dim> shape_4th_derivative (const unsigned int i,
const Point<dim> &p) const;
/**
* Return the fourth 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<4,dim> shape_4th_derivative_component (const unsigned int i,
const Point<dim> &p,
const unsigned int component) const;
protected:
/*
* NOTE: The following function has its definition inlined into the class declaration
* because we otherwise run into a compiler error with MS Visual Studio.
*/
virtual
typename FiniteElement<dim,spacedim>::InternalDataBase *
get_data(const UpdateFlags update_flags,
const Mapping<dim,spacedim> &/*mapping*/,
const Quadrature<dim> &quadrature,
dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const
{
// generate a new data object and
// initialize some fields
InternalData *data = new InternalData;
data->update_each = requires_update_flags(update_flags);
const unsigned int n_q_points = quadrature.size();
// initialize some scratch arrays. we need them for the underlying
// polynomial to put the values and derivatives of shape functions
// to put there, depending on what the user requested
std::vector<double> values(update_flags & update_values ?
this->dofs_per_cell : 0);
std::vector<Tensor<1,dim> > grads(update_flags & update_gradients ?
this->dofs_per_cell : 0);
std::vector<Tensor<2,dim> > grad_grads(update_flags & update_hessians ?
this->dofs_per_cell : 0);
std::vector<Tensor<3,dim> > third_derivatives(update_flags & update_3rd_derivatives ?
this->dofs_per_cell : 0);
std::vector<Tensor<4,dim> > fourth_derivatives; // won't be needed, so leave empty
// now also initialize fields the fields of this class's own
// temporary storage, depending on what we need for the given
// update flags.
//
// there is one exception from the rule: if we are dealing with
// cells (i.e., if this function is not called via
// get_(sub)face_data()), then we can already store things in the
// final location where FEValues::reinit() later wants to see
// things. we then don't need the intermediate space. we determine
// whether we are on a cell by asking whether the number of
// elements in the output array equals the number of quadrature
// points (yes, it's a cell) or not (because in that case the
// number of quadrature points we use here equals the number of
// quadrature points summed over *all* faces or subfaces, whereas
// the number of output slots equals the number of quadrature
// points on only *one* face)
if ((update_flags & update_values)
&&
!((output_data.shape_values.n_rows() > 0)
&&
(output_data.shape_values.n_cols() == n_q_points)))
data->shape_values.reinit (this->dofs_per_cell, n_q_points);
if (update_flags & update_gradients)
data->shape_gradients.reinit (this->dofs_per_cell, n_q_points);
if (update_flags & update_hessians)
data->shape_hessians.reinit (this->dofs_per_cell, n_q_points);
if (update_flags & update_3rd_derivatives)
data->shape_3rd_derivatives.reinit (this->dofs_per_cell, n_q_points);
// next already fill those fields of which we have information by
// now. note that the shape gradients are only those on the unit
// cell, and need to be transformed when visiting an actual cell
if (update_flags & (update_values | update_gradients
| update_hessians | update_3rd_derivatives) )
for (unsigned int i=0; i<n_q_points; ++i)
{
poly_space.compute(quadrature.point(i),
values, grads, grad_grads,
third_derivatives,
fourth_derivatives);
// the values of shape functions at quadrature points don't change.
// consequently, write these values right into the output array if
// we can, i.e., if the output array has the correct size. this is
// the case on cells. on faces, we already precompute data on *all*
// faces and subfaces, but we later on copy only a portion of it
// into the output object; in that case, copy the data from all
// faces into the scratch object
if (update_flags & update_values)
if (output_data.shape_values.n_rows() > 0)
{
if (output_data.shape_values.n_cols() == n_q_points)
for (unsigned int k=0; k<this->dofs_per_cell; ++k)
output_data.shape_values[k][i] = values[k];
else
for (unsigned int k=0; k<this->dofs_per_cell; ++k)
data->shape_values[k][i] = values[k];
}
// for everything else, derivatives need to be transformed,
// so we write them into our scratch space and only later
// copy stuff into where FEValues wants it
if (update_flags & update_gradients)
for (unsigned int k=0; k<this->dofs_per_cell; ++k)
data->shape_gradients[k][i] = grads[k];
if (update_flags & update_hessians)
for (unsigned int k=0; k<this->dofs_per_cell; ++k)
data->shape_hessians[k][i] = grad_grads[k];
if (update_flags & update_3rd_derivatives)
for (unsigned int k=0; k<this->dofs_per_cell; ++k)
data->shape_3rd_derivatives[k][i] = third_derivatives[k];
}
return data;
}
virtual
void
fill_fe_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const CellSimilarity::Similarity cell_similarity,
const Quadrature<dim> &quadrature,
const Mapping<dim,spacedim> &mapping,
const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
const dealii::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const;
virtual
void
fill_fe_face_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int face_no,
const Quadrature<dim-1> &quadrature,
const Mapping<dim,spacedim> &mapping,
const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
const dealii::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) const;
virtual
void
fill_fe_subface_values (const typename Triangulation<dim,spacedim>::cell_iterator &cell,
const unsigned int face_no,
const unsigned int sub_no,
const Quadrature<dim-1> &quadrature,
const Mapping<dim,spacedim> &mapping,
const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
const dealii::internal::FEValues::MappingRelatedData<dim, spacedim> &mapping_data,
const typename FiniteElement<dim,spacedim>::InternalDataBase &fe_internal,
dealii::internal::FEValues::FiniteElementRelatedData<dim, spacedim> &output_data) 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.
*/
Table<2,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.
*/
Table<2,Tensor<1,dim> > shape_gradients;
/**
* Array with shape function hessians in quadrature points. There is one
* row for each shape function, containing values for each quadrature
* point.
*
* We store the hessians in the quadrature points on the unit cell. We
* then only have to apply the transformation when visiting an actual
* cell.
*/
Table<2,Tensor<2,dim> > shape_hessians;
/**
* Array with shape function third derivatives in quadrature points. There
* is one row for each shape function, containing values for each
* quadrature point.
*
* We store the third derivatives in the quadrature points on the unit
* cell. We then only have to apply the transformation when visiting an
* actual cell.
*/
Table<2,Tensor<3,dim> > shape_3rd_derivatives;
};
/**
* Correct the shape third derivatives by subtracting the terms
* corresponding to the Jacobian pushed forward gradient and second
* derivative.
*
* Before the correction, the third derivatives would be given by
* @f[
* D_{ijkl} = \frac{d^3\phi_i}{d \hat x_J d \hat x_K d \hat x_L} (J_{jJ})^{-1} (J_{kK})^{-1} (J_{lL})^{-1},
* @f]
* where $J_{iI}=\frac{d x_i}{d \hat x_I}$. After the correction, the
* correct third derivative would be given by
* @f[
* \frac{d^3\phi_i}{d x_j d x_k d x_l} = D_{ijkl} - H_{mjl} \frac{d^2 \phi_i}{d x_k d x_m}
* - H_{mkl} \frac{d^2 \phi_i}{d x_j d x_m} - H_{mjk} \frac{d^2 \phi_i}{d x_l d x_m}
* - K_{mjkl} \frac{d \phi_i}{d x_m},
* @f]
* where $H_{ijk}$ is the Jacobian pushed-forward derivative and $K_{ijkl}$
* is the Jacobian pushed-forward second derivative.
*/
void
correct_third_derivatives (internal::FEValues::FiniteElementRelatedData<dim,spacedim> &output_data,
const internal::FEValues::MappingRelatedData<dim,spacedim> &mapping_data,
const unsigned int n_q_points,
const unsigned int dof) const;
/**
* The polynomial space. Its type is given by the template parameter
* PolynomialType.
*/
PolynomialType poly_space;
};
/*@}*/
DEAL_II_NAMESPACE_CLOSE
#endif
|