boundaryscheme package#
Submodules#
boundaryscheme.boundaries module#
This file defines the boundaries classes
- class boundaryscheme.boundaries.Boundary(B)[source]#
Bases:
object
This is a class to represent the boundary condition.
- class boundaryscheme.boundaries.DDJ(d, kd, dx=0.1, a=1, gderivatives=None)[source]#
Bases:
Boundary
This is a class to represent the boundary Reconstruction procedure explained in [Boutin, Le Barbenchon, Seguin, 2023 : Stability of finite difference schemes for the hyperbolic initial boundary value problem by winding number computations] and in [Dakin Desprès Jaouen, 2018 : Inverse Lax–Wendroff boundary treatment for compressible Lagrange-remap hydrodynamics on cartesian grids]
Warning
Warning d represents the order and kd represente the truncature, the notation is not in the same order than SILW, it is to respect the order of [DakinDespresJaouen18]
- Parameters:
kd (int) – Index of the beginning of extrapolation using
d (int) – Order of consistency of the boundary condition
dx (float, optional) – Space discretization, defaults to 1/10
a (float, optional) – Velocity, defaults to 1
gderivatives (list, optional) – List of the boundary data and its derivatives, defaults to None
- class boundaryscheme.boundaries.Dirichlet[source]#
Bases:
Boundary
This is a class to represent the homogeneous dirichlet boundary condition.
- Parameters:
d (int) – Order of consistency
- class boundaryscheme.boundaries.SILW(kd, d, dx=0.1, a=1, gderivatives=None)[source]#
Bases:
Boundary
This is a class to represent the Simplified Inverse Lax Wendroff boundary procedure, see [Vilar,Shu, 2015 : Development and stability analysis of the inverse Lax Wendroff boundary treatment for central compact schemes]
- Parameters:
kd (int) – Index of the beginning of extrapolation using
d (int) – Order of consistency of the boundary condition
dx (float, optional) – Space discretization, defaults to 1/10
a (float, optional) – Velocity, defaults to 1
gderivatives (list, optional) – List of the boundary data and its derivatives, defaults to None
boundaryscheme.complex_winding_number module#
This file implements a routine to perform the “point in polygon” inclusion test
- boundaryscheme.complex_winding_number.Indice(L)[source]#
take a list of (x_i,x_j) representing a polygon return the winding number of the origin with respect to the curve
boundaryscheme.pyplot module#
This file defines utils function specific to the library
- boundaryscheme.pyplot.detKLplot(schem, left_bound=<boundaryscheme.boundaries.Dirichlet object>, lamb=None, sigma=0, lambdacursor=False, sigmacursor=False, nparam=300, parametrization_bool=True, fig_size=(10, 8))[source]#
Computes the Kreiss-Lopatinskii determinant curve for different lambdas and different sigmas or with cursor(s)
- Parameters:
schem (class:Scheme) – Scheme depending on a parameter lambda
left_bound (class:Boundary, optional) – Left boundary of the class Boundary, defaults to Dirichlet()
lamb (int or float or list or numpy.ndarray, optional) – A value lambda or a numpy.ndarray of several values of lambda (lambda is the Courant number a.dt/dx), defaults to None
sigma (int or float or list or numpy.ndarray, optional) – Gap between the mesh and the boundary condition, defaults to 0
lambdacursor (bool, optional) – Boolean to indicate the use of a cursor for lambda’s moving among the lamb values, defaults to False
sigmacursor (bool, optional) – Boolean to indicate the use of a cursor for sigma’s moving among the sigma values, defaults to False
n_param (int, optional) – Size of the discretization of the unit circle, defaults to 300
parametrization_bool (bool, optional) – True to have a refinment close to the origin, False to have a uniform discretization of the unit circle, defaults to True
fig_size ((int,int), optional) – Size of the figure, defaults to (6,4)
- Returns:
the Kreiss-Lopatinskii curve for the scheme
- Return type:
plot object
- boundaryscheme.pyplot.detKLplotsimple(schem, n_param=300, parametrization_bool=True)[source]#
Computes the Kreiss-Lopatinskii determinant curve for a given scheme
- Parameters:
schem (class:Scheme) – Scheme considered
n_param (int, optional) – Size of the discretization of the unit circle, defaults to 300
parametrization_bool (bool, optional) – True to have a refinment close to the origin, False to have a uniform discretization of the unit circle, defaults to True
- Returns:
the Kreiss-Lopatinskii curve
- Return type:
plot object
- boundaryscheme.pyplot.nbrzerosdetKL(schem, left_bound=<boundaryscheme.boundaries.Dirichlet object>, lamb=None, sigma=0, nparam=100, nlambda=200, nsigma=30, parametrization_bool=True, fig_size=(15, 7), fontsize=18)[source]#
Computes the number of zeros of the Kreiss-Lopatinskii determinant for different boundary conditions, for different lambdas or/and for different sigmas
See also
for more details and examples, see [Boutin, Le Barbenchon, Seguin, 2023 : Stability of finite difference schemes for the hyperbolic initial boundary value problem by winding number computations]
- Parameters:
schem (class:Scheme) – Scheme depending on a parameter lambda
left_bound (class:Boundary or list of class:Boundary, optional) – Left boundary of the class Boundary, defaults to Dirichlet()
lamb (int or float or list or numpy.ndarray, optional) – A value lambda or a numpy.ndarray of several values of lambda (lambda is the Courant number a.dt/dx), defaults to None
sigma (int or float or list or numpy.ndarray, optional) – Gap between the mesh and the boundary condition, defaults to 0
nparam (int, optional) – Size of the discretization of the unit circle for the Kreiss-Lopatinskii curve, defaults to 100
nlambda (int, optional) – Size of the discretization of the lambda’s, defaults to 200
nsigma (int, optional) – Size of the discretization of the sigma’s, defaults to 30
parametrization_bool (bool, optional) – True to have a refinment close to the origin, False to have a uniform discretization of the unit circle, defaults to True
fig_size ((int,int), optional) – Size of the figure, defaults to (15,7)
fontsize (int , optional) – Size of the font of the text in the plot, defaults to 18
- Returns:
the number of zeros of the Kreiss-Lopatinskii determinant for the scheme
- Return type:
plot object
- boundaryscheme.pyplot.symbolplot(schem, lamb=None, lambdacursor=False, nparam=300, fig_size=(10, 8))[source]#
Computes the Kreiss-Lopatinskii determinant curve for different lambdas and different sigmas or with cursor(s)
- Parameters:
schem (class:Scheme) – Scheme depending on a parameter lambda
lamb (int or float or list or numpy.ndarray, optional) – A value lambda or a numpy.ndarray of several values of lambda (lambda is the Courant number a.dt/dx), defaults to None
lambdacursor (bool, optional) – Boolean to indicate the use of a cursor for lambda’s moving among the lamb values, defaults to False
nparam (int, optional) – Size of the discretization of the unit circle, defaults to 300
fig_size ((int,int), optional) – Size of the figure, defaults to (10,8)
- Returns:
the symbol curve of the scheme
- Return type:
plot object
boundaryscheme.schemes module#
This file defines the schemes classes
- class boundaryscheme.schemes.BB(lamb, boundary=<boundaryscheme.boundaries.Dirichlet object>, sigma=0, **kwargs)[source]#
Bases:
Scheme
This is a class to represent an example of numerical scheme.
- Parameters:
lamb (float) – The Courant number, i.e a.dt/dx where “a” is the velocity, “dt” the time discretization and “dx” the space discretization
boundary (class:Boundary, optional) – Boundary condition, defaults to Dirichlet()
sigma (float, optional) – Gap between the mesh and the boundary condition, defaults to 0
- class boundaryscheme.schemes.BeamWarming(lamb, boundary=<boundaryscheme.boundaries.Dirichlet object>, sigma=0, **kwargs)[source]#
Bases:
SchemeP0
This is a class to represent the Beam-Warming scheme (which is totally upwind).
- Parameters:
lamb (float) – The Courant number, i.e a.dt/dx where “a” is the velocity, “dt” the time discretization and “dx” the space discretization
boundary (class:Boundary, optional) – Boundary condition, defaults to Dirichlet()
sigma (float, optional) – Gap between the mesh and the boundary condition, defaults to 0
- class boundaryscheme.schemes.Dissipatif(lamb, boundary=<boundaryscheme.boundaries.Dirichlet object>, sigma=0, D=0, **kwargs)[source]#
Bases:
Scheme
This is a class to represent a dissipative scheme.
- Parameters:
lamb (float) – The Courant number, i.e a.dt/dx where “a” is the velocity, “dt” the time discretization and “dx” the space discretization
boundary (class:Boundary, optional) – Boundary condition, defaults to Dirichlet()
sigma (float, optional) – Gap between the mesh and the boundary condition, defaults to 0
D (float, optional) – A parameter of the dissipativity, defaults to 0
- class boundaryscheme.schemes.LaxFriedrichs(lamb, boundary=<boundaryscheme.boundaries.Dirichlet object>, sigma=0, **kwargs)[source]#
Bases:
Scheme
This is a class to represent the Lax-Friedrichs scheme.
- Parameters:
lamb (float) – The Courant number, i.e a.dt/dx where “a” is the velocity, “dt” the time discretization and “dx” the space discretization
boundary (class:Boundary, optional) – Boundary condition, defaults to Dirichlet()
sigma (float, optional) – Gap between the mesh and the boundary condition, defaults to 0
- boundaryscheme.schemes.LaxWendroff(order=2)[source]#
This is a function which create the class of LaxWendroff scheme at order ‘order’ :param order: Order of the Lax-Wendroff scheme, defaults to 2 :type order: int, optional
- class boundaryscheme.schemes.Scheme(inter, center, boundary=<boundaryscheme.boundaries.Dirichlet object>, sigma=0, **kwargs)[source]#
Bases:
object
This is a class to represent a numerical scheme.
- Parameters:
inter (list) – List of float coefficient that represent the interior scheme
center (int) – Index of the central coefficient of the scheme
boundary (class:Boundary, optional) – Boundary condition, defaults to Dirichlet()
sigma (float, optional) – Gap between the mesh and the boundary condition, defaults to 0
See also
[Boutin, Le Barbenchon, Seguin, 2023 : Stability of finite difference schemes for the hyperbolic initial boundary value problem by winding number computations]
- Kappa(z0)[source]#
Selection of kappas that come from the inside of the unit disk :param z0: Parameter z0 of the characteristic polynomial :type z0: complex :return: Roots of the characteristic polynomial for z0 that are inside the unit disk or are coming from the inside of the unit disk :rtype: list
- count_root(eta, eps, z0, kappa)[source]#
Counts the number of roots from kappa(z0) that comes from the inside of the unit disk :param eta: Gap between z and z0 :type eta: float :param eps: Threshold :type eps: float :param z0: Parameter z0 :type z0: complex :param kappa: Multiple root kappa linked to z0 :type kappa: complex :return: Number of roots coming from the inside of the unit disk :rtype: int
- detKL(n_param, parametrization_bool)[source]#
Computes the Kreiss-Lopatinskii determinant on the unit circle :param n_param: Number of discretization of the unit circle :type n_param: int :param parametrization_bool: If True, the discretization is refined close to the origin. If False, the discretization is uniform on the unit circle. :type parametrization_bool: bool :return: The parametrization of the unit circle and the value of the Kreiss-Lopatinskii determinant on the unit circle :rtype: (numpy.ndarray, numpy.ndarray)
- pol(z)[source]#
Returns the characteristic polynomial of the scheme :param z: Parameter z of the characteristic polynomial :type z: complex :return: Characteristic polynomial of the scheme :rtype: class: Polynomial
- roots(z)[source]#
Returns the roots of the characteristic polynomial :param z: Parameter z of the characteristic polynomial :type z: complex :return: Roots of the characteristic polynomial sorting by utils.sort function :rtype: list
- scheme()[source]#
Returns the parameters of the scheme : the interior coefficients et the center :return: Parameters of the scheme :rtype: (list,int)
- symbol(n=300)[source]#
To draw the symbol curve of the scheme :param n: Number of points of the parametrisation of the unit circle, defaults to 300 :type n: int, optional :return: Abscissa and Ordinate of the symbol curve of the scheme :rtype: (list,list)
- toep(J, right_bound=array([], shape=(1, 0), dtype=float64))[source]#
Returns the Toeplitz matrix of the scheme with two boundaries :param J: Size of the Toeplitz matrix :type J: int :param right_bound: Right boundary of the scheme :type right_bound: class:Boundary :return: Toeplitz matrix of the scheme :rtype: numpy.ndarray
- class boundaryscheme.schemes.SchemeP0(inter, center, boundary, sigma=0)[source]#
Bases:
Scheme
This is a class to represent a totally upwind scheme (with p = 0).
- Parameters:
inter (list) – List of float coefficient that represent the interior scheme
center (int) – Index of the central coefficient of the scheme
boundary (class:Boundary, optional) – Boundary condition, defaults to Dirichlet()
sigma (float, optional) – Gap between the mesh and the boundary condition, defaults to 0
See also
[Boutin, Le Barbenchon, Seguin, 2023 : On the stability of totally upwind schemes for the hyperbolic initial boundary value problem]
- detKL(n_param, parametrization_bool)[source]#
Redefinition of the detKL method of the :class: Scheme but in the particular case of totally upwind scheme :param n_param: Number of discretization of the unit circle :type n_param: int :param parametrization_bool: If True, the discretization is refined close to the origin. If False, the discretization is uniform on the unit circle. :type parametrization_bool: bool :return: The parametrization of the unit circle and the value of the Kreiss-Lopatinskii determinant on the unit circle :rtype: (numpy.ndarray, numpy.ndarray)
- class boundaryscheme.schemes.ThirdOrder(lamb, boundary=<boundaryscheme.boundaries.Dirichlet object>, sigma=0, **kwargs)[source]#
Bases:
Scheme
This is a class to represent the O3-scheme.
- Parameters:
lamb (float) – The Courant number, i.e a.dt/dx where “a” is the velocity, “dt” the time discretization and “dx” the space discretization
boundary (class:Boundary, optional) – Boundary condition, defaults to Dirichlet()
sigma (float, optional) – Gap between the mesh and the boundary condition, defaults to 0
- class boundaryscheme.schemes.Upwind(lamb, boundary=<boundaryscheme.boundaries.Dirichlet object>, sigma=0, **kwargs)[source]#
Bases:
SchemeP0
This is a class to represent the Upwind scheme (which is totally upwind).
- Parameters:
lamb (float) – The Courant number, i.e a.dt/dx where “a” is the velocity, “dt” the time discretization and “dx” the space discretization
boundary (class:Boundary, optional) – Boundary condition, defaults to Dirichlet()
sigma (float, optional) – Gap between the mesh and the boundary condition, defaults to 0
boundaryscheme.utils module#
This file defines utils function specific to the library
- boundaryscheme.utils.Lwconstructor(order)[source]#
Constructor of Lax-Wendroff scheme.
- Parameters:
order (int) – Order of Lax Wendroff scheme
- Returns:
Function depending on the Courant number lambda and maps to the parameters of the Lax-Wendroff scheme
- Return type:
function
- boundaryscheme.utils.boundary_quasi_toep_to_boundary(quasi_toep_boundary, bn, inter, center)[source]#
Transformation from a boundary condition written as the extraction of the r first rows of the Quasi-Toeplitz matrix to a boundary condition written with ghost points. :param quasi_toep_boundary: Boundary condition BB written as the extraction of the r first rows of a Quasi-Toeplitz matrix :type quasi_toep_boundary: numpy.ndarray :param bn: Boundary data bn such that (U_0^{n+1},…,U_{r-1}^{n+1}) = BB (U_0^n,…,U_{m-1}^n) + b_n(t^n) :type bn: function :param inter: List of float coefficient that represent the interior scheme :type inter: list :param center: Index of the central coefficient of the scheme :type center: int :return: Boundary condition B written as U_{-r} = …, …, U_{-1} =… and the function g_n such that (U_{-r}^n,…,U_{-1}^n) = B(U_0^n,…,U_{m-1}^n) + g_n(t^n) :rtype: (numpy.ndarray,function)
- boundaryscheme.utils.boundary_to_boundary_quasi_toep(boundary_condition, gn, inter, center)[source]#
Transformation from a boundary condition written with ghost points to a boundary condition written as the extraction of the r first rows of the Quasi-Toeplitz matrix
- Parameters:
boundary_condition (numpy.ndarray) – Boundary condition written as U_{-r} = …, …,U_{-1} =…, etc
gn (function) – Boundary data
inter (list) – List of float coefficient that represent the interior scheme
center (int) – Index of the central coefficient of the scheme
- Returns:
Extraction of the r first rows of the Quasi-Toeplitz matrix and the function b_n(t) such that (U_0^{n+1},…,U_{r-1}^{n+1}) = T_J (U_0^n,…,U_{m-1}^n) + b_n(t^n)
- Return type:
(numpy.ndarray,function)
- boundaryscheme.utils.coefBinomial(n, k)[source]#
Computes the binomial coefficient
- Parameters:
n (int) – Number of elements
k (int) – Number of choosen elements
- Returns:
Binomial coefficient “n choose k”
- Return type:
int
- boundaryscheme.utils.epsilon(L)[source]#
Computes the half of the smallest distance between two distincts elements of the list L
- Parameters:
L (list) – List of complex numbers (z_i)
- Returns:
min |z_i-z_j| / 2 (when z_i != z_j)
- Return type:
float
- boundaryscheme.utils.eta_func(eps, kappa0, N, polynom, r)[source]#
Computes the distance eta .. note:: for more details, see [Boutin, Le Barbenchon, Seguin, 2023 : Stability of finite difference schemes for the hyperbolic initial boundary value problem by winding number computations]
- Parameters:
eps (float) – Radius
kappa0 (complex) – Center
N (int) – Number of discretization of the circle centered in kappa0 of radius eps
polynom (Polynomial) – Polynomial
r (int) – Integer
- Returns:
min |polynom(kappa)|/(1+eps)^r for kappa on the circle centered in kappa0 of radius eps
- Return type:
float
- boundaryscheme.utils.neighboor(sector1, z2)[source]#
Checks if z2 is in a neighboring sector of sector1
- Parameters:
sector1 (int) – Integer between 0 and 7 and represent the angular sector [sector1*pi/4, (sector1+1)*pi/4[
z2 (complex) – Complex number
- Returns:
True if and only if the complex z2 is in a neighboring sector of [sector1*pi/4, (sector1+1)*pi/4[
- Return type:
bool
- boundaryscheme.utils.parametrization(n_param, curve_formula)[source]#
Finds the discretization of the curve refine if it is needed (as [ZapataMartin2014] procedure).
- Parameters:
n_param (int) – Integer for the default discretization of [0,2pi]
curve_formula – A fonction : z in the unit circle mapsto a complex and represent a curve
- Returns:
The discretization of the unit circle and the values of the curve for this discretization
- Return type:
(list,list)
- boundaryscheme.utils.recflux(N, lamb)[source]#
Computes the coefficient of the Lax-Wendroff scheme at any order. Recursive function.
- Parameters:
N (int) – Order of Lax Wendroff scheme
lamb (float) – The Courant number, i.e a.dt/dx where “a” is the velocity, “dt” the time discretization and “dx” the space discretization
- Returns:
Parameters of the Lax-Wendroff scheme : interior coefficient and the center.
- Return type:
(np.ndarray,int)
- boundaryscheme.utils.schemscal(schem, scal)[source]#
Multiplication of a scheme by a scalar
- Parameters:
schem ((list,int)) – Representation of a numerical scheme with the interior coefficients and the center
scal (float) – Scalar element
- Returns:
Multiplication of schem by scal
- Return type:
(list,int)
- boundaryscheme.utils.schemsum(schem1, schem2)[source]#
Addition of two numerical schemes
- Parameters:
schem1 ((list,int)) – Representation of a numerical scheme with the interior coefficients and the center
schem2 ((list,int)) – Representation of a numerical scheme with the interior coefficients and the center
- Returns:
Sum of schem1 and schem2
- Return type:
(list,int)
- boundaryscheme.utils.sector(z)[source]#
Finds the sector of z
Warning
the argument is between -pi and pi and here “a” is between 0 and 7
- Parameters:
z (complex) – Complex number
- Returns:
The sector of z, it is an integer “a” between 0 and 7 which correspond to the angular sector [a*pi/4, (a+1)*pi/4[
- Return type:
int