Python API¶
Référence générée automatiquement depuis le module compilé poisson_cpp.
Grilles¶
- class poisson_cpp.Grid1D¶
Bases:
pybind11_object1D node-centered grid covering
[0, L]withNequispaced nodes.- property L¶
Domain length.
- property N¶
Number of nodes.
- property dx¶
Spacing
L / (N - 1).
- x(self: poisson_cpp._core.Grid1D, arg0: SupportsInt | SupportsIndex) float¶
Return the node coordinates as a length-
Nnumpy array.
Solveurs 1D¶
- poisson_cpp.thomas(a: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]'], b: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]'], c: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]'], d: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']¶
Solve a tridiagonal linear system
A x = dby the Thomas algorithm.The matrix
Ahasaon the sub-diagonal,bon the main diagonal, andcon the super-diagonal. Runs inO(N).- Parameters:
a (numpy.ndarray of shape (N,)) – Sub-diagonal;
a[0]is unused.b (numpy.ndarray of shape (N,)) – Main diagonal.
c (numpy.ndarray of shape (N,)) – Super-diagonal;
c[N-1]is unused.d (numpy.ndarray of shape (N,)) – Right-hand side.
- Returns:
Solution vector
x.- Return type:
numpy.ndarray of shape (N,)
- Raises:
RuntimeError – If a zero pivot is encountered (matrix not diagonally dominant).
- poisson_cpp.solve_poisson_1d(rho: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], uL: SupportsFloat | SupportsIndex, uR: SupportsFloat | SupportsIndex, grid: poisson_cpp._core.Grid1D) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']¶
Finite-volume 1D Poisson with Dirichlet BC, uniform permittivity.
Solves
-V''(x) = rho(x) / eps0on a node-centered grid withV(0) = uLandV(L) = uR.- Parameters:
rho (numpy.ndarray of shape (N,)) – Charge density at each node.
uL (float) – Dirichlet values at
x = 0andx = L.uR (float) – Dirichlet values at
x = 0andx = L.grid (Grid1D) – Discretization grid.
- Returns:
Potential
Vat each node.- Return type:
numpy.ndarray of shape (N,)
- poisson_cpp.solve_poisson_1d_dielectric(rho: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], eps_r: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], uL: SupportsFloat | SupportsIndex, uR: SupportsFloat | SupportsIndex, grid: poisson_cpp._core.Grid1D, eps0: SupportsFloat | SupportsIndex = 1.0) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']¶
Finite-volume 1D Poisson with spatially-varying permittivity.
Solves
eps0 * d/dx[eps_r(x) dV/dx] = -rhoon a node-centered grid withV(0) = uLandV(L) = uR. The face permittivity is the harmonic mean of the two adjacent cell values, which preserves the normal component ofD = eps0 eps_r grad Vacross a dielectric interface.- Parameters:
rho (numpy.ndarray of shape (N,)) – Charge density at each node.
eps_r (numpy.ndarray of shape (N,)) – Relative permittivity at each node. Must be strictly positive.
uL (float) – Dirichlet values at
x = 0andx = L.uR (float) – Dirichlet values at
x = 0andx = L.grid (Grid1D) – Discretization grid.
eps0 (float, optional) – Vacuum permittivity. Default
1.0.
- Returns:
Potential
Vat each node.- Return type:
numpy.ndarray of shape (N,)
Solveur 2D itératif (SOR)¶
- class poisson_cpp.Solver2D¶
Bases:
pybind11_object2D finite-volume Poisson solver using SOR red-black with optimal omega.
Boundary conditions: Dirichlet in x (
uLat left,uRat right), homogeneous Neumann in y. Permittivity may be constant or spatial.- Parameters:
grid (Grid2D) – Cell-centered grid.
eps (float or numpy.ndarray) – Constant permittivity, or per-cell array of shape
(Nx, Ny).uL (float) – Dirichlet values on the x-boundaries.
uR (float) – Dirichlet values on the x-boundaries.
- solve(self: poisson_cpp._core.Solver2D, rho: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, n]'], omega: SupportsFloat | SupportsIndex = -1.0, tol: SupportsFloat | SupportsIndex = 1e-08, max_iter: SupportsInt | SupportsIndex = 20000) tuple[Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]'], poisson_cpp._core.SORReport]¶
Run SOR from a zero initial guess.
- Parameters:
rho (numpy.ndarray of shape (Nx, Ny)) – Charge density per cell.
omega (float, optional) – Relaxation factor. Default
-1selectsomega_opt = 2 / (1 + sin(pi/N)).tol (float, optional) – Convergence threshold on max-norm residual. Default
1e-8.max_iter (int, optional) – Iteration cap. Default
20_000.
- Returns:
V (numpy.ndarray of shape (Nx, Ny)) – Potential field.
report (SORReport) – Iteration count and final residual.
- solve_inplace(self: poisson_cpp._core.Solver2D, V: Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]', 'flags.writeable', 'flags.f_contiguous'], rho: Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]', 'flags.f_contiguous'], omega: SupportsFloat | SupportsIndex = -1.0, tol: SupportsFloat | SupportsIndex = 1e-08, max_iter: SupportsInt | SupportsIndex = 20000) poisson_cpp._core.SORReport¶
In-place variant of
solve().Vandrhomust be Fortran-ordered numpy arrays (order='F') matching Eigen’s column-major layout, otherwise pybind11 will copy them.- Parameters:
V (numpy.ndarray of shape (Nx, Ny), order='F') – Initial guess; overwritten with the solution.
rho (numpy.ndarray of shape (Nx, Ny), order='F') – Charge density per cell.
omega – Same meaning as
solve().tol – Same meaning as
solve().max_iter – Same meaning as
solve().
- Returns:
Iteration count and final residual.
- Return type:
Gradient Conjugué¶
- poisson_cpp.solve_poisson_cg(V: Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]', 'flags.writeable', 'flags.f_contiguous'], rho: Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]', 'flags.f_contiguous'], grid: poisson_cpp._core.Grid2D, eps: SupportsFloat | SupportsIndex = 1.0, uL: SupportsFloat | SupportsIndex = 0.0, uR: SupportsFloat | SupportsIndex = 0.0, tol: SupportsFloat | SupportsIndex = 1e-08, max_iter: SupportsInt | SupportsIndex = 10000, use_preconditioner: bool = False, record_history: bool = False) tuple[poisson_cpp._core.CGReport, list[float]]¶
Solve
-eps Laplacian V = rhoby Conjugate Gradient.Cell-centered FV discretization with Dirichlet in x (
uLleft,uRright) and homogeneous Neumann in y. Boundary values are folded into the right-hand side so the operator stays SPD.- Parameters:
V (numpy.ndarray of shape (Nx, Ny), order='F') – Initial guess; updated in place with the solution.
rho (numpy.ndarray of shape (Nx, Ny), order='F') – Charge density per cell.
grid (Grid2D) – Cell-centered grid.
eps (float, optional) – Constant permittivity. Default
1.0.uL (float, optional) – Dirichlet values on the x-boundaries. Default
0.0.uR (float, optional) – Dirichlet values on the x-boundaries. Default
0.0.tol (float, optional) – Stopping criterion on relative residual. Default
1e-8.max_iter (int, optional) – Iteration cap. Default
10_000.use_preconditioner (bool, optional) – If True, use Jacobi PCG. Marginal gain unless
epsvaries spatially. Default False.record_history (bool, optional) – If True, record
||r|| / ||b||at every iteration. Default False.
- Returns:
report (CGReport) – Iteration count and final relative residual.
history (list of float) – Per-iteration relative residuals when
record_history=True, otherwise an empty list.
Spectral (FFTW)¶
- class poisson_cpp.DSTSolver1D¶
Bases:
pybind11_objectSpectral 1D Poisson solver via DST-I (FFTW), homogeneous Dirichlet BC.
O(N log N): forward DST, divide by eigenvalues, backward DST. The FFTW plan is built once at construction.- Parameters:
- solve(self: poisson_cpp._core.DSTSolver1D, rho: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']¶
Solve
eps0 V'' = -rhoon theNinterior nodes.- Parameters:
rho (numpy.ndarray of shape (N,)) – Charge density at the interior nodes.
- Returns:
Potential
Vat the interior nodes (boundary values are zero).- Return type:
numpy.ndarray of shape (N,)
- class poisson_cpp.DSTSolver2D¶
Bases:
pybind11_objectSpectral 2D Poisson solver via DST-I (FFTW), homogeneous Dirichlet on all four faces.
O(N^2 log N)via tensor product of 1D DSTs.- Parameters:
- solve(self: poisson_cpp._core.DSTSolver2D, rho: Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]', 'flags.f_contiguous']) Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]']¶
Solve
-eps0 Laplacian V = rhoon the interior grid.- Parameters:
rho (numpy.ndarray of shape (Nx, Ny)) – Charge density at the interior nodes.
- Returns:
Potential
Vat the interior nodes (boundary values are zero).- Return type:
numpy.ndarray of shape (Nx, Ny)
AMR quadtree¶
- class poisson_cpp.Quadtree¶
Bases:
pybind11_objectCell-centered quadtree on the square domain
[0, L] x [0, L]with Morton-encoded keys. Only leaves are stored.- Parameters:
- L(self: poisson_cpp._core.Quadtree) float¶
Domain side length.
- at(self: poisson_cpp._core.Quadtree, key: SupportsInt | SupportsIndex) poisson_cpp._core.Cell¶
Copy of the leaf’s
Celldata.
- balance_2to1(self: poisson_cpp._core.Quadtree) None¶
Enforce the 2:1 balance constraint on the current tree.
- build(self: poisson_cpp._core.Quadtree, predicate: collections.abc.Callable[[SupportsInt | SupportsIndex], bool], level_max: SupportsInt | SupportsIndex, rho_func: collections.abc.Callable[[SupportsFloat | SupportsIndex, SupportsFloat | SupportsIndex], float]) None¶
Refine until every leaf either fails
predicateor hitslevel_max, enforce 2:1 balance, then evaluaterho_func(x, y)at every leaf center and store it inCell.rho.
- cell_center(self: poisson_cpp._core.Quadtree, key: SupportsInt | SupportsIndex) tuple[float, float]¶
Geometric center
(x, y)of the given cell key.
- cell_size(self: poisson_cpp._core.Quadtree, level: SupportsInt | SupportsIndex) float¶
Side length of a cell at the given level.
- is_leaf(self: poisson_cpp._core.Quadtree, key: SupportsInt | SupportsIndex) bool¶
Return True if the cell key is currently a leaf.
- leaves(self: poisson_cpp._core.Quadtree) dict[int, poisson_cpp._core.Cell]¶
Return
{key: Cell}for all leaves (copy).
- level_min(self: poisson_cpp._core.Quadtree) int¶
Initial uniform refinement level.
- neighbour_leaves(self: poisson_cpp._core.Quadtree, key: SupportsInt | SupportsIndex, dir: poisson_cpp._core.Direction) list[int]¶
Return up to 2 leaf neighbours of
keyacross facedir.
- num_leaves(self: poisson_cpp._core.Quadtree) int¶
Total number of leaves in the tree.
- refine(self: poisson_cpp._core.Quadtree, key: SupportsInt | SupportsIndex) None¶
Subdivide the leaf
keyinto its 4 children.keymust currently be a leaf.
- class poisson_cpp.Cell¶
Bases:
pybind11_objectPer-leaf data: potential
Vand charge densityrho.- property V¶
Potential at the leaf.
- property rho¶
Charge density at the leaf.
- class poisson_cpp.Direction¶
Bases:
pybind11_objectGeometric direction across a face of a quadtree cell.
Members:
N
S
E
W
- E = <Direction.E: 2>¶
- N = <Direction.N: 0>¶
- S = <Direction.S: 1>¶
- W = <Direction.W: 3>¶
- Direction.name -> str
- property value¶
- class poisson_cpp.AMRArrays¶
Bases:
pybind11_objectFlat array view of a Quadtree, suitable for fast iterative solves.
Built by
extract_arrays().Vandrhoare mutable; modify them then callwriteback()to push results back into the tree.- property V¶
Potential per leaf.
- property Vc¶
Diagonal coefficient of the FV stencil per leaf.
- property h¶
Cell size per leaf.
- property keys¶
Morton key of each leaf, ordered to match V/rho/h.
- property nb0¶
First neighbour index per face (-1 if absent).
- property nb1¶
Second neighbour index per face (-1 if absent).
- property rho¶
Charge density per leaf.
- property w0¶
Off-diagonal weight for
nb0per face.
- property w1¶
Off-diagonal weight for
nb1per face.
- poisson_cpp.extract_arrays(tree: poisson_cpp._core.Quadtree) poisson_cpp._core.AMRArrays¶
Build an
AMRArraysview from a fully builtQuadtree.The returned view stores precomputed FV stencil weights handling 2:1 coarse/fine interfaces and Dirichlet (V=0) ghost cells at the boundary.
- poisson_cpp.writeback(tree: poisson_cpp._core.Quadtree, keys: collections.abc.Sequence[SupportsInt | SupportsIndex], V: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']) None¶
Push the flat
Varray back into the tree leaves.- Parameters:
tree (Quadtree) – Tree to update in place.
keys (list[int]) – Leaf keys, must match the order of
V(usearr.keys).V (numpy.ndarray of shape (n_leaves,))
- poisson_cpp.amr_sor(arr: poisson_cpp._core.AMRArrays, omega: SupportsFloat | SupportsIndex = 1.85, tol: SupportsFloat | SupportsIndex = 1e-08, max_iter: SupportsInt | SupportsIndex = 20000, eps0: SupportsFloat | SupportsIndex = 1.0) poisson_cpp._core.AMRSORReport¶
Run Gauss-Seidel-like SOR on the AMR arrays, in place.
- Parameters:
arr (AMRArrays) – Flat view from
extract_arrays().arr.Vis updated in place.omega (float, optional) – Relaxation factor. Default
1.85.tol (float, optional) – Convergence threshold. Default
1e-8.max_iter (int, optional) – Iteration cap. Default
20_000.eps0 (float, optional) – Permittivity factor. Default
1.0.
- Return type:
- poisson_cpp.amr_residual(arr: poisson_cpp._core.AMRArrays, eps0: SupportsFloat | SupportsIndex = 1.0) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']¶
FV residual at every leaf:
r_i = sum w * V_neigh + h_i^2 rho_i / eps0 - Vc_i V_i.- Parameters:
- Return type:
numpy.ndarray of shape (n_leaves,)
- class poisson_cpp.AMRSORParams¶
Bases:
pybind11_objectTuning parameters for the AMR SOR solver.
- property eps0¶
Permittivity factor in
rhs = h^2 rho / eps0.
- property max_iter¶
Iteration cap.
- property omega¶
Relaxation factor.
- property tol¶
Stopping criterion on max-norm of the per-iteration update.
- class poisson_cpp.AMRSORReport¶
Bases:
pybind11_objectResult of an AMR SOR run.
- property iterations¶
Number of sweeps performed.
- property residual¶
Final
||V_new - V||_inf.
- poisson_cpp.make_key(level: SupportsInt | SupportsIndex, i: SupportsInt | SupportsIndex, j: SupportsInt | SupportsIndex) int¶
Encode
(level, i, j)into a Morton CellKey.
- poisson_cpp.level_of(key: SupportsInt | SupportsIndex) int¶
Refinement level encoded in the CellKey.
- poisson_cpp.i_of(key: SupportsInt | SupportsIndex) int¶
x-coordinate (cell index at the cell’s level) encoded in the key.
- poisson_cpp.j_of(key: SupportsInt | SupportsIndex) int¶
y-coordinate (cell index at the cell’s level) encoded in the key.
Multigrille¶
- poisson_cpp.gs_smooth(V: Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]', 'flags.writeable', 'flags.f_contiguous'], rho: Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]', 'flags.f_contiguous'], h: SupportsFloat | SupportsIndex, n_iter: SupportsInt | SupportsIndex) None¶
Gauss-Seidel red-black smoothing sweeps on the uniform cell-centered FV Poisson operator with Dirichlet
V = 0on the full boundary.- Parameters:
V (numpy.ndarray of shape (N, N), order='F') – Initial guess; updated in place after
n_itersweeps.rho (numpy.ndarray of shape (N, N), order='F') – Right-hand side.
h (float) – Cell size.
n_iter (int) – Number of sweeps.
- poisson_cpp.laplacian_fv(V: Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]', 'flags.f_contiguous'], h: SupportsFloat | SupportsIndex) Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]']¶
Apply the uniform FV Laplacian: returns
A VwhereAis the 5-point stencil with DirichletV = 0on the boundary.- Parameters:
V (numpy.ndarray of shape (N, N))
h (float) – Cell size.
- Return type:
numpy.ndarray of shape (N, N)
- poisson_cpp.restrict_avg(r: Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]', 'flags.f_contiguous']) Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]']¶
4-cell averaging restriction:
(N, N) -> (N/2, N/2).Nmust be even.- Parameters:
r (numpy.ndarray of shape (N, N))
- Return type:
numpy.ndarray of shape (N/2, N/2)
- poisson_cpp.prolongate_const(c: Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]', 'flags.f_contiguous']) Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]']¶
Piecewise-constant prolongation (order 0):
(N/2, N/2) -> (N, N).- Parameters:
c (numpy.ndarray of shape (M, M))
- Return type:
numpy.ndarray of shape (2M, 2M)
- poisson_cpp.prolongate_bilinear(c: Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]', 'flags.f_contiguous']) Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]']¶
Bilinear prolongation (order 2):
(M, M) -> (2M, 2M)for cell-centered FV. Each fine cell gets a weighted combination of the enclosing coarse cell and its two neighbours along the fine cell’s offset.- Parameters:
c (numpy.ndarray of shape (M, M))
- Return type:
numpy.ndarray of shape (2M, 2M)
- poisson_cpp.vcycle_uniform(V: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, n]'], rho: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, n]'], h: SupportsFloat | SupportsIndex, n_pre: SupportsInt | SupportsIndex = 2, n_post: SupportsInt | SupportsIndex = 2, n_min: SupportsInt | SupportsIndex = 4) Annotated[numpy.typing.NDArray[numpy.float64], '[m, n]']¶
Recursive V-cycle multigrid on a uniform grid with Dirichlet
V = 0on the boundary. Coarsens by factor 2 down to size<= n_min.- Parameters:
V (numpy.ndarray of shape (N, N)) – Initial guess.
rho (numpy.ndarray of shape (N, N)) – Right-hand side.
h (float) – Cell size at the finest level.
n_pre (int, optional) – Pre/post Gauss-Seidel sweeps at each level. Default
2.n_post (int, optional) – Pre/post Gauss-Seidel sweeps at each level. Default
2.n_min (int, optional) – Coarsest grid size at which to stop recursing. Default
4.
- Returns:
Updated potential.
- Return type:
numpy.ndarray of shape (N, N)
- poisson_cpp.vcycle_amr_composite(arr: poisson_cpp._core.AMRArrays, tree: poisson_cpp._core.Quadtree, params: poisson_cpp._core.CompositeParams = <poisson_cpp._core.CompositeParams object at 0x7f708b946330>) None¶
One composite 2-grid V-cycle on an AMR tree.
- Steps:
SOR pre-smoothing on AMR leaves.
AMR residual
r.Volume-weighted restriction
r -> r_con a uniform2^level_mingrid.n_coarse_cyclesuniform V-cycles to approximately solveA delta = r_c.Bilinear prolongation of
deltaback to AMR leaves (V += delta).SOR post-smoothing on AMR leaves.
The AMR array state is updated in place.
- Parameters:
arr (AMRArrays) – Flat view, updated in place.
tree (Quadtree) – The tree from which
arrwas extracted.params (CompositeParams, optional) – Tuning parameters; defaults match the C++ defaults.
- class poisson_cpp.CompositeParams¶
Bases:
pybind11_objectTuning parameters for the 2-grid composite V-cycle on AMR.
- property eps0¶
Permittivity factor.
- property n_coarse_cycles¶
Number of uniform V-cycles on the coarse problem.
- property n_post¶
Post-smoothing SOR sweeps on AMR.
- property n_pre¶
Pre-smoothing SOR sweeps on AMR.
- property omega¶
SOR omega for AMR smoothing.
Utilitaires¶
- poisson_cpp.harmonic_mean(a: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]'], b: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']¶
Element-wise harmonic mean
2 a b / (a + b).Used to compute the effective permittivity at the face between two cells with different
eps_r(preserves continuity of the normal component ofD = eps0 eps_r grad V).- Parameters:
a (numpy.ndarray of shape (N,)) – Per-element values, must be strictly positive.
b (numpy.ndarray of shape (N,)) – Per-element values, must be strictly positive.
- Return type:
numpy.ndarray of shape (N,)
- poisson_cpp.fftw_install_hint()[source]¶
Return the FFTW install command for the current platform.
DSTSolver1D/2D need FFTW3 at build time. After installing FFTW with the command below, reinstall poisson-cpp from source so CMake picks it up.
- Return type:
- poisson_cpp.dump_amr_snapshot(tree, path, extra=None)[source]¶
Serialize an AMR
Quadtreeto a JSON file.The schema matches the one written by the C++ CLI
poisson_demo --problem amr --output ...: a list of leaves with their geometry (key,level,x,y,h) and per-leaf scientific data (V,rho).- Parameters:
- Returns:
Absolute path of the written file.
- Return type:
- poisson_cpp.has_fftw3 = True¶
bool(x) -> bool
Returns True when the argument x is true, False otherwise. The builtins True and False are the only two instances of the class bool. The class bool is a subclass of the class int, and cannot be subclassed.
- poisson_cpp.__version__ = '0.2.0'¶
str(object=’’) -> str str(bytes_or_buffer[, encoding[, errors]]) -> str
Create a new string object from the given object. If encoding or errors is specified, then the object must expose a data buffer that will be decoded using the given encoding and error handler. Otherwise, returns the result of object.__str__() (if defined) or repr(object). encoding defaults to sys.getdefaultencoding(). errors defaults to ‘strict’.