|
numerics 0.1.0
|
Functions | |
| void | diffusion_step_2d (Vector &u, int N, double coeff, Backend b=best_backend) |
| void | diffusion_step_2d_dirichlet (Vector &u, int N, double coeff, Backend b=best_backend) |
| void | diffusion_step_2d_4th_dirichlet (Vector &u, int N, double coeff, Backend b=best_backend) |
| void | diffusion_step_2d_dirichlet (ScalarField2D &g, double coeff, Backend b=best_backend) |
| void | diffusion_step_2d_4th_dirichlet (ScalarField2D &g, double coeff, Backend b=best_backend) |
| SparseMatrix | laplacian_sparse_2d (int N) |
| SparseMatrix | backward_euler_matrix (int N, double coeff) |
| SparseMatrix | backward_euler_matrix (const Grid2D &grid, double coeff) |
| Grid2D overload – N is read from the grid. | |
| LinearSolver | make_cg_solver (const SparseMatrix &A, real tol=1e-6) |
| Matrix | poisson2d_fd (const Matrix &f, int N) |
| Solve -Delta u = f on (0,1)^2 (Dirichlet) via DST with FD eigenvalues. | |
| Matrix | poisson2d (const Matrix &f, int N) |
| Solve -Delta u = f on (0,1)^2 (Dirichlet) via DST with exact eigenvalues. | |
|
inline |
Grid2D overload – N is read from the grid.
Definition at line 109 of file diffusion.hpp.
References backward_euler_matrix(), and num::Grid2D::N.
|
inline |
Build the backward Euler system matrix A = I - coeff*L. Solve A * u^{n+1} = u^n at each time step.
Definition at line 90 of file diffusion.hpp.
References num::SparseMatrix::from_triplets().
Referenced by backward_euler_matrix().
|
inline |
Explicit Euler diffusion step on a periodic 2D grid: u += coeff * Lap_periodic(u)
Definition at line 31 of file diffusion.hpp.
References num::axpy(), num::laplacian_stencil_2d_periodic(), and num::BasicVector< T >::size().
|
inline |
Definition at line 63 of file diffusion.hpp.
References diffusion_step_2d_4th_dirichlet(), num::ScalarField2D::N(), and num::ScalarField2D::vec().
|
inline |
Explicit Euler diffusion step with Dirichlet BCs and 4th-order Laplacian. Stability requires coeff <= ~0.08 (vs 0.25 for 2nd order).
Definition at line 49 of file diffusion.hpp.
References num::axpy(), num::laplacian_stencil_2d_4th(), and num::BasicVector< T >::size().
Referenced by diffusion_step_2d_4th_dirichlet().
|
inline |
Definition at line 58 of file diffusion.hpp.
References diffusion_step_2d_dirichlet(), num::ScalarField2D::N(), and num::ScalarField2D::vec().
|
inline |
Explicit Euler diffusion step with Dirichlet (zero) BCs: u += coeff * Lap_dirichlet(u)
Definition at line 40 of file diffusion.hpp.
References num::axpy(), num::laplacian_stencil_2d(), and num::BasicVector< T >::size().
Referenced by diffusion_step_2d_dirichlet().
|
inline |
Build the N^2 x N^2 sparse 5-point Laplacian matrix (Dirichlet BCs). Entry (k,k) = -4, off-diagonals = +1 for each interior neighbor.
Definition at line 70 of file diffusion.hpp.
References num::SparseMatrix::from_triplets().
|
inline |
Returns a LinearSolver that solves A*x = rhs using matrix-free CG. Captures A by reference – A must outlive the returned callable.
Definition at line 115 of file diffusion.hpp.
References num::cg_matfree(), make_cg_solver(), num::matvec(), and num::sparse_matvec().
Referenced by make_cg_solver().
Solve -Delta u = f on (0,1)^2 (Dirichlet) via DST with exact eigenvalues.
Replaces the FD eigenvalue D_k = 2(1 - cos(k*pi/(N+1))) with the exact eigenvalue (k*pi)^2 of the continuous operator -d^2/dx^2. The error is determined by the quadrature approximation of the DST projection rather than the FD truncation of the eigenvalue, giving spectral accuracy.
| f | N x N matrix of RHS values. |
| N | Grid dimension. Must satisfy N = 2^p - 1. |
Definition at line 145 of file poisson.cpp.
References num::pi.
Solve -Delta u = f on (0,1)^2 (Dirichlet) via DST with FD eigenvalues.
| f | N x N matrix of RHS values on the interior grid (row i, col j corresponds to the grid point ((i+1)*h, (j+1)*h)). |
| N | Grid dimension. Must satisfy N = 2^p - 1. |
Definition at line 106 of file poisson.cpp.
References num::pi.