""" 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, # Domain parameters Lx=20.0, # Domain length (centered around 0) Nx=1024, # Number of grid points # Physical parameters hbar=1.0, # Reduced Planck constant mass=1.0, # Particle mass omega=1.0, # Harmonic oscillator frequency # Solver parameters dealias=3/2, # Dealiasing factor stop_sim_time=5.0, # Final simulation time timestep=1e-3, # Time step size timestepper=d3.RK443, # Time integration scheme 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__() # Basic parameter validation 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") # Store domain and grid parameters self.Lx = Lx # Physical domain size [-Lx/2, Lx/2] self.Nx = Nx # Number of Fourier modes for spatial resolution # Store physical parameters (in natural units where convenient) self.hbar = hbar # Reduced Planck constant - sets quantum scale self.mass = mass # Particle mass - affects kinetic energy self.omega = omega # Oscillator frequency - sets energy and time scales # Store numerical solver parameters self.dealias = dealias # Prevents aliasing in spectral methods self.stop_sim_time = stop_sim_time # Total evolution time self.timestep = timestep # Time step (must satisfy CFL condition) self.timestepper = timestepper # Runge-Kutta integration scheme self.dtype = dtype # Floating point precision # Setup Dedalus spectral method infrastructure self.xcoord = d3.Coordinate("x") # Define spatial coordinate self.dist = d3.Distributor(self.xcoord, dtype=dtype) # Handle parallel distribution # Create Fourier basis: periodic functions on [-Lx/2, Lx/2] # RealFourier uses cosines/sines, ideal for real-valued problems self.xbasis = d3.RealFourier(self.xcoord, size=Nx, bounds=(-Lx/2, Lx/2), dealias=dealias) # Generate physical grid points for visualization and analysis 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: # Generate random physical parameters for wave packet diversity # Center position: avoid boundaries to prevent edge effects x0 = np.random.uniform(-self.Lx/4, self.Lx/4) # Wave packet width: balance between localization and numerical stability sigma = np.random.uniform(0.5, 2.0) # Initial momentum: determines kinetic energy and motion direction k0 = np.random.uniform(-2.0, 2.0) # Initial amplitude: will be normalized, but affects relative scales amplitude = np.random.uniform(0.5, 2.0) # Generate complex Gaussian wave packet: ψ(x) = A*exp[-(x-x₀)²/2σ²]*exp(ik₀x) # This represents a localized particle with definite momentum # Use full grid from Dedalus (includes dealiasing padding) x_grid = self.x.ravel() # Get full dealiased grid for initialization psi_init = generate_gaussian_wave_packet( x_grid, x0, sigma, k0, amplitude ) # Normalize wavefunction to satisfy quantum probability condition ∫|ψ|²dx = 1 # This ensures proper quantum state throughout evolution norm = np.sqrt(np.trapezoid(np.abs(psi_init)**2, x_grid)) psi_init /= norm # Solve time-dependent Schrödinger equation and return full solution data # Pass the complex initial condition (splitting happens in solve method) 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. """ # Setup PDE fields for real and imaginary parts of wavefunction # These will store ψᵣ(x,t) and ψᵢ(x,t) at each time step psi_r = self.dist.Field(name="psi_r", bases=self.xbasis) # Real part ψᵣ psi_i = self.dist.Field(name="psi_i", bases=self.xbasis) # Imaginary part ψᵢ # Define spatial derivative operator for kinetic energy term # d/dx operator using spectral differentiation (exact for polynomials) dx = lambda A: d3.Differentiate(A, self.xcoord) # Create harmonic oscillator potential V(x) = ½mω²x² # This confines the quantum particle to oscillate around x=0 V = self.dist.Field(name="V", bases=self.xbasis) V['g'] = 0.5 * self.mass * (self.omega * self.x)**2 # Set grid values # Setup coupled PDE system for real and imaginary parts # # Starting from the Schrödinger equation: iℏ ∂ψ/∂t = Ĥψ # where Ĥ = -ℏ²/2m ∇² + V(x) and ψ = ψᵣ + iψᵢ # # Substituting ψ = ψᵣ + iψᵢ and separating real/imaginary parts: # Real part: ∂ψᵣ/∂t = (ℏ/2m)∇²ψᵢ - V(x)ψᵢ/ℏ # Imaginary part: ∂ψᵢ/∂t = -(ℏ/2m)∇²ψᵣ + V(x)ψᵣ/ℏ # # These equations are coupled: ψᵣ evolution depends on ψᵢ and vice versa # This preserves the unitary evolution and probability conservation # Create namespace with all variables for Dedalus equation parser hbar = self.hbar mass = self.mass namespace = locals() # Include all local variables (hbar, mass, dx, etc.) # Define the initial value problem with two coupled fields problem = d3.IVP([psi_r, psi_i], namespace=namespace) # Add the coupled evolution equations # Note: dx(dx(field)) computes the second derivative ∇²field 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") # Split complex initial condition into real and imaginary parts # This converts ψ = ψᵣ + iψᵢ to two real-valued fields for the PDE solver psi_r_init = np.real(initial_condition) # Real component ψᵣ(x,0) psi_i_init = np.imag(initial_condition) # Imaginary component ψᵢ(x,0) # Set initial conditions in Dedalus fields # ['g'] accessor sets grid point values directly psi_r['g'] = psi_r_init # Initialize real part ψᵣ(x,0) psi_i['g'] = psi_i_init # Initialize imaginary part ψᵢ(x,0) # Build time integrator from the PDE system # This compiles the equations and sets up the Runge-Kutta stepper solver = problem.build_solver(self.timestepper) solver.stop_sim_time = self.stop_sim_time # Set maximum evolution time # Initialize storage arrays for solution snapshots # ['g', 1] gets grid values in proper layout for post-processing psi_r_list = [psi_r['g', 1].copy()] # Store ψᵣ(x,t) snapshots psi_i_list = [psi_i['g', 1].copy()] # Store ψᵢ(x,t) snapshots t_list = [solver.sim_time] # Store corresponding time values energy_list = [] # Store energy for conservation check # Calculate initial total energy E = ⟨ψ|Ĥ|ψ⟩ # This should be conserved throughout the evolution (Hamiltonian is Hermitian) x_1d = self.x.ravel() # Get full grid for energy calculation dx = x_1d[1] - x_1d[0] # Uniform grid spacing for integration energy_initial = self._compute_energy(psi_r_init, psi_i_init, V, dx) energy_list.append(energy_initial) # Main time evolution loop # Save snapshots periodically to avoid memory issues while capturing dynamics save_frequency = max(1, int(0.05 / self.timestep)) # Save every ~0.05 time units while solver.proceed: # Continue until stop_sim_time reached # Advance solution by one time step using Runge-Kutta solver.step(self.timestep) # Periodically store snapshots for analysis/visualization if solver.iteration % save_frequency == 0: # Extract current wavefunction components from Dedalus fields psi_r_current = psi_r['g', 1].copy() # Current real part psi_i_current = psi_i['g', 1].copy() # Current imaginary part # Store snapshot data psi_r_list.append(psi_r_current) psi_i_list.append(psi_i_current) t_list.append(solver.sim_time) # Monitor energy conservation (should remain constant) energy = self._compute_energy(psi_r_current, psi_i_current, V, dx) energy_list.append(energy) # Convert to numpy arrays 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) # Compute probability density |ψ|² prob_density = psi_r_trajectory**2 + psi_i_trajectory**2 # Create combined state vector (concatenate real and imaginary parts) state_trajectory = np.concatenate([psi_r_trajectory, psi_i_trajectory], axis=1) # Compute potential energy at grid points V_values = V['g', 1].copy() return { # Coordinates "spatial_coordinates": self.x.ravel(), "time_coordinates": time_coordinates, # Initial conditions "psi_r_initial": psi_r_init, "psi_i_initial": psi_i_init, # Solution trajectories "psi_r_trajectory": psi_r_trajectory, "psi_i_trajectory": psi_i_trajectory, "state_trajectory": state_trajectory, # Combined (2*Nx,) for ML "probability_density": prob_density, # |ψ|² # Physical quantities "potential": V_values, # Harmonic oscillator potential "total_energy": energy_trajectory, # Energy conservation check # Physical parameters "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. """ # Compute kinetic energy: T = -ℏ²/2m ⟨ψ|∇²ψ⟩ # For complex ψ = ψᵣ + iψᵢ: T = -ℏ²/2m ∫[ψᵣ∇²ψᵣ + ψᵢ∇²ψᵢ] dx # Use second-order finite differences for ∇² approximation psi_r_xx = np.gradient(np.gradient(psi_r, dx), dx) # ∂²ψᵣ/∂x² psi_i_xx = np.gradient(np.gradient(psi_i, dx), dx) # ∂²ψᵢ/∂x² # Integrate kinetic energy density over space kinetic = -(self.hbar**2)/(2*self.mass) * np.trapezoid( psi_r * psi_r_xx + psi_i * psi_i_xx, dx=dx ) # Compute potential energy: V = ⟨ψ|V|ψ⟩ = ∫ |ψ(x)|² V(x) dx # Probability density |ψ|² = ψᵣ² + ψᵢ² weighted by potential V_values = V_field['g', 1] # Extract potential values from Dedalus field potential = np.trapezoid((psi_r**2 + psi_i**2) * V_values, dx=dx) # Total energy = kinetic + potential (should be conserved) return kinetic + potential