Tabulation Project basix
moments.h
1 // Copyright (c) 2020 Chris Richardson & Matthew Scroggs
2 // FEniCS Project
3 // SPDX-License-Identifier: MIT
4 
5 #pragma once
6 
7 #include "cell.h"
8 #include <Eigen/Dense>
9 
10 namespace basix
11 {
12 
13 class FiniteElement;
14 
18 namespace moments
19 {
34 Eigen::MatrixXd make_integral_moments(const FiniteElement& moment_space,
35  const cell::type celltype,
36  const int value_size, const int poly_deg,
37  const int q_deg);
38 
47 std::pair<Eigen::ArrayXXd, Eigen::MatrixXd> make_integral_moments_interpolation(
48  const FiniteElement& moment_space, const cell::type celltype,
49  const int value_size, const int poly_deg, const int q_deg);
50 
65 Eigen::MatrixXd make_dot_integral_moments(const FiniteElement& moment_space,
66  const cell::type celltype,
67  const int value_size,
68  const int poly_deg, const int q_deg);
69 
78 std::pair<Eigen::ArrayXXd, Eigen::MatrixXd>
80  const cell::type celltype,
81  const int value_size,
82  const int poly_deg, const int q_deg);
83 
96 Eigen::MatrixXd make_tangent_integral_moments(const FiniteElement& moment_space,
97  const cell::type celltype,
98  const int value_size,
99  const int poly_deg,
100  const int q_deg);
101 
114 std::pair<Eigen::ArrayXXd, Eigen::MatrixXd>
116  const cell::type celltype,
117  const int value_size,
118  const int poly_deg,
119  const int q_deg);
120 
133 Eigen::MatrixXd make_normal_integral_moments(const FiniteElement& moment_space,
134  const cell::type celltype,
135  const int value_size,
136  const int poly_deg,
137  const int q_deg);
138 
151 std::pair<Eigen::ArrayXXd, Eigen::MatrixXd>
153  const cell::type celltype,
154  const int value_size,
155  const int poly_deg, const int q_deg);
156 
157 }; // namespace moments
158 } // namespace basix
Definition: finite-element.h:148
type
Cell type.
Definition: cell.h:21
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)
Definition: moments.cpp:175
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)
Definition: moments.cpp:34
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)
Definition: moments.cpp:248
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)
Definition: moments.cpp:379
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)
Definition: moments.cpp:315
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)
Definition: moments.cpp:530
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)
Definition: moments.cpp:106
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)
Definition: moments.cpp:442
basix
Definition: basix.h:7
int value_size(int handle)
Value size.
Definition: basix.cpp:56