|
|
""" |
|
|
1D Time-Dependent Schrödinger Equation Dataset with Harmonic Oscillator Potential. |
|
|
|
|
|
Solves the time-dependent Schrödinger equation: |
|
|
iℏ ∂ψ/∂t = Ĥψ |
|
|
where Ĥ = -ℏ²/2m ∇² + V(x) and V(x) = ½mω²x² |
|
|
|
|
|
The complex wavefunction ψ = ψᵣ + iψᵢ is split into real and imaginary parts, |
|
|
which are evolved separately and concatenated into a 2×Nx state vector. |
|
|
|
|
|
Physical applications: quantum mechanics, atom optics, Bose-Einstein condensates. |
|
|
""" |
|
|
|
|
|
import numpy as np |
|
|
from torch.utils.data import IterableDataset |
|
|
|
|
|
import dedalus.public as d3 |
|
|
import logging |
|
|
|
|
|
logger = logging.getLogger(__name__) |
|
|
|
|
|
|
|
|
def generate_gaussian_wave_packet(x, x0, sigma, k0, amplitude=1.0): |
|
|
""" |
|
|
Generate a Gaussian wave packet initial condition for the Schrödinger equation. |
|
|
|
|
|
Creates a localized wave packet: ψ(x) = A * exp[-(x-x₀)²/(2σ²)] * exp(ik₀x) |
|
|
This represents a particle localized around position x₀ with momentum k₀. |
|
|
|
|
|
Args: |
|
|
x: Spatial coordinates array |
|
|
x0: Center position of the wave packet |
|
|
sigma: Width parameter (standard deviation of Gaussian envelope) |
|
|
k0: Initial momentum (wave number) |
|
|
amplitude: Amplitude scaling factor |
|
|
|
|
|
Returns: |
|
|
Complex wavefunction ψ = ψᵣ + iψᵢ as numpy array |
|
|
|
|
|
Notes: |
|
|
The wave packet satisfies the normalization ∫|ψ|²dx = 1 after normalization |
|
|
in the calling code. |
|
|
""" |
|
|
envelope = amplitude * np.exp(-(x - x0)**2 / (2 * sigma**2)) |
|
|
phase = np.exp(1j * k0 * x) |
|
|
return envelope * phase |
|
|
|
|
|
|
|
|
class SchrodingerDataset(IterableDataset): |
|
|
""" |
|
|
Dataset for 1D time-dependent Schrödinger equation with harmonic oscillator potential. |
|
|
|
|
|
This dataset generates solutions to the quantum harmonic oscillator by solving: |
|
|
iℏ ∂ψ/∂t = Ĥψ |
|
|
where the Hamiltonian is: |
|
|
Ĥ = -ℏ²/2m ∇² + V(x) |
|
|
and the potential is: |
|
|
V(x) = ½mω²x² |
|
|
|
|
|
The complex wavefunction ψ = ψᵣ + iψᵢ is split into real and imaginary parts |
|
|
that are evolved separately using coupled PDEs: |
|
|
∂ψᵣ/∂t = (ℏ/2m)∇²ψᵢ - V(x)ψᵢ/ℏ |
|
|
∂ψᵢ/∂t = -(ℏ/2m)∇²ψᵣ + V(x)ψᵣ/ℏ |
|
|
|
|
|
For machine learning applications, the solution is provided as: |
|
|
- Individual trajectories: ψᵣ(x,t) and ψᵢ(x,t) |
|
|
- Combined state vector: [ψᵣ, ψᵢ] with shape (2×Nx,) |
|
|
- Probability density: |ψ|² = ψᵣ² + ψᵢ² |
|
|
- Energy conservation: Total energy ⟨ψ|Ĥ|ψ⟩ |
|
|
|
|
|
Physical Applications: |
|
|
- Quantum mechanics education and visualization |
|
|
- Trapped atom dynamics in optical/magnetic traps |
|
|
- Bose-Einstein condensate mean-field dynamics |
|
|
- Neural operator learning for quantum systems |
|
|
|
|
|
Mathematical Properties: |
|
|
- Unitary evolution (probability conservation) |
|
|
- Energy conservation in the absence of dissipation |
|
|
- Spectral accuracy via Fourier pseudospectral methods |
|
|
""" |
|
|
def __init__( |
|
|
self, |
|
|
|
|
|
Lx=20.0, |
|
|
Nx=1024, |
|
|
|
|
|
|
|
|
hbar=1.0, |
|
|
mass=1.0, |
|
|
omega=1.0, |
|
|
|
|
|
|
|
|
dealias=3/2, |
|
|
stop_sim_time=5.0, |
|
|
timestep=1e-3, |
|
|
timestepper=d3.RK443, |
|
|
dtype=np.float64, |
|
|
): |
|
|
""" |
|
|
Initialize Schrödinger equation dataset with harmonic oscillator potential. |
|
|
|
|
|
Sets up the numerical infrastructure to solve the time-dependent Schrödinger |
|
|
equation using Dedalus spectral methods. Creates Fourier basis functions |
|
|
for spatial derivatives and configures the coupled PDE system for real |
|
|
and imaginary parts of the wavefunction. |
|
|
|
|
|
Mathematical Setup: |
|
|
Domain: x ∈ [-Lx/2, Lx/2] with periodic boundary conditions |
|
|
Grid: Nx Fourier modes with dealiasing for nonlinear terms |
|
|
Time Integration: High-order Runge-Kutta schemes (RK443 recommended) |
|
|
|
|
|
Args: |
|
|
Lx: Spatial domain length. Domain extends from -Lx/2 to +Lx/2. |
|
|
Should be large enough to contain the wave packet throughout |
|
|
the simulation without significant boundary effects. |
|
|
Nx: Number of Fourier modes. Higher values give better spatial |
|
|
resolution. Recommended: 256-1024 for typical simulations. |
|
|
hbar: Reduced Planck constant (ℏ). In natural units, often set to 1. |
|
|
mass: Particle mass (m). Controls kinetic energy scale. |
|
|
omega: Harmonic oscillator frequency (ω). Sets potential energy scale |
|
|
and oscillation period T = 2π/ω. |
|
|
dealias: Dealiasing factor for spectral methods. 3/2 removes aliasing |
|
|
for quadratic nonlinearities (standard for Schrödinger eq). |
|
|
stop_sim_time: Maximum simulation time. Should be several oscillation |
|
|
periods (>> 2π/ω) to see full dynamics. |
|
|
timestep: Time step size. Should satisfy CFL condition for stability. |
|
|
Recommended: dt ≤ 0.01 * (2π/ω) for accuracy. |
|
|
timestepper: Dedalus time integration scheme. RK443 (4th order) |
|
|
provides good accuracy/stability balance. |
|
|
dtype: Floating point precision. np.float64 recommended for accuracy. |
|
|
|
|
|
Raises: |
|
|
ValueError: If domain or grid parameters are invalid. |
|
|
ImportError: If Dedalus is not properly installed. |
|
|
""" |
|
|
super().__init__() |
|
|
|
|
|
|
|
|
if Lx <= 0: |
|
|
raise ValueError("Domain length Lx must be positive") |
|
|
if Nx <= 0 or not isinstance(Nx, int): |
|
|
raise ValueError("Grid points Nx must be a positive integer") |
|
|
if hbar <= 0: |
|
|
raise ValueError("Reduced Planck constant hbar must be positive") |
|
|
if mass <= 0: |
|
|
raise ValueError("Particle mass must be positive") |
|
|
if omega <= 0: |
|
|
raise ValueError("Oscillator frequency omega must be positive") |
|
|
if timestep <= 0: |
|
|
raise ValueError("Time step must be positive") |
|
|
if stop_sim_time <= 0: |
|
|
raise ValueError("Simulation time must be positive") |
|
|
|
|
|
|
|
|
self.Lx = Lx |
|
|
self.Nx = Nx |
|
|
|
|
|
|
|
|
self.hbar = hbar |
|
|
self.mass = mass |
|
|
self.omega = omega |
|
|
|
|
|
|
|
|
self.dealias = dealias |
|
|
self.stop_sim_time = stop_sim_time |
|
|
self.timestep = timestep |
|
|
self.timestepper = timestepper |
|
|
self.dtype = dtype |
|
|
|
|
|
|
|
|
self.xcoord = d3.Coordinate("x") |
|
|
self.dist = d3.Distributor(self.xcoord, dtype=dtype) |
|
|
|
|
|
|
|
|
|
|
|
self.xbasis = d3.RealFourier(self.xcoord, size=Nx, bounds=(-Lx/2, Lx/2), dealias=dealias) |
|
|
|
|
|
|
|
|
self.x = self.dist.local_grid(self.xbasis) |
|
|
|
|
|
def __iter__(self): |
|
|
""" |
|
|
Generate infinite samples from the dataset. |
|
|
|
|
|
Creates diverse quantum wave packet initial conditions by randomly sampling |
|
|
Gaussian wave packet parameters. Each sample represents a different physical |
|
|
scenario with varying localization, momentum, and energy. |
|
|
|
|
|
Initial Condition Generation: |
|
|
1. Random wave packet center x₀ ~ U[-Lx/4, Lx/4] |
|
|
2. Random width σ ~ U[0.5, 2.0] (controls localization) |
|
|
3. Random momentum k₀ ~ U[-2.0, 2.0] (initial velocity) |
|
|
4. Random amplitude A ~ U[0.5, 2.0] (before normalization) |
|
|
5. Normalize to ensure ∫|ψ|²dx = 1 (probability conservation) |
|
|
|
|
|
Physics: |
|
|
- Narrow wave packets (small σ) are highly localized but spread quickly |
|
|
- Wide wave packets (large σ) are delocalized but maintain shape longer |
|
|
- Higher |k₀| gives faster initial motion and higher kinetic energy |
|
|
- All initial conditions are proper quantum states (normalized) |
|
|
|
|
|
Yields: |
|
|
Dict containing complete solution trajectory and physical quantities |
|
|
""" |
|
|
while True: |
|
|
|
|
|
|
|
|
x0 = np.random.uniform(-self.Lx/4, self.Lx/4) |
|
|
|
|
|
|
|
|
sigma = np.random.uniform(0.5, 2.0) |
|
|
|
|
|
|
|
|
k0 = np.random.uniform(-2.0, 2.0) |
|
|
|
|
|
|
|
|
amplitude = np.random.uniform(0.5, 2.0) |
|
|
|
|
|
|
|
|
|
|
|
|
|
|
x_grid = self.x.ravel() |
|
|
psi_init = generate_gaussian_wave_packet( |
|
|
x_grid, x0, sigma, k0, amplitude |
|
|
) |
|
|
|
|
|
|
|
|
|
|
|
norm = np.sqrt(np.trapezoid(np.abs(psi_init)**2, x_grid)) |
|
|
psi_init /= norm |
|
|
|
|
|
|
|
|
|
|
|
yield self.solve(psi_init) |
|
|
|
|
|
def solve(self, initial_condition): |
|
|
""" |
|
|
Solve the time-dependent Schrödinger equation using spectral methods. |
|
|
|
|
|
Integrates the coupled PDE system forward in time using high-order |
|
|
Runge-Kutta methods, storing snapshots of the solution and computing |
|
|
physical quantities like energy for validation. |
|
|
|
|
|
Numerical Method: |
|
|
- Spectral differentiation for spatial derivatives (machine precision) |
|
|
- High-order Runge-Kutta time stepping (4th order RK443) |
|
|
- Adaptive time step control via Dedalus CFL condition |
|
|
- Energy monitoring for conservation verification |
|
|
|
|
|
Args: |
|
|
initial_condition: Complex initial wavefunction ψ(x,0) = ψᵣ + iψᵢ |
|
|
|
|
|
Returns: |
|
|
Dictionary containing: |
|
|
- Solution trajectories: ψᵣ(x,t), ψᵢ(x,t), |ψ(x,t)|² |
|
|
- Coordinate arrays: spatial grid x, time points t |
|
|
- ML-ready state vector: concatenated [ψᵣ, ψᵢ] |
|
|
- Physical quantities: total energy, potential V(x) |
|
|
- System parameters: ℏ, m, ω for reproducibility |
|
|
|
|
|
Note: |
|
|
Energy should be conserved to within numerical precision. |
|
|
Large energy drift indicates time step is too large or |
|
|
simulation time exceeds numerical stability limits. |
|
|
""" |
|
|
|
|
|
|
|
|
psi_r = self.dist.Field(name="psi_r", bases=self.xbasis) |
|
|
psi_i = self.dist.Field(name="psi_i", bases=self.xbasis) |
|
|
|
|
|
|
|
|
|
|
|
dx = lambda A: d3.Differentiate(A, self.xcoord) |
|
|
|
|
|
|
|
|
|
|
|
V = self.dist.Field(name="V", bases=self.xbasis) |
|
|
V['g'] = 0.5 * self.mass * (self.omega * self.x)**2 |
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
hbar = self.hbar |
|
|
mass = self.mass |
|
|
namespace = locals() |
|
|
|
|
|
|
|
|
problem = d3.IVP([psi_r, psi_i], namespace=namespace) |
|
|
|
|
|
|
|
|
|
|
|
problem.add_equation("dt(psi_r) - (hbar/(2*mass))*dx(dx(psi_i)) = - V*psi_i/hbar") |
|
|
problem.add_equation("dt(psi_i) + (hbar/(2*mass))*dx(dx(psi_r)) = V*psi_r/hbar") |
|
|
|
|
|
|
|
|
|
|
|
psi_r_init = np.real(initial_condition) |
|
|
psi_i_init = np.imag(initial_condition) |
|
|
|
|
|
|
|
|
|
|
|
psi_r['g'] = psi_r_init |
|
|
psi_i['g'] = psi_i_init |
|
|
|
|
|
|
|
|
|
|
|
solver = problem.build_solver(self.timestepper) |
|
|
solver.stop_sim_time = self.stop_sim_time |
|
|
|
|
|
|
|
|
|
|
|
psi_r_list = [psi_r['g', 1].copy()] |
|
|
psi_i_list = [psi_i['g', 1].copy()] |
|
|
t_list = [solver.sim_time] |
|
|
energy_list = [] |
|
|
|
|
|
|
|
|
|
|
|
x_1d = self.x.ravel() |
|
|
dx = x_1d[1] - x_1d[0] |
|
|
energy_initial = self._compute_energy(psi_r_init, psi_i_init, V, dx) |
|
|
energy_list.append(energy_initial) |
|
|
|
|
|
|
|
|
|
|
|
save_frequency = max(1, int(0.05 / self.timestep)) |
|
|
|
|
|
while solver.proceed: |
|
|
|
|
|
solver.step(self.timestep) |
|
|
|
|
|
|
|
|
if solver.iteration % save_frequency == 0: |
|
|
|
|
|
psi_r_current = psi_r['g', 1].copy() |
|
|
psi_i_current = psi_i['g', 1].copy() |
|
|
|
|
|
|
|
|
psi_r_list.append(psi_r_current) |
|
|
psi_i_list.append(psi_i_current) |
|
|
t_list.append(solver.sim_time) |
|
|
|
|
|
|
|
|
energy = self._compute_energy(psi_r_current, psi_i_current, V, dx) |
|
|
energy_list.append(energy) |
|
|
|
|
|
|
|
|
psi_r_trajectory = np.array(psi_r_list) |
|
|
psi_i_trajectory = np.array(psi_i_list) |
|
|
time_coordinates = np.array(t_list) |
|
|
energy_trajectory = np.array(energy_list) |
|
|
|
|
|
|
|
|
prob_density = psi_r_trajectory**2 + psi_i_trajectory**2 |
|
|
|
|
|
|
|
|
state_trajectory = np.concatenate([psi_r_trajectory, psi_i_trajectory], axis=1) |
|
|
|
|
|
|
|
|
V_values = V['g', 1].copy() |
|
|
|
|
|
return { |
|
|
|
|
|
"spatial_coordinates": self.x.ravel(), |
|
|
"time_coordinates": time_coordinates, |
|
|
|
|
|
|
|
|
"psi_r_initial": psi_r_init, |
|
|
"psi_i_initial": psi_i_init, |
|
|
|
|
|
|
|
|
"psi_r_trajectory": psi_r_trajectory, |
|
|
"psi_i_trajectory": psi_i_trajectory, |
|
|
"state_trajectory": state_trajectory, |
|
|
"probability_density": prob_density, |
|
|
|
|
|
|
|
|
"potential": V_values, |
|
|
"total_energy": energy_trajectory, |
|
|
|
|
|
|
|
|
"hbar": self.hbar, |
|
|
"mass": self.mass, |
|
|
"omega": self.omega, |
|
|
} |
|
|
|
|
|
def _compute_energy(self, psi_r, psi_i, V_field, dx): |
|
|
""" |
|
|
Compute total energy of the quantum system using wavefunction components. |
|
|
|
|
|
Calculates the expectation value of the Hamiltonian: E = ⟨ψ|Ĥ|ψ⟩ |
|
|
where Ĥ = -ℏ²/2m ∇² + V(x) is the quantum harmonic oscillator Hamiltonian. |
|
|
|
|
|
The energy should be conserved during unitary evolution, making this |
|
|
a crucial diagnostic for numerical accuracy and stability. |
|
|
|
|
|
Mathematical Details: |
|
|
Total Energy: E = T + V |
|
|
Kinetic Energy: T = -ℏ²/2m ∫ ψ*(x) ∇²ψ(x) dx |
|
|
Potential Energy: V = ∫ |ψ(x)|² V(x) dx |
|
|
|
|
|
For split wavefunction ψ = ψᵣ + iψᵢ: |
|
|
T = -ℏ²/2m ∫ [ψᵣ ∇²ψᵣ + ψᵢ ∇²ψᵢ] dx |
|
|
V = ∫ [ψᵣ² + ψᵢ²] V(x) dx |
|
|
|
|
|
Args: |
|
|
psi_r: Real part of wavefunction at current time |
|
|
psi_i: Imaginary part of wavefunction at current time |
|
|
V_field: Potential field from Dedalus |
|
|
dx: Spatial grid spacing for numerical integration |
|
|
|
|
|
Returns: |
|
|
Total energy (float): Should be conserved throughout evolution |
|
|
|
|
|
Note: |
|
|
Uses finite differences for spatial derivatives. For high accuracy, |
|
|
the Dedalus spectral derivatives could be used instead, but this |
|
|
simple approach is sufficient for energy monitoring. |
|
|
""" |
|
|
|
|
|
|
|
|
|
|
|
psi_r_xx = np.gradient(np.gradient(psi_r, dx), dx) |
|
|
psi_i_xx = np.gradient(np.gradient(psi_i, dx), dx) |
|
|
|
|
|
|
|
|
kinetic = -(self.hbar**2)/(2*self.mass) * np.trapezoid( |
|
|
psi_r * psi_r_xx + psi_i * psi_i_xx, dx=dx |
|
|
) |
|
|
|
|
|
|
|
|
|
|
|
V_values = V_field['g', 1] |
|
|
potential = np.trapezoid((psi_r**2 + psi_i**2) * V_values, dx=dx) |
|
|
|
|
|
|
|
|
return kinetic + potential |
|
|
|