Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
FE_NedelecSZ< dim, spacedim >::InternalData Class Reference

#include <deal.II/fe/fe_nedelec_sz.h>

Inheritance diagram for FE_NedelecSZ< dim, spacedim >::InternalData:
[legend]

Public Member Functions

virtual std::size_t memory_consumption () const
 

Public Attributes

std::vector< std::vector< Tensor< 1, dim > > > shape_values
 
std::vector< std::vector< DerivativeForm< 1, dim, dim > > > shape_grads
 
std::vector< std::vector< DerivativeForm< 2, dim, dim > > > shape_hessians
 
std::vector< std::vector< std::vector< double > > > sigma_imj_values
 
std::vector< std::vector< std::vector< double > > > sigma_imj_grads
 
std::vector< std::vector< double > > edge_sigma_values
 
std::vector< std::vector< double > > edge_sigma_grads
 
std::vector< std::vector< double > > edge_lambda_values
 
std::vector< std::vector< double > > edge_lambda_grads_2d
 
std::vector< std::vector< std::vector< double > > > edge_lambda_grads_3d
 
std::vector< std::vector< std::vector< double > > > edge_lambda_gradgrads_3d
 
std::vector< std::vector< double > > face_lambda_values
 
std::vector< std::vector< double > > face_lambda_grads
 
UpdateFlags update_each
 

Detailed Description

template<int dim, int spacedim = dim>
class FE_NedelecSZ< dim, spacedim >::InternalData

Derived Internal data which is used to store cell-independent data. Note that due to the nature of this element, a number of useful pre-computed quantities are stored for the computation of cell-dependent shape functions.

The main quantities which are stored are associated with edge and face parameterizations. These are:

The definitions of these functionals, as well as the edge and face parameterizations and edge and face extension parameters, can be found on page 82 of Zaglmayr's thesis. The details of the definition of the globally-defined edge and face orientations can be found on page 67.

Definition at line 242 of file fe_nedelec_sz.h.

Member Function Documentation

◆ memory_consumption()

virtual std::size_t FiniteElement< dim, spacedim >::InternalDataBase::memory_consumption ( ) const
virtualinherited

Return an estimate (in bytes) for the memory consumption of this object.

Member Data Documentation

◆ shape_values

template<int dim, int spacedim = dim>
std::vector<std::vector<Tensor<1, dim> > > FE_NedelecSZ< dim, spacedim >::InternalData::shape_values
mutable

Storage for shape functions on the reference element. We only pre-compute cell-based DoFs, as the edge- and face-based DoFs depend on the cell.

Due to the cell-dependent DoFs, this variable is declared mutable.

Definition at line 251 of file fe_nedelec_sz.h.

◆ shape_grads

template<int dim, int spacedim = dim>
std::vector<std::vector<DerivativeForm<1, dim, dim> > > FE_NedelecSZ< dim, spacedim >::InternalData::shape_grads
mutable

Storage for shape function gradients on the reference element. We only pre-compute cell-based DoFs, as the edge- and face-based DoFs depend on the cell.

Due to the cell-dependent DoFs, this variable is declared mutable.

Definition at line 260 of file fe_nedelec_sz.h.

◆ shape_hessians

template<int dim, int spacedim = dim>
std::vector<std::vector<DerivativeForm<2, dim, dim> > > FE_NedelecSZ< dim, spacedim >::InternalData::shape_hessians
mutable

Storage for shape function hessians on the reference element. We only pre-compute cell-based DoFs, as the edge- and face-based DoFs depend on the cell.

Due to the cell-dependent DoFs, this variable is declared mutable.

Definition at line 270 of file fe_nedelec_sz.h.

◆ sigma_imj_values

template<int dim, int spacedim = dim>
std::vector<std::vector<std::vector<double> > > FE_NedelecSZ< dim, spacedim >::InternalData::sigma_imj_values

Storage for all possible edge parameterization between vertices. These are required in the computation of edge- and face-based DoFs, which are cell-dependent.

The edge parameterization of an edge, E, starting at vertex i and ending at vertex \(j\) is given by \(\sigma_{E} = \sigma_{i} - \sigma{j}\).

sigma_imj_values[q][i][j] stores the value of the edge parameterization connected by vertices \(i\) and \(j\) at the q-th quadrature point.

Note that not all of the \(i\) and \(j\) combinations result in valid edges on the hexahedral cell, but they are computed in this fashion for use with non-standard edge and face orientations.

Definition at line 287 of file fe_nedelec_sz.h.

◆ sigma_imj_grads

template<int dim, int spacedim = dim>
std::vector<std::vector<std::vector<double> > > FE_NedelecSZ< dim, spacedim >::InternalData::sigma_imj_grads

Storage for gradients of all possible edge parameterizations between vertices. These are required in the computation of edge- and face-based DoFs, which are cell-dependent. Note that the components of the gradient are constant.

The edge parameterization of an edge, \(E\), starting at vertex \(i\) and ending at vertex \(j\) is given by \(\sigma_{E} = \sigma_{i} - \sigma{j}\).

sigma_imj_grads[i][j][d] stores the gradient of the edge parameterization connected by vertices \(i\) and \(j\) in component \(d\).

Note that the gradient of the edge parameterization is constant on an edge, so we do not need to store it at every quadrature point.

Definition at line 304 of file fe_nedelec_sz.h.

◆ edge_sigma_values

template<int dim, int spacedim = dim>
std::vector<std::vector<double> > FE_NedelecSZ< dim, spacedim >::InternalData::edge_sigma_values

Storage for values of edge parameterizations at quadrature points. These are stored for the 12 edges such that the global vertex numbering would follow the order defined by the "standard" deal.II cell.

