qmat.nodes
Implement the base NodesGenerator class that can be used to generate
quadrature nodes for many kind of distribution and types.
Its implementation is mostly based on the book of
[Gautschi, 2004].
Examples
>>> from qmat.nodes import NodesGenerator
>>> gen = NodesGenerator(nodeType="CHEBY-1", quadType="RADAU-RIGHT")
>>> nodes = gen.getNodes(10)
>>> coarse = gen.getNodes(5)
Attributes
Available types (distributions) for nodes |
|
Available quadrature type for each node distribution |
Classes
Class that can be used to generate generic distribution of nodes derived from Gauss quadrature rule. |
Module Contents
- NODE_TYPES = ['EQUID', 'LEGENDRE', 'CHEBY-1', 'CHEBY-2', 'CHEBY-3', 'CHEBY-4']
Available types (distributions) for nodes
- QUAD_TYPES = ['GAUSS', 'RADAU-LEFT', 'RADAU-RIGHT', 'LOBATTO']
Available quadrature type for each node distribution
- class NodesGenerator(nodeType='LEGENDRE', quadType='LOBATTO')[source]
Class that can be used to generate generic distribution of nodes derived from Gauss quadrature rule. Its implementation is fully inspired from [Gautschi, 2004].
- Parameters:
nodeType (str, optional) – The type of node distribution, can be
EQUID : equidistant nodes
LEGENDRE : node distribution from Legendre polynomials
CHEBY-1 : node distribution from Chebychev polynomials (1st kind)
CHEBY-2 : node distribution from Chebychev polynomials (2nd kind)
CHEBY-3 : node distribution from Chebychev polynomials (3rd kind)
CHEBY-4 : node distribution from Chebychev polynomials (4th kind)
The default is ‘LEGENDRE’.
quadType (str, optional) – The quadrature type, can be
GAUSS : inner point only, no node at boundary
RADAU-LEFT : only left boundary as node
RADAU-RIGHT : only right boundary as node
LOBATTO : left and right boundary as node
The default is ‘LOBATTO’.
- nodeType = 'LEGENDRE'
Type of the node distribution
- quadType = 'LOBATTO'
Quadrature type
- getNodes(nNodes)[source]
Computes a given number of quadrature nodes.
- Parameters:
nNodes (int) – Number of nodes to compute.
- Returns:
nodes – Nodes located in [-1, 1], in increasing order.
- Return type:
np.1darray
- getOrthogPolyCoefficients(num_coeff)[source]
Produces a given number of analytic three-term recurrence coefficients.
- Parameters:
num_coeff (int) – Number of coefficients to compute.
- Returns:
alpha (np.1darray) – The alpha coefficients of the three-term recurrence.
beta (np.1darray) – The beta coefficients of the three-term recurrence.
- evalOrthogPoly(t, alpha, beta)[source]
Evaluates the two higher order orthogonal polynomials corresponding to the given (alpha,beta) coefficients.
- Parameters:
t (float or np.1darray) – The point where to evaluate the orthogonal polynomials.
alpha (np.1darray) – The alpha coefficients of the three-term recurrence.
beta (np.1darray) – The beta coefficients of the three-term recurrence.
- Returns:
pi[0] (float or np.1darray) – The second higher order orthogonal polynomial evaluation.
pi[1] (float or np.1darray) – The higher oder orthogonal polynomial evaluation.
- getTridiagCoefficients(nNodes)[source]
Computes recurrence coefficients for the tridiagonal Jacobian matrix, taking into account the quadrature type.
- Parameters:
nNodes (int) – Number of nodes that should be computed from those coefficients.
- Returns:
alpha (np.1darray) – The modified alpha coefficients of the three-term recurrence.
beta (np.1darray) – The modified beta coefficients of the three-term recurrence.