Reference documentation for deal.II version GIT relicensing-493-g2fa0c96cf2 2024-04-29 17:00:02+00:00
\(\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 | Private Attributes | List of all members
FECouplingValues< dim1, dim2, spacedim > Class Template Reference

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

Public Member Functions

 FECouplingValues ()
 
 FECouplingValues (const FEValuesBase< dim1, spacedim > &fe_values_1, const FEValuesBase< dim2, spacedim > &fe_values_2, const DoFCouplingType &dof_coupling_type=DoFCouplingType::independent, const QuadratureCouplingType &quadrature_coupling_type=QuadratureCouplingType::tensor_product)
 
void reinit (const FEValuesBase< dim1, spacedim > &fe_values_1, const FEValuesBase< dim2, spacedim > &fe_values_2, const DoFCouplingType &dof_coupling_type=DoFCouplingType::independent, const QuadratureCouplingType &quadrature_coupling_type=QuadratureCouplingType::tensor_product)
 
template<typename Extractor >
FEValuesExtractors::FirstCoupling< Extractor > get_first_extractor (const Extractor &extractor) const
 
template<typename Extractor >
FEValuesExtractors::SecondCoupling< Extractor > get_second_extractor (const Extractor &extractor) const
 
double JxW (const unsigned int quadrature_point) const
 
std::pair< Point< spacedim >, Point< spacedim > > quadrature_point (const unsigned int quadrature_point) const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intquadrature_point_indices () const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intfirst_dof_indices () const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intsecond_dof_indices () const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intcoupling_dof_indices () const
 
std::vector< types::global_dof_indexget_coupling_dof_indices (const std::vector< types::global_dof_index > &dof_indices_1, const std::vector< types::global_dof_index > &dof_indices_2, const types::global_dof_index dofs_offset_1=0, const types::global_dof_index dofs_offset_2=0) const
 
std::pair< unsigned int, unsigned intcoupling_dof_to_dof_indices (const unsigned int coupling_dof_index) const
 
std::pair< unsigned int, unsigned intcoupling_quadrature_to_quadrature_indices (const unsigned int quadrature_point) const
 
unsigned int n_coupling_dofs () const
 
unsigned int n_first_dofs () const
 
unsigned int n_second_dofs () const
 
unsigned int n_quadrature_points () const
 
Extractors Methods to extract individual components
template<typename Extractor >
const FEValuesViews::RenumberedView< FEValuesViews::View< dim1, spacedim, Extractor > > operator[] (const FEValuesExtractors::FirstCoupling< Extractor > &extractor) const
 
template<typename Extractor >
const FEValuesViews::RenumberedView< FEValuesViews::View< dim2, spacedim, Extractor > > operator[] (const FEValuesExtractors::SecondCoupling< Extractor > &extractor) const
 

Private Attributes

DoFCouplingType dof_coupling_type
 
QuadratureCouplingType quadrature_coupling_type
 
SmartPointer< const FEValuesBase< dim1, spacedim > > first_fe_values
 
SmartPointer< const FEValuesBase< dim2, spacedim > > second_fe_values
 
std::unique_ptr< const FEValuesViews::RenumberingDatafirst_renumbering_data
 
std::unique_ptr< const FEValuesViews::RenumberingDatasecond_renumbering_data
 
unsigned int n_quadrature_points_
 
unsigned int n_coupling_dofs_
 

Detailed Description

template<int dim1, int dim2 = dim1, int spacedim = dim1>
class FECouplingValues< dim1, dim2, spacedim >

FECouplingValues is a class that facilitates the integration of finite element data between two different finite element objects, possibly living on different grids, and with possibly different topological dimensions (i.e., cells, faces, edges, and any combination thereof).

This class provides a way to simplify the implementation of the following abstract operation:

\[ \int_{T_1} \int{T_2} K(x_1, x_2) \phi^1_i(x_1) \phi^2_j(x_2) dT_1 dT_2 \]

for three different types of Kernels \(K\):

For the first case, one may think that the only natural way to proceed is to compute the double integral by simply nesting two loops:

\[ \int_{T_1} \int{T_2} K(x_1, x_2) \phi^1_i(x_1) \phi^2_j(x_2) dT_1 dT_2 \approx \sum_{q_1} \sum_{q_2} K(x_1^{q_1}, x_2^{q_2}) \phi^1_i(x_1^{q_1}) \phi^2_j(x_2^{q_2}) w_1^{q_1} w_2^{q_2}, \]

where \(x_1^{q_1}\) and \(x_2^{q_2}\) are the quadrature points in \(T_1\) and \(T_2\) respectively, and \(w_1^{q_1}\) and \(w_2^{q_2}\) are the corresponding quadrature weights.

This, however is not the only way to proceed. In fact, such an integral can be rewritten as a single loop over corresponding elements of two arrays of points with the same length that can be thought of as a single quadrature rule on the set \(T_1\times T_2\). For singular kernels, for example, this is often the only way to proceed, since the quadrature formula on \(T_1\times T_2\) is usually not written as a tensor product quadrature formula, and one needs to build a custom quadrature formula for this purpose.

This class allows one to treat the three cases above in the same way, and to approximate the integral as follows:

\[ \int_{T_1} \int{T_2} K(x_1, x_2) \phi^1_i(x_1) \phi^2_j(x_2) dT_1 dT_2 \approx \sum_{i=1}^{N_q} K(x_1^{i}, x_2^{i}) \phi^1_i(x_1^{i}) \phi^2_j(x_2^{i}) w_1^{i} w_2^i, \]

Since the triple of objects \((\{q\}, \{w\}, \{\phi\})\) is usually provided by a class derived from the FEValuesBase class, this is the type that the class needs at construction time. \(T_1\) and \(T_2\) can be two arbitrary cells, faces, or edges belonging to possibly different meshes (or to meshes with different topological dimensions), \(\phi^1_i\) and \(\phi^2_j\) are basis functions defined on \(T_1\) and \(T_2\), respectively.

The case of the Dirac distribution is when \(T_1\) and \(T_2\) correspond to the common face of two neighboring cells. In this case, this class provides a functionality which is similar to the FEInterfaceValues class, and gives you a way to access values of basis functions on the neighboring cells, as well as their gradients and Hessians, in a unified fashion, on the face.

Similarly, this class can be used to couple bulk and surface meshes across the faces of the bulk mesh. In this case, the two FEValuesBase objects will have different topological dimension (i.e., one will be a cell in a co-dimension one triangulation, and the other a face of a bulk grid with co-dimension zero), and the QuadratureCouplingType argument is usually chosen to be QuadratureCouplingType::reorder, since the quadrature points of the two different FEValuesBase objects are not necessarily generated with the same ordering.

The type of integral to compute is controlled by the QuadratureCouplingType argument (see the documentation of that enum class for more details), while the type degrees of freedom coupling is controlled by the DoFCouplingType argument (see the documentation of that enum class for more details).

As an example usage of this class, consider the a bilinear form of the form:

\[ \int_{T_1} \int{T_2} K_1(x_1, x_2) v_i(x_1) u_j(x_2) dT_1 dT_2 + \int_{T_1} \int{T_2} K_2(x_1, x_2) p_i(x_1) q_j(x_2) dT_1 dT_2 \]

where the finite dimensional space has two scalar components. We indicate with \(v_i\) and \(p_i\) the trial functions, and with \(u_j\) and \(q_j\) the corresponding test functions. \(K_1\) and \(K_2\) are coupling kernels: such a formulation is used, for example, to write the bilinear forms of Galerkin boundary element methods.

The corresponding implementation would look like the following:

... // double loop over cells that yields cell_1 and cell_2
fe_values_1.reinit(cell_1);
fe_values_2.reinit(cell_2);
CouplingFEValues<dim, dim, spacedim> cfv(fe_values1, fe_values2,
FullMatrix<double> local_matrix(fe_values1.dofs_per_cell,
fe_values2.dofs_per_cell);
// Extractor on first cell
const auto v1 = cfv.get_first_extractor(FEValuesExtractor::Scalar(0));
const auto p1 = cfv.get_first_extractor(FEValuesExtractor::Scalar(1));
// Extractor on second cell
const auto u2 = cfv.get_second_extractor(FEValuesExtractor::Scalar(0));
const auto q2 = cfv.get_second_extractor(FEValuesExtractor::Scalar(1));
...
for (const unsigned int q : cfv.quadrature_point_indices()) {
const auto &[x_q,y_q] = cfv.quadrature_point(q);
for (const unsigned int i : cfv.first_dof_indices())
{
const auto &v_i = cfv[v1].value(i, q);
const auto &p_i = cfv[p1].value(i, q);
for (const unsigned int j : cfv.second_dof_indices())
{
const auto &u_j = cfv[u2].value(j, q);
const auto &q_j = cfv[q2].value(j, q);
local_matrix(i, j) += (K1(x_q, y_q) * v_i * u_j +
K2(x_q, y_q) * p_i * q_j) *
cfv[0].JxW(q) *
cfv[1].JxW(q);
}
}
}
std_cxx20::ranges::iota_view< unsigned int, unsigned int > first_dof_indices() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
double JxW(const unsigned int quadrature_point) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > second_dof_indices() const
const unsigned int v1

In the above loop, the quadrature points of the two FEValuesBase objects are grouped together in a single loop, while the dof indices are kept separate.

Internally, this class provides an abstraction to organize coupling terms between two arbitrary FEValuesBase objects, and provides a unified way to access the values of the two basis, and of the two sets of quadrature points and weights. The underlying data is stored in the two FEValuesBase objects, and this class provides access to two RenumberedView objects, which can be obtained by calling the operator[] method of this class with the arguments constructed by calling the get_first_extractor() and get_second_extractor() methods of this class with a standard FEValuesExtractors object, i.e.,

// extract the first scalar component of the basis
FEValuesExtractor scalar(0);
...
FECouplingValues<dim> fe_coupling(fev1, fev1);
// Extractors for the two FEValuesBase objects are returned by the
// get_first_extractor() and
// get_second_extractor() methods
const auto first = fe_coupling.get_first_extractor(scalar);
const auto second = fe_coupling.get_second_extractor(scalar);
// Now we can use the augmented extractors to access the values of the two
// FEValuesBase objects
const auto & first_vi = fe_coupling[first].value(i, q);
const auto & second_vi = fe_coupling[second].value(i, q);
Point< 2 > second
Definition grid_out.cc:4624
Point< 2 > first
Definition grid_out.cc:4623

Definition at line 731 of file fe_coupling_values.h.

Constructor & Destructor Documentation

◆ FECouplingValues() [1/2]

template<int dim1, int dim2 = dim1, int spacedim = dim1>
FECouplingValues< dim1, dim2, spacedim >::FECouplingValues ( )

Construct the FECouplingValues in an invalid state. Before you can use this object, you must call the reinit() function.

◆ FECouplingValues() [2/2]

template<int dim1, int dim2 = dim1, int spacedim = dim1>
FECouplingValues< dim1, dim2, spacedim >::FECouplingValues ( const FEValuesBase< dim1, spacedim > &  fe_values_1,
const FEValuesBase< dim2, spacedim > &  fe_values_2,
const DoFCouplingType dof_coupling_type = DoFCouplingType::independent,
const QuadratureCouplingType quadrature_coupling_type = QuadratureCouplingType::tensor_product 
)

Construct the FECouplingValues with two arbitrary FEValuesBase objects. This class assumes that the FEValuesBase objects that are given at construction time are initialized and ready to use (i.e., that you have called the reinit() function on them before calling this constructor).

Notice that the actual renumbering of the degrees of freedom and quadrature points is done at construction time, or upon calling the reinit() function. If you change the underlying FEValuesBase objects after construction, you must call the reinit() function to update the renumbering. This may or may not be necessary, depending on the type of coupling that you are using.

This really depends on the application and on the specific type of coupling. For example, for volume/volume coupling, i.e., for operators with non-local and non-singular kernels of type

\[ \int_K \int_T f(\phi_i(x)-\phi_j(y), x-y) dx dy \]

you may initialize FECouplingValues once, and just reinit the underlying FEValuesBase objects on different cells K and T, without the need to recompute the coupling (i.e., the numbering is always the same, and nothing differs from what happened in the first call).

For cell/surface coupling, the same cell may couple with different faces, so the renumbering must be really computed from scratch for each pair of FEValuesBase objects, so reinitializing the underlying cells and faces will make the renumbering itself invalid, and FECouplingValues must be reinitialized (o constructed from scratch) after calling fe_values_1.reinit(K) and fe_values_1.reinit(T).

Parameters
fe_values_1The first FEValuesBase object.
fe_values_2The second FEValuesBase object.
dof_coupling_typeThe type of dof coupling to use.
quadrature_coupling_typeThe type of quadrature to use for the coupling.

Member Function Documentation

◆ reinit()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
void FECouplingValues< dim1, dim2, spacedim >::reinit ( const FEValuesBase< dim1, spacedim > &  fe_values_1,
const FEValuesBase< dim2, spacedim > &  fe_values_2,
const DoFCouplingType dof_coupling_type = DoFCouplingType::independent,
const QuadratureCouplingType quadrature_coupling_type = QuadratureCouplingType::tensor_product 
)

Reinitialize the FECouplingValues with two arbitrary FEValuesBase objects. The FEValuesBase objects must be initialized and ready to use, i.e., you must have called the reinit() function on them before calling this method.

This method computes the actual renumbering of the degrees of freedom and quadrature points. If you change the underlying FEValuesBase objects after calling this method, you mey need to call the reinit() function to update the renumbering. This may or may not be necessary, depending on the type of coupling that you are using.

This really depends on the application and on the specific type of coupling. For example, for volume/volume coupling, i.e., for operators with non-local and non-singular kernels of type

\[ \int_K \int_T f(\phi_i(x)-\phi_j(y), x-y) dx dy \]

you may initialize FECouplingValues once, and just reinit the underlying FEValuesBase objects on different cells K and T, without the need to recompute the coupling (i.e., the numbering is always the same, and nothing differs from what happened in the first call).

For cell/surface coupling, the same cell may couple with different faces, so the renumbering must be really computed from scratch for each pair of FEValuesBase objects, so reinitializing the underlying cells and faces will make the renumbering itself invalid, and FECouplingValues must be reinitialized (o constructed from scratch) after calling fe_values_1.reinit(K) and fe_values_1.reinit(T).

Parameters
fe_values_1The first FEValuesBase object.
fe_values_2The second FEValuesBase object.
dof_coupling_typeThe type of dof coupling to use.
quadrature_coupling_typeThe type of quadrature to use for the coupling.

◆ get_first_extractor()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
template<typename Extractor >
FEValuesExtractors::FirstCoupling< Extractor > FECouplingValues< dim1, dim2, spacedim >::get_first_extractor ( const Extractor &  extractor) const

Return a FirstCoupling object that can be used to extract values from the first FEValuesBase object.

◆ get_second_extractor()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
template<typename Extractor >
FEValuesExtractors::SecondCoupling< Extractor > FECouplingValues< dim1, dim2, spacedim >::get_second_extractor ( const Extractor &  extractor) const

Return a SecondCoupling object that can be used to extract values from the second FEValuesBase object.

◆ JxW()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
double FECouplingValues< dim1, dim2, spacedim >::JxW ( const unsigned int  quadrature_point) const

Return the value of JxW at the given quadrature point.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_JxW_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ quadrature_point()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
std::pair< Point< spacedim >, Point< spacedim > > FECouplingValues< dim1, dim2, spacedim >::quadrature_point ( const unsigned int  quadrature_point) const

Return the two quadrature points in real space at the given quadrature point index, corresponding to a quadrature point in the set \(T_1\times T_2\).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_quadrature_points flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ quadrature_point_indices()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
std_cxx20::ranges::iota_view< unsigned int, unsigned int > FECouplingValues< dim1, dim2, spacedim >::quadrature_point_indices ( ) const

Return an object that can be thought of as an array containing all indices from zero to n_quadrature_points. This allows to write code using range-based for loops.

See also
deal.II and Modern C++ standards

◆ first_dof_indices()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
std_cxx20::ranges::iota_view< unsigned int, unsigned int > FECouplingValues< dim1, dim2, spacedim >::first_dof_indices ( ) const

Return an object that can be thought of as an array containing all indices from zero (inclusive) to n_first_dofs() (exclusive). This allows one to write code using range-based for loops.

◆ second_dof_indices()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
std_cxx20::ranges::iota_view< unsigned int, unsigned int > FECouplingValues< dim1, dim2, spacedim >::second_dof_indices ( ) const

Return an object that can be thought of as an array containing all indices from zero (inclusive) to n_second_dofs() (exclusive). This allows one to write code using range-based for loops.

◆ coupling_dof_indices()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
std_cxx20::ranges::iota_view< unsigned int, unsigned int > FECouplingValues< dim1, dim2, spacedim >::coupling_dof_indices ( ) const

Return an object that can be thought of as an array containing all indices from zero (inclusive) to n_coupling_dof_indices() (exclusive). This allows one to write code using range-based for loops.

This function only makes sense when the dof coupling type is DoFCouplingType::contiguous, and, in this case, n_coupling_dof_indices() is the sum of n_first_dofs() and n_second_dofs(). If DoFCouplingType::independent is used, you should call first_dof_indices() or second_dof_indices() instead, and calling this function will throw an exception.

◆ get_coupling_dof_indices()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
std::vector< types::global_dof_index > FECouplingValues< dim1, dim2, spacedim >::get_coupling_dof_indices ( const std::vector< types::global_dof_index > &  dof_indices_1,
const std::vector< types::global_dof_index > &  dof_indices_2,
const types::global_dof_index  dofs_offset_1 = 0,
const types::global_dof_index  dofs_offset_2 = 0 
) const

Return the set of joint DoF indices, possibly offset by the given values.

For coupled problems, the DoF indices of the two FEValuesBase objects are generally constructed from different DoFHandler objects. According to how the user wants to construct the coupling matrices, the two dof indices may map local dof indices to global dof indices starting from zero (i.e., we are assembling a single rectangular matrix, and the first dof indices run on rows of the matrix, while the second dof indices run on columns of the matrix), or the user may be wanting to map the coupling matrix to a block of a larger block matrix. In this second case the dof indices of either or both blocks may need to be shifted by a fixed amount, to fall within the right range of the target block matrix. This is usually required if the user wants to keep their solution vectors in a single block vector, where each block corresponds to a different DoFHandler object. In this case, all dofs referring to the first FEValuesBase object are offset by the value of dofs_offset_1, while all dofs referring to the second FEValuesBase object are offset by the value of dofs_offset_2, i.e., the following relationship is satisfied:

fe_values_1.get_dof_indices(dof_indices_1);
fe_values_2.get_dof_indices(dof_indices_2);
auto coupling_dof_indices = fe_coupling_values.
get_coupling_dof_indices(dof_indices_1, dof_indices_2,
dofs_offset_1, dofs_offset_2);
for(const auto &i : fe_coupling_values.coupling_dof_indices()) {
auto [id1, id2] = fe_coupling_values.coupling_dof_to_dof_indices(i);
// always true
coupling_dof_indices[i] == dof_indices_1[id1] + dofs_offset_1;
// always true
coupling_dof_indices[i] == dof_indices_2[id2] + dofs_offset_2;
}
std::vector< types::global_dof_index > get_coupling_dof_indices(const std::vector< types::global_dof_index > &dof_indices_1, const std::vector< types::global_dof_index > &dof_indices_2, const types::global_dof_index dofs_offset_1=0, const types::global_dof_index dofs_offset_2=0) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > coupling_dof_indices() const
static const unsigned int invalid_unsigned_int
Definition types.h:220

This method only makes sense when the dof coupling type is DoFCouplingType::contiguous. If DoFCouplingType::independent is used, calling this function will result in an exception.

◆ coupling_dof_to_dof_indices()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
std::pair< unsigned int, unsigned int > FECouplingValues< dim1, dim2, spacedim >::coupling_dof_to_dof_indices ( const unsigned int  coupling_dof_index) const

Convert a coupling dof index into the corresponding local DoF indices of the two FEValuesObjects. If a DoF is only active on one of the FEValuesObjects, the other index will be numbers::invalid_unsigned_int.

◆ coupling_quadrature_to_quadrature_indices()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
std::pair< unsigned int, unsigned int > FECouplingValues< dim1, dim2, spacedim >::coupling_quadrature_to_quadrature_indices ( const unsigned int  quadrature_point) const

Convert a quadrature index into the corresponding local quadrature indices of the two FEValuesObjects. Both indices are guaranteed to be valid within the corresponding FEValuesObject.

◆ operator[]() [1/2]

template<int dim1, int dim2 = dim1, int spacedim = dim1>
template<typename Extractor >
const FEValuesViews::RenumberedView< FEValuesViews::View< dim1, spacedim, Extractor > > FECouplingValues< dim1, dim2, spacedim >::operator[] ( const FEValuesExtractors::FirstCoupling< Extractor > &  extractor) const

Create a combined view of the first FECouplingValues object that represents a view of the possibly vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews.

◆ operator[]() [2/2]

template<int dim1, int dim2 = dim1, int spacedim = dim1>
template<typename Extractor >
const FEValuesViews::RenumberedView< FEValuesViews::View< dim2, spacedim, Extractor > > FECouplingValues< dim1, dim2, spacedim >::operator[] ( const FEValuesExtractors::SecondCoupling< Extractor > &  extractor) const

Create a combined view of the second FECouplingValues object that represents a view of the possibly vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews.

◆ n_coupling_dofs()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
unsigned int FECouplingValues< dim1, dim2, spacedim >::n_coupling_dofs ( ) const

Return the number of coupling DoF indices. If DoFCouplingType::independent is used, this function returns numbers::invalid_unsigned_int, otherwise it returns the sum of n_first_dofs() and n_second_dofs().

◆ n_first_dofs()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
unsigned int FECouplingValues< dim1, dim2, spacedim >::n_first_dofs ( ) const

Return the number of first DoF indices. This generally coincides with n_dofs_per_cell of the first FEValuesBase object.

◆ n_second_dofs()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
unsigned int FECouplingValues< dim1, dim2, spacedim >::n_second_dofs ( ) const

Return the number of second DoF indices. This generally coincides with n_dofs_per_cell of the second FEValuesBase object.

◆ n_quadrature_points()

template<int dim1, int dim2 = dim1, int spacedim = dim1>
unsigned int FECouplingValues< dim1, dim2, spacedim >::n_quadrature_points ( ) const

Return the number of quadrature points.

Member Data Documentation

◆ dof_coupling_type

template<int dim1, int dim2 = dim1, int spacedim = dim1>
DoFCouplingType FECouplingValues< dim1, dim2, spacedim >::dof_coupling_type
private

The dof coupling type used by this object.

Definition at line 1027 of file fe_coupling_values.h.

◆ quadrature_coupling_type

template<int dim1, int dim2 = dim1, int spacedim = dim1>
QuadratureCouplingType FECouplingValues< dim1, dim2, spacedim >::quadrature_coupling_type
private

The quadrature coupling type used by this object.

Definition at line 1032 of file fe_coupling_values.h.

◆ first_fe_values

template<int dim1, int dim2 = dim1, int spacedim = dim1>
SmartPointer<const FEValuesBase<dim1, spacedim> > FECouplingValues< dim1, dim2, spacedim >::first_fe_values
private

Pointer to first FEValuesBase object.

Definition at line 1037 of file fe_coupling_values.h.

◆ second_fe_values

template<int dim1, int dim2 = dim1, int spacedim = dim1>
SmartPointer<const FEValuesBase<dim2, spacedim> > FECouplingValues< dim1, dim2, spacedim >::second_fe_values
private

Pointer to second FEValuesBase object.

Definition at line 1042 of file fe_coupling_values.h.

◆ first_renumbering_data

template<int dim1, int dim2 = dim1, int spacedim = dim1>
std::unique_ptr<const FEValuesViews::RenumberingData> FECouplingValues< dim1, dim2, spacedim >::first_renumbering_data
private

Renumbering data for the first FEValuesBase object.

Definition at line 1047 of file fe_coupling_values.h.

◆ second_renumbering_data

template<int dim1, int dim2 = dim1, int spacedim = dim1>
std::unique_ptr<const FEValuesViews::RenumberingData> FECouplingValues< dim1, dim2, spacedim >::second_renumbering_data
private

Renumbering data for the second FEValuesBase object.

Definition at line 1052 of file fe_coupling_values.h.

◆ n_quadrature_points_

template<int dim1, int dim2 = dim1, int spacedim = dim1>
unsigned int FECouplingValues< dim1, dim2, spacedim >::n_quadrature_points_
private

Number of quadrature points.

Definition at line 1057 of file fe_coupling_values.h.

◆ n_coupling_dofs_

template<int dim1, int dim2 = dim1, int spacedim = dim1>
unsigned int FECouplingValues< dim1, dim2, spacedim >::n_coupling_dofs_
private

Number of coupling DoF indices. If DoFCouplingType::independent is used, this is numbers::invalid_unsigned_int, while it is n_first_dofs() + n_second_dofs() in the the case of DoFCouplingType::contiguous.

Definition at line 1064 of file fe_coupling_values.h.


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