|
numerics
|
#include <tdse_solver.hpp>
Public Member Functions | |
| TDSESolver (int N, double L, double dt) | |
| void | set_potential (Potential p, double param=0.0) |
| Build potential and recompute Thomas factorisations. | |
| void | init_gaussian (double x0, double y0, double kx, double ky, double sigma) |
| Gaussian wavepacket: psi = A*exp(-(r-r0)^2/sigma^2)*exp(i*k*r), then normalised. | |
| void | set_mode (int k) |
| Set psi to the k-th eigenmode (modes must be computed first). | |
| void | step () |
| Advance one full time step (Strang splitting). | |
| void | compute_modes (int k=8) |
| int | at (int i, int j) const |
| double | prob (int i, int j) const |
| double | phase_ang (int i, int j) const |
| double | compute_norm () const |
| integral|psi|^2 dA (Gauss-Legendre over each grid strip). | |
| double | compute_energy () const |
| <psi|H|psi> using 5-point Laplacian and the potential. | |
| void | renormalize () |
Public Attributes | |
| int | N |
| Interior grid points per axis. | |
| double | L |
| Domain length ([0,L]x[0,L]) | |
| double | h |
| Grid spacing = L/(N+1) | |
| double | dt |
| Time step. | |
| num::CVector | psi |
| NxN wavefunction (complex, length N^2) | |
| num::Vector | V |
| NxN potential (real, length N^2) | |
| std::vector< EigenMode > | modes |
| bool | modes_ready = false |
| Stats | stats |
Definition at line 70 of file tdse_solver.hpp.
| tdse::TDSESolver::TDSESolver | ( | int | N, |
| double | L, | ||
| double | dt | ||
| ) |
Definition at line 29 of file tdse_solver.cpp.
References init_gaussian().
|
inline |
Definition at line 107 of file tdse_solver.hpp.
References N.
Referenced by init_gaussian(), main(), phase_ang(), prob(), and set_potential().
| double tdse::TDSESolver::compute_energy | ( | ) | const |
<psi|H|psi> using 5-point Laplacian and the potential.
Definition at line 173 of file tdse_solver.cpp.
References h, num::laplacian_stencil_2d(), N, psi, num::BasicVector< T >::size(), and V.
Referenced by main().
| void tdse::TDSESolver::compute_modes | ( | int | k = 8 | ) |
Run Lanczos to find the k lowest eigenstates of H. For CircularWell also fills EigenMode::exact_energy via Bessel zero search.
Definition at line 195 of file tdse_solver.cpp.
References num::bessel_j(), num::brent(), tdse::CircularWell, tdse::EigenMode::energy, h, L, num::lanczos(), num::laplacian_stencil_2d(), modes, modes_ready, N, tdse::Stats::n_modes, tdse::EigenMode::phi, stats, and V.
Referenced by main().
| double tdse::TDSESolver::compute_norm | ( | ) | const |
integral|psi|^2 dA (Gauss-Legendre over each grid strip).
Definition at line 168 of file tdse_solver.cpp.
References h, num::norm(), and psi.
Referenced by main().
| void tdse::TDSESolver::init_gaussian | ( | double | x0, |
| double | y0, | ||
| double | kx, | ||
| double | ky, | ||
| double | sigma | ||
| ) |
Gaussian wavepacket: psi = A*exp(-(r-r0)^2/sigma^2)*exp(i*k*r), then normalised.
Definition at line 116 of file tdse_solver.cpp.
References at(), h, N, psi, and renormalize().
Referenced by main(), and TDSESolver().
|
inline |
|
inline |
| void tdse::TDSESolver::renormalize | ( | ) |
Definition at line 188 of file tdse_solver.cpp.
References h, num::norm(), psi, and num::scale().
Referenced by init_gaussian().
| void tdse::TDSESolver::set_mode | ( | int | k | ) |
Set psi to the k-th eigenmode (modes must be computed first).
Definition at line 130 of file tdse_solver.cpp.
References modes, modes_ready, N, and psi.
Referenced by main().
| void tdse::TDSESolver::set_potential | ( | Potential | p, |
| double | param = 0.0 |
||
| ) |
Build potential and recompute Thomas factorisations.
Definition at line 44 of file tdse_solver.cpp.
References at(), tdse::Barrier, tdse::CircularWell, tdse::DoubleSlit, h, tdse::Harmonic, L, modes_ready, N, num::scale(), and V.
Referenced by main().
| void tdse::TDSESolver::step | ( | ) |
Advance one full time step (Strang splitting).
Definition at line 153 of file tdse_solver.cpp.
References dt, stats, and tdse::Stats::step_ms.
Referenced by main().
| double tdse::TDSESolver::dt |
| double tdse::TDSESolver::h |
Grid spacing = L/(N+1)
Definition at line 74 of file tdse_solver.hpp.
Referenced by compute_energy(), compute_modes(), compute_norm(), init_gaussian(), renormalize(), and set_potential().
| double tdse::TDSESolver::L |
Domain length ([0,L]x[0,L])
Definition at line 73 of file tdse_solver.hpp.
Referenced by compute_modes(), and set_potential().
| std::vector<EigenMode> tdse::TDSESolver::modes |
Definition at line 80 of file tdse_solver.hpp.
Referenced by compute_modes(), main(), and set_mode().
| bool tdse::TDSESolver::modes_ready = false |
Definition at line 81 of file tdse_solver.hpp.
Referenced by compute_modes(), main(), set_mode(), and set_potential().
| int tdse::TDSESolver::N |
Interior grid points per axis.
Definition at line 72 of file tdse_solver.hpp.
Referenced by at(), compute_energy(), compute_modes(), init_gaussian(), set_mode(), and set_potential().
| num::CVector tdse::TDSESolver::psi |
NxN wavefunction (complex, length N^2)
Definition at line 77 of file tdse_solver.hpp.
Referenced by compute_energy(), compute_norm(), init_gaussian(), phase_ang(), prob(), renormalize(), and set_mode().
| Stats tdse::TDSESolver::stats |
Definition at line 83 of file tdse_solver.hpp.
Referenced by compute_modes(), main(), and step().
| num::Vector tdse::TDSESolver::V |
NxN potential (real, length N^2)
Definition at line 78 of file tdse_solver.hpp.
Referenced by compute_energy(), compute_modes(), main(), and set_potential().