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 {
29 std::vector<Eigen::MatrixXd>
31 
37 std::vector<Eigen::MatrixXd>
39 
45 std::vector<Eigen::MatrixXd>
47 
53 std::vector<Eigen::MatrixXd>
55 
70 std::pair<Eigen::ArrayXXd, Eigen::MatrixXd>
71 make_integral_moments(const FiniteElement& moment_space,
72  const cell::type celltype, const int value_size,
73  const int poly_deg, const int q_deg);
74 
89 std::pair<Eigen::ArrayXXd, Eigen::MatrixXd>
90 make_dot_integral_moments(const FiniteElement& moment_space,
91  const cell::type celltype, const int value_size,
92  const int poly_deg, const int q_deg);
93 
106 std::pair<Eigen::ArrayXXd, Eigen::MatrixXd>
107 make_tangent_integral_moments(const FiniteElement& moment_space,
108  const cell::type celltype, const int value_size,
109  const int poly_deg, const int q_deg);
110 
123 std::pair<Eigen::ArrayXXd, Eigen::MatrixXd>
124 make_normal_integral_moments(const FiniteElement& moment_space,
125  const cell::type celltype, const int value_size,
126  const int poly_deg, const int q_deg);
127 
128 }; // namespace moments
129 } // namespace basix
Definition: finite-element.h:181
type
Cell type.
Definition: cell.h:21
std::pair< Eigen::ArrayXXd, 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:505
std::vector< Eigen::MatrixXd > create_normal_moment_dof_transformations(const FiniteElement &moment_space)
Definition: moments.cpp:257
std::pair< Eigen::ArrayXXd, 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:437
std::vector< Eigen::MatrixXd > create_dot_moment_dof_transformations(const FiniteElement &moment_space)
Definition: moments.cpp:39
std::vector< Eigen::MatrixXd > create_tangent_moment_dof_transformations(const FiniteElement &moment_space)
Definition: moments.cpp:270
std::pair< Eigen::ArrayXXd, 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:363
std::pair< Eigen::ArrayXXd, 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:284
std::vector< Eigen::MatrixXd > create_moment_dof_transformations(const FiniteElement &moment_space)
Definition: moments.cpp:211
Placeholder.
Definition: basix.h:10