Basix 0.5.1

Home     Installation     Demos     C++ docs     Python docs

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 <array>
9 #include <tuple>
10 #include <utility>
11 #include <vector>
12 
13 namespace basix
14 {
15 
16 class FiniteElement;
17 
19 namespace moments
20 {
21 
41 std::tuple<std::vector<std::vector<double>>, std::array<std::size_t, 2>,
42  std::vector<std::vector<double>>, std::array<std::size_t, 4>>
43 make_integral_moments(const FiniteElement& moment_space, cell::type celltype,
44  std::size_t value_size, int q_deg);
45 
68 std::tuple<std::vector<std::vector<double>>, std::array<std::size_t, 2>,
69  std::vector<std::vector<double>>, std::array<std::size_t, 4>>
71  std::size_t value_size, int q_deg);
72 
90 std::tuple<std::vector<std::vector<double>>, std::array<std::size_t, 2>,
91  std::vector<std::vector<double>>, std::array<std::size_t, 4>>
93  std::size_t value_size, int q_deg);
94 
111 std::tuple<std::vector<std::vector<double>>, std::array<std::size_t, 2>,
112  std::vector<std::vector<double>>, std::array<std::size_t, 4>>
114  std::size_t value_size, int q_deg);
115 
116 } // namespace moments
117 } // namespace basix
basix::cell::type
type
Cell type.
Definition: cell.h:19
basix::moments::make_normal_integral_moments
std::tuple< std::vector< std::vector< double > >, std::array< std::size_t, 2 >, std::vector< std::vector< double > >, std::array< std::size_t, 4 > > make_normal_integral_moments(const FiniteElement &V, cell::type celltype, std::size_t value_size, int q_deg)
Compute interpolation points and weights for normal integral moments.
Definition: moments.cpp:333
basix::moments::make_tangent_integral_moments
std::tuple< std::vector< std::vector< double > >, std::array< std::size_t, 2 >, std::vector< std::vector< double > >, std::array< std::size_t, 4 > > make_tangent_integral_moments(const FiniteElement &V, cell::type celltype, std::size_t value_size, int q_deg)
Make interpolation points and weights for tangent integral moments.
Definition: moments.cpp:262
basix::FiniteElement
A finite element.
Definition: finite-element.h:143
basix::moments::make_integral_moments
std::tuple< std::vector< std::vector< double > >, std::array< std::size_t, 2 >, std::vector< std::vector< double > >, std::array< std::size_t, 4 > > make_integral_moments(const FiniteElement &moment_space, cell::type celltype, std::size_t value_size, int q_deg)
Make interpolation points and weights for simple integral moments.
Definition: moments.cpp:103
basix
Basix: FEniCS runtime basis evaluation library.
Definition: cell.h:15
basix::moments::make_dot_integral_moments
std::tuple< std::vector< std::vector< double > >, std::array< std::size_t, 2 >, std::vector< std::vector< double > >, std::array< std::size_t, 4 > > make_dot_integral_moments(const FiniteElement &V, cell::type celltype, std::size_t value_size, int q_deg)
Make interpolation points and weights for dot product integral moments.
Definition: moments.cpp:190