edge_sigma_values[m][q] stores the edge parameterization value at the q-th quadrature point on edge m.

These values change with the orientation of the edges of a physical cell and so must take the "sign" into account when used for computation.

Definition at line 317 of file fe_nedelec_sz.h.

◆ edge_sigma_grads

template<int dim, int spacedim = dim>
std::vector<std::vector<double> > FE_NedelecSZ< dim, spacedim >::InternalData::edge_sigma_grads

Storage for gradients of edge parameterization at quadrature points. These are stored for the 12 edges such that the global vertex numbering would follow the order defined by the "standard" deal.II cell.

edge_sigma_grads[m][d] stores the gradient of the edge parameterization for component d on edge m.

These values change with the orientation of the edges of a physical cell and so must take the "sign" into account when used for computation.

Definition at line 330 of file fe_nedelec_sz.h.

◆ edge_lambda_values

template<int dim, int spacedim = dim>
std::vector<std::vector<double> > FE_NedelecSZ< dim, spacedim >::InternalData::edge_lambda_values

Storage for edge extension parameters at quadrature points. These are stored for the 12 edges such that the global vertex numbering would follow the order defined by the "standard" deal.II cell.

The edge extension parameter of an edge, \(E\), starting at vertex \(i\) and ending at vertex \(j\) is given by \(\lambda_{E} = \lambda_{i} + \lambda_{j}\).

Note that under this definition, the values of \(\lambda_{E}\) do not change with the orientation of the edge.

edge_lambda_values[m][q] stores the edge extension parameter value at the \(q\)-th quadrature point on edge \(m\).

Definition at line 347 of file fe_nedelec_sz.h.

◆ edge_lambda_grads_2d

template<int dim, int spacedim = dim>
std::vector<std::vector<double> > FE_NedelecSZ< dim, spacedim >::InternalData::edge_lambda_grads_2d

Storage for gradients of edge extension parameters in 2d. In this case they are constant. These are stored for the 12 edges such that the global vertex numbering* would follow the order defined by the "standard" deal.II cell.

edge_lambda_grads_2d[m][d] stores the gradient of the edge extension parameter for component \(d\) on edge \(m\).

Definition at line 358 of file fe_nedelec_sz.h.

◆ edge_lambda_grads_3d

template<int dim, int spacedim = dim>
std::vector<std::vector<std::vector<double> > > FE_NedelecSZ< dim, spacedim >::InternalData::edge_lambda_grads_3d

Storage for gradients of edge extension parameters in 3d. In this case they are non-constant. These are stored for the 12 edges such that the global vertex numbering* would follow the order defined by the "standard" deal.II cell.

edge_lambda_grads_3d[m][q][d] stores the gradient of the edge extension parameter for component \(d\) at the \(q\)-th quadrature point on edge m.

Definition at line 369 of file fe_nedelec_sz.h.

◆ edge_lambda_gradgrads_3d

template<int dim, int spacedim = dim>
std::vector<std::vector<std::vector<double> > > FE_NedelecSZ< dim, spacedim >::InternalData::edge_lambda_gradgrads_3d

Storage for 2nd derivatives of edge extension parameters in 3d, which are constant across the cell. These are stored for the 12 edges such that the global vertex numbering* would follow the order defined by the "standard" deal.II cell.

edge_lambda_gradgrads_3d[m][d1][d2] stores the 2nd derivatives of the edge extension parameters with respect to components d1 and d2 on edge \(m\).

Definition at line 381 of file fe_nedelec_sz.h.

◆ face_lambda_values

template<int dim, int spacedim = dim>
std::vector<std::vector<double> > FE_NedelecSZ< dim, spacedim >::InternalData::face_lambda_values

Storage for the face extension parameters. These are stored for the 6 faces such that the global vertex numbering would follow the order defined by the "standard" deal.II cell.

The face extension parameter of a face, F, defined by the vertices v1, v2, v3, v4 is given by \(\lambda_{F} = \lambda_{v1} + \lambda_{v2} + \lambda_{v3} + \lambda_{v4}\).

Note that under this definition, the values of \(\lambda_{F}\) do not change with the orientation of the face.

face_lambda_values[m][q] stores the face extension parameter value at the \(q\)-th quadrature point on face \(m\).

Definition at line 399 of file fe_nedelec_sz.h.

◆ face_lambda_grads

template<int dim, int spacedim = dim>
std::vector<std::vector<double> > FE_NedelecSZ< dim, spacedim >::InternalData::face_lambda_grads

Storage for gradients of face extension parameters. These are stored for the 6 faces such that the global vertex numbering would follow the order defined by the "standard" deal.II cell.

face_lambda_grads[m][d] stores the gradient of the face extension parameters for component \(d\) on face \(m\).

Definition at line 409 of file fe_nedelec_sz.h.

◆ update_each

UpdateFlags FiniteElement< dim, spacedim >::InternalDataBase::update_each
inherited

A set of update flags specifying the kind of information that an implementation of the FiniteElement interface needs to compute on each cell or face, i.e., in FiniteElement::fill_fe_values() and friends.

This set of flags is stored here by implementations of FiniteElement::get_data(), FiniteElement::get_face_data(), or FiniteElement::get_subface_data(), and is that subset of the update flags passed to those functions that require re-computation on every cell. (The subset of the flags corresponding to information that can be computed once and for all already at the time of the call to FiniteElement::get_data() – or an implementation of that interface – need not be stored here because it has already been taken care of.)

Definition at line 720 of file fe.h.


The documentation for this class was generated from the following file: