/usr/include/deal.II/fe/fe_bernstein.h is in libdeal.ii-dev 8.5.1-3.
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 | // ---------------------------------------------------------------------
//
// Copyright (C) 2000 - 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_bernstein_h
#define dealii__fe_bernstein_h
#include <deal.II/base/config.h>
#include <deal.II/base/tensor_product_polynomials.h>
#include <deal.II/fe/fe_q_base.h>
DEAL_II_NAMESPACE_OPEN
/*!@addtogroup fe */
/*@{*/
/**
* Implementation of a scalar Bernstein finite element @p that we call
* FE_Bernstein in analogy with FE_Q that yields the finite element space of
* continuous, piecewise Bernstein polynomials of degree @p p in each
* coordinate direction. This class is realized using tensor product
* polynomials of Bernstein basis polynomials.
*
*
* The standard constructor of this class takes the degree @p p of this finite
* element.
*
* For more information about the <tt>spacedim</tt> template parameter check
* the documentation of FiniteElement or the one of Triangulation.
*
* <h3>Implementation</h3>
*
* The constructor creates a TensorProductPolynomials object that includes the
* tensor product of @p Bernstein polynomials of degree @p p. This @p
* TensorProductPolynomials object provides all values and derivatives of the
* shape functions.
*
* <h3>Numbering of the degrees of freedom (DoFs)</h3>
*
* The original ordering of the shape functions represented by the
* TensorProductPolynomials is a tensor product numbering. However, the shape
* functions on a cell are renumbered beginning with the shape functions whose
* support points are at the vertices, then on the line, on the quads, and
* finally (for 3d) on the hexes. See the documentation of FE_Q for more
* details.
*
*
* @author Marco Tezzele, Luca Heltai
* @date 2013, 2015
*/
template <int dim, int spacedim=dim>
class FE_Bernstein : public FE_Q_Base<TensorProductPolynomials<dim>,dim,spacedim>
{
public:
/**
* Constructor for tensor product polynomials of degree @p p.
*/
FE_Bernstein (const unsigned int p);
/**
* FE_Bernstein is not interpolatory in the element interior, which prevents
* this element from defining an interpolation matrix. An exception will be
* thrown.
*
* Overrides the implementation from FE_Q_Base.
*/
virtual void
get_interpolation_matrix (const FiniteElement<dim,spacedim> &source,
FullMatrix<double> &matrix) const;
/**
* FE_Bernstein is not interpolatory in the element interior, which prevents
* this element from defining a restriction matrix. An exception will be
* thrown.
*
* Overrides the implementation from FE_Q_Base.
*/
virtual const FullMatrix<double> &
get_restriction_matrix (const unsigned int child,
const RefinementCase<dim> &refinement_case=RefinementCase<dim>::isotropic_refinement) const;
/**
* FE_Bernstein is not interpolatory in the element interior, which prevents
* this element from defining a prolongation matrix. An exception will be
* thrown.
*
* Overrides the implementation from FE_Q_Base.
*/
virtual const FullMatrix<double> &
get_prolongation_matrix (const unsigned int child,
const RefinementCase<dim> &refinement_case=RefinementCase<dim>::isotropic_refinement) const;
/**
* Return the matrix interpolating from a face of one element to the face of
* the neighboring element. The size of the matrix is then
* <tt>source.dofs_per_face</tt> times <tt>this->dofs_per_face</tt>. The
* FE_Bernstein element family only provides interpolation matrices for
* elements of the same type and FE_Nothing. For all other elements, an
* exception of type
* FiniteElement<dim,spacedim>::ExcInterpolationNotImplemented is thrown.
*/
virtual void
get_face_interpolation_matrix (const FiniteElement<dim,spacedim> &source,
FullMatrix<double> &matrix) const;
/**
* Return the matrix interpolating from a face of one element to the face of
* the neighboring element. The size of the matrix is then
* <tt>source.dofs_per_face</tt> times <tt>this->dofs_per_face</tt>. The
* FE_Bernstein element family only provides interpolation matrices for
* elements of the same type and FE_Nothing. For all other elements, an
* exception of type
* FiniteElement<dim,spacedim>::ExcInterpolationNotImplemented is thrown.
*/
virtual void
get_subface_interpolation_matrix (const FiniteElement<dim,spacedim> &source,
const unsigned int subface,
FullMatrix<double> &matrix) const;
/**
* Return whether this element implements its hanging node constraints in
* the new way, which has to be used to make elements "hp compatible".
*/
virtual bool hp_constraints_are_implemented () const;
/**
* If, on a vertex, several finite elements are active, the hp code first
* assigns the degrees of freedom of each of these FEs different global
* indices. It then calls this function to find out which of them should get
* identical values, and consequently can receive the same global DoF index.
* This function therefore returns a list of identities between DoFs of the
* present finite element object with the DoFs of @p fe_other, which is a
* reference to a finite element object representing one of the other finite
* elements active on this particular vertex. The function computes which of
* the degrees of freedom of the two finite element objects are equivalent,
* both numbered between zero and the corresponding value of dofs_per_vertex
* of the two finite elements. The first index of each pair denotes one of
* the vertex dofs of the present element, whereas the second is the
* corresponding index of the other finite element.
*/
virtual
std::vector<std::pair<unsigned int, unsigned int> >
hp_vertex_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
/**
* Same as hp_vertex_dof_indices(), except that the function treats degrees
* of freedom on lines.
*/
virtual
std::vector<std::pair<unsigned int, unsigned int> >
hp_line_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
/**
* Same as hp_vertex_dof_indices(), except that the function treats degrees
* of freedom on quads.
*/
virtual
std::vector<std::pair<unsigned int, unsigned int> >
hp_quad_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
/**
* Return whether this element dominates the one given as argument when they
* meet at a common face, whether it is the other way around, whether
* neither dominates, or if either could dominate.
*
* For a definition of domination, see FiniteElementDomination::Domination
* and in particular the
* @ref hp_paper "hp paper".
*/
virtual
FiniteElementDomination::Domination
compare_for_face_domination (const FiniteElement<dim,spacedim> &fe_other) const;
/**
* Return a string that uniquely identifies a finite element. This class
* returns <tt>FE_Bernstein<dim>(degree)</tt>, with @p dim and @p degree
* replaced by appropriate values.
*/
virtual std::string get_name () const;
protected:
/**
* @p clone function instead of a copy constructor.
*
* This function is needed by the constructors of @p FESystem.
*/
virtual FiniteElement<dim,spacedim> *clone() const;
/**
* Only for internal use. Its full name is @p get_dofs_per_object_vector
* function and it creates the @p dofs_per_object vector that is needed
* within the constructor to be passed to the constructor of @p
* FiniteElementData.
*/
static std::vector<unsigned int> get_dpo_vector(const unsigned int degree);
/**
* This function renumbers Bernstein basis functions from hierarchic to
* lexicographic numbering.
*/
TensorProductPolynomials<dim> renumber_bases(const unsigned int degree);
};
/*@}*/
DEAL_II_NAMESPACE_CLOSE
#endif
|