core.Collocation module

class CollBase(num_nodes=None, tleft=0, tright=1, node_type='LEGENDRE', quad_type=None, **kwargs)[source]

Bases: object

Generic collocation class, that contains everything to do integration over intervals and between nodes. It can be used to produce many kind of quadrature nodes from various distribution (awesome!).

It is based on the two main parameters that define the nodes :

  • node_type : the node distribution used for the collocation method

  • quad_type : the type of quadrature used (inclusion of not of boundary)

Current implementation provides the following available parameter values for node_type :

  • EQUID : equidistant node distribution

  • LEGENDRE : distribution from Legendre polynomials

  • CHEBY-{1,2,3,4} : distribution from Chebyshev polynomials of a given kind

The type of quadrature cann be GAUSS (only inner nodes), RADAU-LEFT (inclusion of the left boundary), RADAU-RIGHT (inclusion of the right boundary) and LOBATTO (inclusion of left and right boundary).

Here is the equivalency table with the (old) original classes implemented in pySDC :

Original Class

node_type

quad_type

Equidistant

EQUID

LOBATTO

EquidistantInner

EQUID

GAUSS

EquidistantNoLeft

EQUID

RADAU-RIGHT

CollGaussLegendre

LEGENDRE

GAUSS

CollGaussLobatto

LEGENDRE

LOBATTO

CollGaussRadau_Left

LEGENDRE

RADAU-LEFT

CollGaussRadau_Right

LEGENDRE

RADAU-RIGHT

num_nodes

number of collocation nodes

Type:

int

tleft

left interval point

Type:

float

tright

right interval point

Type:

float

nodes

array of quadrature nodes

Type:

numpy.ndarray

weights

array of quadrature weights for the full interval

Type:

numpy.ndarray

Qmat

matrix containing the weights for tleft to node

Type:

numpy.ndarray

Smat

matrix containing the weights for node to node

Type:

numpy.ndarray

delta_m

array of distances between nodes

Type:

numpy.ndarray

right_is_node

flag to indicate whether right point is collocation node

Type:

bool

left_is_node

flag to indicate whether left point is collocation node

Type:

bool

static evaluate(weights, data)[source]

Evaluates the quadrature over the full interval

Parameters:
  • weights (numpy.ndarray) – array of quadrature weights for the full interval

  • data (numpy.ndarray) – f(x) to be integrated

Returns:

integral over f(x) between tleft and tright

Return type:

numpy.ndarray