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
|
// Copyright (c) 2020 Chris Richardson & Matthew Scroggs
// FEniCS Project
// SPDX-License-Identifier: MIT
#pragma once
#include "cell.h"
#include <Eigen/Dense>
namespace basix
{
class FiniteElement;
/// ## Integral moments
/// These functions generate dual set matrices for integral moments
/// against spaces on a subentity of the cell
namespace moments
{
/// Make simple integral moments
///
/// This will compute the integral of each function in the moment space
/// over each sub entity of the moment space's cell type in a cell with
/// the given type. For example, if the input cell type is a triangle,
/// and the moment space is a P1 space on an edge, this will perform two
/// integrals for each of the 3 edges of the triangle.
///
/// @param moment_space The space to compute the integral moments against
/// @param celltype The cell type of the cell on which the space is being
/// defined
/// @param value_size The value size of the space being defined
/// @param poly_deg The polynomial degree of the poly set that defines the space
/// @param q_deg The quadrature degree used for the integrals
Eigen::MatrixXd make_integral_moments(const FiniteElement& moment_space,
const cell::type celltype,
const int value_size, const int poly_deg,
const int q_deg);
/// Make interpolation points and weights for simple integral moments
///
/// @param moment_space The space to compute the integral moments against
/// @param celltype The cell type of the cell on which the space is being
/// defined
/// @param value_size The value size of the space being defined
/// @param poly_deg The polynomial degree of the poly set that defines the space
/// @param q_deg The quadrature degree used for the integrals
std::pair<Eigen::ArrayXXd, Eigen::MatrixXd> make_integral_moments_interpolation(
const FiniteElement& moment_space, const cell::type celltype,
const int value_size, const int poly_deg, const int q_deg);
/// Make dot product integral moments
///
/// This will compute the integral of each function in the moment space over
/// each sub entity of the moment space's cell type in a cell with the given
/// type. For example, if the input cell type is a triangle, and the moment
/// space is a P1 space on an edge, this will perform two integrals for each of
/// the 3 edges of the triangle.
///
/// @param moment_space The space to compute the integral moments against
/// @param celltype The cell type of the cell on which the space is being
/// defined
/// @param value_size The value size of the space being defined
/// @param poly_deg The polynomial degree of the poly set that defines the space
/// @param q_deg The quadrature degree used for the integrals
Eigen::MatrixXd make_dot_integral_moments(const FiniteElement& moment_space,
const cell::type celltype,
const int value_size,
const int poly_deg, const int q_deg);
/// Make interpolation points and weights for dot product integral moments
///
/// @param moment_space The space to compute the integral moments against
/// @param celltype The cell type of the cell on which the space is being
/// defined
/// @param value_size The value size of the space being defined
/// @param poly_deg The polynomial degree of the poly set that defines the space
/// @param q_deg The quadrature degree used for the integrals
std::pair<Eigen::ArrayXXd, Eigen::MatrixXd>
make_dot_integral_moments_interpolation(const FiniteElement& moment_space,
const cell::type celltype,
const int value_size,
const int poly_deg, const int q_deg);
/// Make tangential integral moments
///
/// These can only be used when the moment space is defined on edges of
/// the cell
///
/// @param moment_space The space to compute the integral moments against
/// @param celltype The cell type of the cell on which the space is
/// being defined
/// @param value_size The value size of the space being defined
/// @param poly_deg The polynomial degree of the poly set that defines
/// the space
/// @param q_deg The quadrature degree used for the integrals
Eigen::MatrixXd make_tangent_integral_moments(const FiniteElement& moment_space,
const cell::type celltype,
const int value_size,
const int poly_deg,
const int q_deg);
/// Make interpolation points and weights for tangent integral moments
///
/// These can only be used when the moment space is defined on edges of
/// the cell
///
/// @param moment_space The space to compute the integral moments against
/// @param celltype The cell type of the cell on which the space is
/// being defined
/// @param value_size The value size of the space being defined
/// @param poly_deg The polynomial degree of the poly set that defines
/// the space
/// @param q_deg The quadrature degree used for the integrals
std::pair<Eigen::ArrayXXd, Eigen::MatrixXd>
make_tangent_integral_moments_interpolation(const FiniteElement& moment_space,
const cell::type celltype,
const int value_size,
const int poly_deg,
const int q_deg);
/// Make normal integral moments
///
/// These can only be used when the moment space is defined on facets of
/// the cell
///
/// @param moment_space The space to compute the integral moments against
/// @param celltype The cell type of the cell on which the space is
/// being defined
/// @param value_size The value size of the space being defined
/// @param poly_deg The polynomial degree of the poly set that defines
/// the space
/// @param q_deg The quadrature degree used for the integrals
Eigen::MatrixXd make_normal_integral_moments(const FiniteElement& moment_space,
const cell::type celltype,
const int value_size,
const int poly_deg,
const int q_deg);
/// Make interpolation points and weights for normal integral moments
///
/// These can only be used when the moment space is defined on facets of
/// the cell
///
/// @param moment_space The space to compute the integral moments against
/// @param celltype The cell type of the cell on which the space is
/// being defined
/// @param value_size The value size of the space being defined
/// @param poly_deg The polynomial degree of the poly set that defines
/// the space
/// @param q_deg The quadrature degree used for the integrals
std::pair<Eigen::ArrayXXd, Eigen::MatrixXd>
make_normal_integral_moments_interpolation(const FiniteElement& moment_space,
const cell::type celltype,
const int value_size,
const int poly_deg, const int q_deg);
}; // namespace moments
} // namespace basix
|