|
AXOM
Axom provides a robust, flexible software infrastructure for the development of multi-physics applications and computational tools.
|
Consists of methods that evaluate scalar-field integrals on curves and regions defined by 2D curves, and vector-field integrals on curves. More...
#include "axom/core.hpp"#include "axom/config.hpp"#include "axom/core/utilities/Utilities.hpp"#include "axom/primal/geometry/CurvedPolygon.hpp"#include "axom/primal/operators/detail/evaluate_integral_curve_impl.hpp"#include <cmath>Namespaces | |
| axom | |
| axom::primal | |
Functions | |
Evaluates scalar-field line integrals for functions f : R^n -> R^m | |
| template<typename CurveType , typename Lambda , typename LambdaRetType = std::invoke_result_t<Lambda, typename CurveType::PointType>> | |
| LambdaRetType | axom::primal::evaluate_line_integral (const primal::CurvedPolygon< CurveType > &cpoly, Lambda &&integrand, int npts) |
| Evaluate a line integral along the boundary of a CurvedPolygon object for a function with an arbitrary return type. More... | |
| template<typename CurveType , typename Lambda , typename LambdaRetType = std::invoke_result_t<Lambda, typename CurveType::PointType>> | |
| LambdaRetType | axom::primal::evaluate_line_integral (const CurveType &c, Lambda &&integrand, int npts) |
| Evaluate a line integral along the boundary of a generic curve for a function with an arbitrary return type. More... | |
| template<typename CurveType , typename Lambda , typename LambdaRetType = std::invoke_result_t<Lambda, typename CurveType::PointType>> | |
| LambdaRetType | axom::primal::evaluate_line_integral (const axom::Array< CurveType > &carray, Lambda &&integrand, int npts) |
| Evaluate a line integral on an array of NURBS curves on a scalar field for a function with an arbitrary return type. More... | |
Evaluates vector-field line integrals for functions f : R^n -> R^n | |
| template<typename CurveType , typename Lambda , typename FuncRetType = typename CurveType::NumericType> | |
| FuncRetType | axom::primal::evaluate_vector_line_integral (const CurvedPolygon< CurveType > &cpoly, Lambda &&vector_integrand, int npts) |
| Evaluate a vector-field line integral along the boundary of a CurvedPolygon object. More... | |
| template<typename CurveType , typename Lambda , typename FuncRetType = typename CurveType::NumericType> | |
| FuncRetType | axom::primal::evaluate_vector_line_integral (const CurveType &c, Lambda &&vector_integrand, int npts) |
| Evaluate a vector-field line integral on a single generic curve. More... | |
| template<typename CurveType , typename Lambda , typename FuncRetType = typename CurveType::NumericType> | |
| FuncRetType | axom::primal::evaluate_vector_line_integral (const axom::Array< CurveType > &carray, Lambda &&vector_integrand, int npts) |
| Evaluate a line integral on an array of generic curves on a vector field. More... | |
Evaluates scalar-field 2D area integrals for functions f : R^2 -> R^m | |
| template<typename CurveType , typename Lambda , typename LambdaRetType = std::invoke_result_t<Lambda, typename CurveType::PointType>> | |
| LambdaRetType | axom::primal::evaluate_area_integral (const primal::CurvedPolygon< CurveType > &cpoly, Lambda &&integrand, int npts_Q, int npts_P=0) |
| Evaluate an integral on the interior of a CurvedPolygon object. More... | |
| template<typename CurveType , typename Lambda , typename LambdaRetType = std::invoke_result_t<Lambda, typename CurveType::PointType>> | |
| LambdaRetType | axom::primal::evaluate_area_integral (const axom::Array< CurveType > &carray, Lambda &&integrand, int npts_Q, int npts_P=0) |
| Evaluate an integral on the interior of a region bound by 2D curves. More... | |
Consists of methods that evaluate scalar-field integrals on curves and regions defined by 2D curves, and vector-field integrals on curves.
All integrals are evaluated numerically with Gauss-Legendre quadrature
Scalar-field line integrals and scalar-field area integrals are of the form int_D f(x) dr, with f : R^n -> R^m, D is a curve or a 2D region bound by curves
Vector-field line integrals are of form int_C f(x) \cdot d\vec{r}, with f : R^n -> R^n, C is a curve
2D area integrals computed with "Spectral Mesh-Free Quadrature for Planar Regions Bounded by Rational Parametric Curves" by D. Gunderman et al. https://doi.org/10.1016/j.cad.2020.102944