BSpline Finite Element Exterior Calculus
Loading...
Searching...
No Matches
m_bspline_basis Module Reference

The B-spline basis module. More...

Data Types

type  bsplinefun
 The B-spline function type. More...
 
interface  bsplinespace
 The B-spline space type. More...
 
interface  evaluate
 Evaluate a B-spline (function) More...
 
interface  evaluate_support
 
interface  get_index_ranges
 Get the index ranges for integrating a B-spline or the product of two B-splines. More...
 
interface  get_petsc
 Convert a B-spline function to a PETSc vector. More...
 
interface  quadrature
 Compute quadrature of a B-spline. More...
 
interface  quadrature_product
 Compute quadrature of the product of two B-splines. More...
 
interface  size
 Get the size of the B-spline space (i.e., the number of linearly independent B-splines) More...
 

Functions/Subroutines

pure type(bsplinespace) function, public derivative (bspline, allow_periodic_degree0)
 Returns the B-spline space containing the derivatives of the given B-spline space.
 
pure type(bsplinespace) function, public anti_derivative (bspline)
 Returns the B-spline space containing the anti-derivatives of the given B-spline space.
 
pure real(wp) function, public evaluate_single (this, x, j, derivative)
 Evaluate the j-th B-spline of the B-spline space 'this' at x.
 
real(wp) function, public find_root_newton (this, x0, rhs, max_iter, tol, success)
 Find a root of the B-spline function 'this' minus rhs using Newton's method.
 
real(wp) function, public find_root_brent (this, x0, x1, rhs, max_iter, tol, success, res0, res1)
 Find a root of the B-spline function 'this' minus rhs using Brent's method.
 
pure integer function, public get_j_minus_i_quad (this, j, i)
 Determine j_actual - i, where the input j may exceed the number of B-splines.
 
pure subroutine, public get_j_internal (this, j, j_internal, bspline_mirrored)
 Get the internal B-spline index and whether the B-spline is mirrored.
 
pure subroutine, public get_internal_and_actual_inds (bspline, i, j_minus_i, j, j_internal, bspline_mirrored)
 Get the internal B-spline index and the actual B-spline index.
 
pure integer function, public get_interval_index (bspline, x)
 Get the interval index for a point in the B-spline space.
 

Detailed Description

The B-spline basis module.

This module provides

  • BSplineSpace: the space of B-spline basis functions
  • BSplineFun: the B-spline function which is represented by a linear combination of the basis functions from a space

The space is defined by the number of uniform intervals, the B-spline degree, as well as the boundary type:

  • The B-spline space can be clamped meaning that only the first B-spline is nonzero at the left boundary and only the last B-spline is nonzero at the right boundary.
  • Or the B-spline space can be periodic meaning that the first 'degree' B-splines are combined with the last 'degree' B-splines, thereby enforcing periodicity.

Function/Subroutine Documentation

◆ anti_derivative()

pure type(bsplinespace) function, public m_bspline_basis::anti_derivative ( type(bsplinespace), intent(in) bspline)

Returns the B-spline space containing the anti-derivatives of the given B-spline space.

Parameters
[in]bsplineThe B-spline space
Returns
The B-spline space containing the anti-derivatives of the given B-spline space
Note
This requires the B-spline space to be scaled

◆ derivative()

pure type(bsplinespace) function, public m_bspline_basis::derivative ( type(bsplinespace), intent(in) bspline,
logical, intent(in), optional allow_periodic_degree0 )

Returns the B-spline space containing the derivatives of the given B-spline space.

Parameters
[in]bsplineThe B-spline space
[in]allow_periodic_degree0Whether to allow periodic degree 0 B-splines in the 0-form space (default: false)
Returns
The B-spline space containing the derivatives of the given B-spline space
Note
This requires the B-spline space to be unscaled

◆ evaluate_single()

pure real(wp) function, public m_bspline_basis::evaluate_single ( type(bsplinespace), intent(in) this,
real(wp), intent(in) x,
integer, intent(in) j,
logical, intent(in), optional derivative )

Evaluate the j-th B-spline of the B-spline space 'this' at x.

Parameters
[in]thisThe B-spline space
[in]xThe evaluation point
[in]jThe B-spline index
[in]derivative_(optional)_ If true, the derivative of the B-spline is evaluated
Returns
The value of the B-spline at x

