/usr/include/deal.II/fe/fe_face.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 452 453 454 455 456 457 458 459 460 461 462 463 | // ---------------------------------------------------------------------
//
// Copyright (C) 2009 - 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_face_h
#define dealii__fe_face_h
#include <deal.II/base/config.h>
#include <deal.II/base/tensor_product_polynomials.h>
#include <deal.II/base/polynomial_space.h>
#include <deal.II/fe/fe_poly_face.h>
DEAL_II_NAMESPACE_OPEN
/**
* A finite element, which is a tensor product polynomial on each face and
* undefined in the interior of the cells. The basis functions on the faces
* are Lagrange polynomials based on the support points of the
* (dim-1)-dimensional Gauss--Lobatto quadrature rule. For element degree one
* and two, the polynomials hence correspond to the usual Lagrange polynomials
* on equidistant points.
*
* Although the name does not give it away, the element is discontinuous at
* locations where faces of cells meet. In particular, this finite element is
* the trace space of FE_RaviartThomas on the faces and serves in hybridized
* methods, e.g. in combination with the FE_DGQ element. Its use is
* demonstrated in the step-51 tutorial program.
*
* @note Since this element is defined only on faces, only FEFaceValues and
* FESubfaceValues will be able to extract reasonable values from any face
* polynomial. In order to make the use of FESystem simpler, using a (cell)
* FEValues object will not fail using this finite element space, but all
* shape function values extracted will be equal to zero.
*
* @ingroup fe
* @author Guido Kanschat, Martin Kronbichler
* @date 2009, 2011, 2013
*/
template <int dim, int spacedim=dim>
class FE_FaceQ : public FE_PolyFace<TensorProductPolynomials<dim-1>, dim, spacedim>
{
public:
/**
* Constructor for tensor product polynomials of degree <tt>p</tt>. The
* shape functions created using this constructor correspond to Lagrange
* polynomials in each coordinate direction.
*/
FE_FaceQ (const unsigned int p);
virtual FiniteElement<dim,spacedim> *clone() const;
/**
* Return a string that uniquely identifies a finite element. This class
* returns <tt>FE_FaceQ<dim>(degree)</tt>, with <tt>dim</tt> and
* <tt>degree</tt> replaced by appropriate values.
*/
virtual std::string get_name () const;
/**
* Return the matrix interpolating from a face of 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>. This
* element 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 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>. This
* element 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;
/**
* This function returns @p true, if the shape function @p shape_index has
* non-zero function values somewhere on the face @p face_index.
*/
virtual bool has_support_on_face (const unsigned int shape_index,
const unsigned int face_index) 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;
/**
* 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 FiniteElementBase::Domination and in
* particular the
* @ref hp_paper "hp paper".
*/
virtual
FiniteElementDomination::Domination
compare_for_face_domination (const FiniteElement<dim,spacedim> &fe_other) const;
/**
* Returns a list of constant modes of the element. For this element, it
* simply returns one row with all entries set to true.
*/
virtual std::pair<Table<2,bool>, std::vector<unsigned int> >
get_constant_modes () const;
private:
/**
* Return vector with dofs per vertex, line, quad, hex.
*/
static std::vector<unsigned int> get_dpo_vector (const unsigned int deg);
};
/**
* Specialization of FE_FaceQ for 1D. In that case, the finite element only
* consists of one degree of freedom in each of the two faces (= vertices) of
* a cell, irrespective of the degree. However, this element still accepts a
* degree in its constructor and also returns that degree. This way,
* dimension-independent programming with trace elements is also possible in
* 1D (even though there is no computational benefit at all from it in 1D).
*
* @ingroup fe
* @author Guido Kanschat, Martin Kronbichler
* @date 2014
*/
template <int spacedim>
class FE_FaceQ<1,spacedim> : public FiniteElement<1,spacedim>
{
public:
/**
* Constructor.
*/
FE_FaceQ (const unsigned int p);
/**
* Clone method.
*/
virtual FiniteElement<1,spacedim> *clone() const;
/**
* Return a string that uniquely identifies a finite element. This class
* returns <tt>FE_FaceQ<dim>(degree)</tt>, with <tt>dim</tt> and
* <tt>degree</tt> replaced by appropriate values.
*/
virtual std::string get_name () const;
// for documentation, see the FiniteElement base class
virtual
UpdateFlags
requires_update_flags (const UpdateFlags update_flags) const;
/**
* Return the matrix interpolating from a face of 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>. This
* element 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<1,spacedim> &source,
FullMatrix<double> &matrix) const;
/**
* Return the matrix interpolating from a face of 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>. This
* element 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<1,spacedim> &source,
const unsigned int subface,
FullMatrix<double> &matrix) const;
/**
* This function returns @p true, if the shape function @p shape_index has
* non-zero function values somewhere on the face @p face_index.
*/
virtual bool has_support_on_face (const unsigned int shape_index,
const unsigned int face_index) 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;
/**
* 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 FiniteElementBase::Domination and in
* particular the
* @ref hp_paper "hp paper".
*/
virtual
FiniteElementDomination::Domination
compare_for_face_domination (const FiniteElement<1,spacedim> &fe_other) const;
/**
* Returns a list of constant modes of the element. For this element, it
* simply returns one row with all entries set to true.
*/
virtual std::pair<Table<2,bool>, std::vector<unsigned int> >
get_constant_modes () const;
protected:
/*
* NOTE: The following functions have their definitions inlined into the class declaration
* because we otherwise run into a compiler error with MS Visual Studio.
*/
virtual
typename FiniteElement<1,spacedim>::InternalDataBase *
get_data (const UpdateFlags /*update_flags*/,
const Mapping<1,spacedim> &/*mapping*/,
const Quadrature<1> &/*quadrature*/,
dealii::internal::FEValues::FiniteElementRelatedData<1, spacedim> &/*output_data*/) const
{
return new typename FiniteElement<1, spacedim>::InternalDataBase;
}
typename FiniteElement<1,spacedim>::InternalDataBase *
get_face_data(const UpdateFlags update_flags,
const Mapping<1,spacedim> &/*mapping*/,
const Quadrature<0> &quadrature,
dealii::internal::FEValues::FiniteElementRelatedData<1, spacedim> &/*output_data*/) const
{
// generate a new data object and initialize some fields
typename FiniteElement<1,spacedim>::InternalDataBase *data =
new typename FiniteElement<1,spacedim>::InternalDataBase;
data->update_each = requires_update_flags(update_flags);
const unsigned int n_q_points = quadrature.size();
AssertDimension(n_q_points, 1);
(void)n_q_points;
// No derivatives of this element are implemented.
if (data->update_each & update_gradients || data->update_each & update_hessians)
{
Assert(false, ExcNotImplemented());
}
return data;
}
typename FiniteElement<1,spacedim>::InternalDataBase *
get_subface_data(const UpdateFlags update_flags,
const Mapping<1,spacedim> &mapping,
const Quadrature<0> &quadrature,
dealii::internal::FEValues::FiniteElementRelatedData<1, spacedim> &output_data) const
{
return get_face_data(update_flags, mapping, quadrature, output_data);
}
virtual
void
fill_fe_values (const typename Triangulation<1,spacedim>::cell_iterator &cell,
const CellSimilarity::Similarity cell_similarity,
const Quadrature<1> &quadrature,
const Mapping<1,spacedim> &mapping,
const typename Mapping<1,spacedim>::InternalDataBase &mapping_internal,
const dealii::internal::FEValues::MappingRelatedData<1, spacedim> &mapping_data,
const typename FiniteElement<1,spacedim>::InternalDataBase &fe_internal,
dealii::internal::FEValues::FiniteElementRelatedData<1, spacedim> &output_data) const;
virtual
void
fill_fe_face_values (const typename Triangulation<1,spacedim>::cell_iterator &cell,
const unsigned int face_no,
const Quadrature<0> &quadrature,
const Mapping<1,spacedim> &mapping,
const typename Mapping<1,spacedim>::InternalDataBase &mapping_internal,
const dealii::internal::FEValues::MappingRelatedData<1, spacedim> &mapping_data,
const typename FiniteElement<1,spacedim>::InternalDataBase &fe_internal,
dealii::internal::FEValues::FiniteElementRelatedData<1, spacedim> &output_data) const;
virtual
void
fill_fe_subface_values (const typename Triangulation<1,spacedim>::cell_iterator &cell,
const unsigned int face_no,
const unsigned int sub_no,
const Quadrature<0> &quadrature,
const Mapping<1,spacedim> &mapping,
const typename Mapping<1,spacedim>::InternalDataBase &mapping_internal,
const dealii::internal::FEValues::MappingRelatedData<1, spacedim> &mapping_data,
const typename FiniteElement<1,spacedim>::InternalDataBase &fe_internal,
dealii::internal::FEValues::FiniteElementRelatedData<1, spacedim> &output_data) const;
private:
/**
* Return vector with dofs per vertex, line, quad, hex.
*/
static
std::vector<unsigned int>
get_dpo_vector (const unsigned int deg);
};
/**
* A finite element, which is a Legendre element of complete polynomials on
* each face (i.e., it is the face equivalent of what FE_DGP is on cells) and
* undefined in the interior of the cells. The basis functions on the faces
* are from Polynomials::Legendre.
*
* Although the name does not give it away, the element is discontinuous at
* locations where faces of cells meet. The element serves in hybridized
* methods, e.g. in combination with the FE_DGP element. An example of
* hybridizes methods can be found in the step-51 tutorial program.
*
* @note Since this element is defined only on faces, only FEFaceValues and
* FESubfaceValues will be able to extract reasonable values from any face
* polynomial. In order to make the use of FESystem simpler, using a (cell)
* FEValues object will not fail using this finite element space, but all
* shape function values extracted will be equal to zero.
*
* @ingroup fe
* @author Martin Kronbichler
* @date 2013
*/
template <int dim, int spacedim=dim>
class FE_FaceP : public FE_PolyFace<PolynomialSpace<dim-1>, dim, spacedim>
{
public:
/**
* Constructor for complete basis of polynomials of degree <tt>p</tt>. The
* shape functions created using this constructor correspond to Legendre
* polynomials in each coordinate direction.
*/
FE_FaceP(unsigned int p);
/**
* Clone method.
*/
virtual FiniteElement<dim,spacedim> *clone() const;
/**
* Return a string that uniquely identifies a finite element. This class
* returns <tt>FE_FaceP<dim>(degree)</tt> , with <tt>dim</tt> and
* <tt>degree</tt> replaced by appropriate values.
*/
virtual std::string get_name () const;
/**
* Return the matrix interpolating from a face of 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>. This
* element 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 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>. This
* element 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;
/**
* This function returns @p true, if the shape function @p shape_index has
* non-zero function values somewhere on the face @p face_index.
*/
virtual bool has_support_on_face (const unsigned int shape_index,
const unsigned int face_index) 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;
/**
* 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 FiniteElementBase::Domination and in
* particular the
* @ref hp_paper "hp paper".
*/
virtual
FiniteElementDomination::Domination
compare_for_face_domination (const FiniteElement<dim,spacedim> &fe_other) const;
/**
* Returns a list of constant modes of the element. For this element, the
* first entry on each face is true, all other are false (as the constant
* function is represented by the first base function of Legendre
* polynomials).
*/
virtual std::pair<Table<2,bool>, std::vector<unsigned int> >
get_constant_modes () const;
private:
/**
* Return vector with dofs per vertex, line, quad, hex.
*/
static std::vector<unsigned int> get_dpo_vector (const unsigned int deg);
};
/**
* FE_FaceP in 1D, i.e., with degrees of freedom on the element vertices.
*/
template <int spacedim>
class FE_FaceP<1,spacedim> : public FE_FaceQ<1,spacedim>
{
public:
/**
* Constructor.
*/
FE_FaceP (const unsigned int p);
/**
* Returns the name of the element
*/
std::string get_name() const;
};
DEAL_II_NAMESPACE_CLOSE
#endif
|