Public API

Integrals

MeshIntegrals.integralFunction
integral(f, geometry[, rule]; FP=Float64)

Numerically integrate a given function f(::Point) over the domain defined by a geometry using a particular numerical integration rule with floating point precision of type FP.

Arguments

  • f: an integrand function with a method f(::Meshes.Point)
  • geometry: some Meshes.Geometry that defines the integration domain
  • rule: optionally, the IntegrationRule used for integration (by default GaussKronrod() in 1D and HAdaptiveCubature() else)
  • FP: optionally, the floating point precision desired (Float64 by default)

Note that reducing FP below Float64 will incur some loss of precision. By contrast, increasing FP to e.g. BigFloat will typically increase precision (at the expense of longer runtimes).

source

Aliases

MeshIntegrals.lineintegralFunction
lineintegral(f, geometry[, rule]; FP=Float64)

Numerically integrate a given function f(::Point) along a line-like geometry using a particular numerical integration rule with floating point precision of type FP.

Rule types available:

  • GaussKronrod (default)
  • GaussLegendre
  • HAdaptiveCubature
source
MeshIntegrals.surfaceintegralFunction
surfaceintegral(f, geometry[, rule]; FP=Float64)

Numerically integrate a given function f(::Point) along a surface geometry using a particular numerical integration rule with floating point precision of type FP.

Algorithm types available:

  • GaussKronrod
  • GaussLegendre
  • HAdaptiveCubature (default)
source
MeshIntegrals.volumeintegralFunction
volumeintegral(f, geometry[, rule]; FP=Float64)

Numerically integrate a given function f(::Point) throughout a volumetric geometry using a particular numerical integration rule with floating point precision of type FP.

Algorithm types available:

  • GaussKronrod
  • GaussLegendre
  • HAdaptiveCubature (default)
source

Integration Rules

MeshIntegrals.GaussKronrodType
GaussKronrod(kwargs...)

The h-adaptive Gauss-Kronrod quadrature rule implemented by QuadGK.jl. All standard QuadGK.quadgk keyword arguments are supported. This rule works natively for one dimensional geometries; some two- and three-dimensional geometries are additionally supported using nested integral solvers with the specified kwarg settings.

source
MeshIntegrals.GaussLegendreType
GaussLegendre(n)

An n'th-order Gauss-Legendre quadrature rule. Nodes and weights are efficiently calculated using FastGaussQuadrature.jl.

So long as the integrand function can be well-approximated by a polynomial of order 2n-1, this method should yield results with 16-digit accuracy in O(n) time. If the function is know to have some periodic content, then n should (at a minimum) be greater than the expected number of periods over the geometry, e.g. length(geometry)/λ.

source
MeshIntegrals.HAdaptiveCubatureType
HAdaptiveCubature(kwargs...)

The h-adaptive cubature rule implemented by HCubature.jl. All standard HCubature.hcubature keyword arguments are supported.

source

Derivatives

MeshIntegrals.jacobianFunction
jacobian(geometry, ts; ε=1e-6)

Calculate the Jacobian of a geometry at some parametric point ts using a simple finite-difference approximation with step size ε.

Arguments

  • geometry: some Meshes.Geometry of N parametric dimensions
  • ts: a parametric point specified as a vector or tuple of length N
  • ε: step size to use for the finite-difference approximation
source