◆ find_root_brent()

real(wp) function, public m_bspline_basis::find_root_brent ( type(bsplinefun), intent(in) this,
real(wp), intent(in) x0,
real(wp), intent(in) x1,
real(wp), intent(in), optional rhs,
integer, intent(in), optional max_iter,
real(wp), intent(in), optional tol,
logical, intent(out), optional success,
real(wp), intent(in), optional res0,
real(wp), intent(in), optional res1 )

Find a root of the B-spline function 'this' minus rhs using Brent's method.

Parameters
[in]thisThe B-spline function
[in]x0The left-hand side of the bracketed interval for the root
[in]x1The right-hand side of the bracketed interval for the root
[in]rhs_(optional)_ The right-hand side value to find the root of (default: 0)
[in]max_iter_(optional)_ The maximum number of iterations for the root-finding method (default: 20)
[in]tol_(optional)_ The tolerance for convergence of the root-finding method (default: 1e-6)
[out]success_(optional)_ A logical flag indicating whether the root-finding method converged successfully
[in]res0_(optional)_ The value of the function minus rhs at x0 (if cheaply available)
[in]res1_(optional)_ The value of the function minus rhs at x1 (if cheaply available)

◆ find_root_newton()

real(wp) function, public m_bspline_basis::find_root_newton ( type(bsplinefun), intent(in) this,
real(wp), intent(in), optional x0,
real(wp), intent(in), optional rhs,
integer, intent(in), optional max_iter,
real(wp), intent(in), optional tol,
logical, intent(out), optional success )

Find a root of the B-spline function 'this' minus rhs using Newton's method.

Parameters
[in]thisThe B-spline function
[in]x0_(optional)_ The initial guess for the root (default: 0.5)
[in]rhs_(optional)_ The right-hand side value to find the root of (default: 0)
[in]max_iter_(optional)_ The maximum number of iterations for the root-finding method (default: 20)
[in]tol_(optional)_ The tolerance for convergence of the root-finding method (default: 1e-8)
[out]success_(optional)_ A logical flag indicating whether the root-finding method converged successfully
Returns
The root of the B-spline function 'this' minus rhs

◆ get_internal_and_actual_inds()

pure subroutine, public m_bspline_basis::get_internal_and_actual_inds ( type(bsplinespace), intent(in) bspline,
integer, intent(in) i,
integer, intent(in) j_minus_i,
integer, intent(out) j,
integer, intent(out), optional j_internal,
logical, intent(out), optional bspline_mirrored )

Get the internal B-spline index and the actual B-spline index.

Parameters
[in]bsplineThe B-spline space
[in]iThe interval index
[in]j_minus_iThe difference between the B-spline index and the interval index
[out]jThe actual B-spline index (0:nr_bsplines-1)
[out]j_internal_(optional)_ The internal B-spline index (0:degree)
[out]bspline_mirrored_(optional)_ Whether the B-spline is mirrored (i.e., the B-spline is on the right boundary of the domain)

◆ get_interval_index()

pure integer function, public m_bspline_basis::get_interval_index ( type(bsplinespace), intent(in) bspline,
real(wp), intent(in) x )

Get the interval index for a point in the B-spline space.

Parameters
[in]bsplineThe B-spline space
[in]xThe point
Returns
The interval index for the point in the B-spline space (0:nr_intervals-1)

◆ get_j_internal()

pure subroutine, public m_bspline_basis::get_j_internal ( type(bsplinespace), intent(in) this,
integer, intent(in) j,
integer, intent(out) j_internal,
logical, intent(out) bspline_mirrored )

Get the internal B-spline index and whether the B-spline is mirrored.

Parameters
[in]thisThe B-spline space
[in]jThe B-spline index
[out]j_internalThe internal B-spline index (0:degree)
[out]bspline_mirrored_(optional)_ Whether the B-spline is mirrored (i.e., the B-spline is on the right boundary of the domain)

◆ get_j_minus_i_quad()

pure integer function, public m_bspline_basis::get_j_minus_i_quad ( type(bsplinespace), intent(in) this,
integer, intent(in) j,
integer, intent(in) i )

Determine j_actual - i, where the input j may exceed the number of B-splines.

Parameters
[in]thisThe B-spline basis
[in]jThe B-spline index
[in]iThe interval index
Returns
The j_actual - i