qmat.lagrange ============= .. py:module:: qmat.lagrange .. autoapi-nested-parse:: Base module for Barycentric Lagrange Approximation, based on `[Berrut & Trefethen, 2004] `_. Allows to easily build integration / interpolation / derivative matrices, from any list of node points. .. rubric:: Examples >>> # Base usage to generate a quadrature matrix >>> from qmat.lagrange import LagrangeApproximation, np >>> >>> grid = np.linspace(0, 1, num=5) >>> approx = LagrangeApproximation(grid) >>> Q = approx.getIntegrationMatrix([(0, tau) for tau in grid]) >>> >>> # Interpolation >>> fGrid = np.linspace(0, 1, num=200) >>> u = np.exp(grid) >>> P = approx.getInterpolationMatrix(fGrid) >>> uFine = P @ u >>> >>> # Alternative interpolation using the object as a function >>> uFine = approx(fGrid) >>> >>> # Derivative >>> D = approx.getDerivativeMatrix() >>> du = D @ u Classes ------- .. autoapisummary:: qmat.lagrange.LagrangeApproximation Functions --------- .. autoapisummary:: qmat.lagrange.computeFejerRule qmat.lagrange.getSparseInterpolationMatrix Module Contents --------------- .. py:function:: computeFejerRule(n) Compute a Fejer rule of the first kind, using DFT `[Waldvogel, 2006] `_. Inspired from quadpy (https://github.com/nschloe/quadpy @Nico_Schlömer) :Parameters: **n** (*int*) -- Number of points for the quadrature rule. :returns: * **nodes** (*np.1darray(n)*) -- The nodes of the quadrature rule * **weights** (*np.1darray(n)*) -- The weights of the quadrature rule. .. py:class:: LagrangeApproximation(points, weightComputation='AUTO', scaleWeights=False, scaleRef='MAX', duplicates='USE_LEFT', fValues=None) Class approximating any function on a given set of points using barycentric Lagrange interpolation. Let note :math:`(t_j)_{0\leq j numpy.ndarray Compute the interpolation matrix for a given set of discrete "time" points. For instance, if we note :math:`\vec{f}` the vector containing the :math:`f_j=f(t_j)` values, and :math:`(\tau_m)_{0\leq m numpy.ndarray Compute the integration matrix for a given set of intervals. For instance, if we note :math:`\vec{f}` the vector containing the :math:`f_j=f(t_j)` values, and :math:`(\tau_{m,left}, \tau_{m,right})_{0\leq m numpy.ndarray Generate derivative matrix of first or second order (or both) based on the Lagrange interpolant. The first order differentiation matrix :math:`D^{(1)}` approximates .. math:: D^{(1)} u \simeq \frac{du}{dx} on the interpolation points. The formula is : .. math:: D^{(1)}_{ij} = \frac{w_j/w_i}{x_i-x_j} for :math:`i \neq j` and .. math:: D^{(1)}_{jj} = -\sum_{i \neq j} D^{(1)}_{ij}` The second order differentiation matrix :math:`D^{(2)}` approximates .. math:: D^{(2)} u \simeq \frac{d^2u}{dx^2} on the interpolation points. The formula is : .. math:: D^{(1)}_{ij} = -2\frac{w_j/w_i}{x_i-x_j}\left[ \frac{1}{x_i-x_j} + \sum_{k \neq i}\frac{w_k/w_i}{x_i-x_k} \right] for :math:`i \neq j` and .. math:: D^{(2)}_{jj} = -\sum_{i \neq j} D^{(2)}_{ij} ⚠️ If you want a derivative matrix with many points (~1000 or more), favor the use of `weightComputation="STABLE"` when initializing the `LagrangeApproximation` object. If not, some (very small) weights could be approximated by zeros, which would make the computation of the derivative matrices fail ... .. note:: There is a typo in the formula for :math:`D^{(2)}` given in the paper of Berrut and Trefethen. The formula above is the correct one. :Parameters: * **order** (*int or str, optional*) -- The order of the derivative matrix, use "ALL" to retrieve both. The default is 1. * **duplicates** (*bool*) -- Wether or not take into account duplicates in the points. This has no impact if all interpolating points are distinct. Default is True. :returns: **D** -- Derivative matrix. If order="ALL", return a tuple containing all derivative matrices in increasing derivative order. :rtype: np.2darray or tuple of np.2darray .. py:method:: getDerivationMatrix(*args, **kwargs) -> numpy.ndarray .. py:function:: getSparseInterpolationMatrix(inPoints, outPoints, order, gridPeriod=-1) Get a sparse interpolation matrix from `inPoints` to `outPoints` of order `order` using barycentric Lagrange interpolation. The matrix will have `order` entries per line, and tends to be banded when both `inPoints` and `outPoints` are equispaced and cover the same interval. :Parameters: * **inPoints** (*np.1darray*) -- The points you want to interpolate from * **outPoints** (*np.1darray*) -- The points you want to interpolate to * **order** (*int*) -- Order of the interpolation * **grid_period** (*float*) -- Period of the grid. Negative values indicate non-periodic grids :returns: **A** -- Sparse interpolation matrix :rtype: scipy.sparse.csc_matrix(len(outPoints), len(inPoints))