Welcome to GASFIR’s documentation!

General Approximator for Strong-Field Ionization Rates

GASFIR is the Python implementation of the ionization-rate retrieval framework introduced in:

Agarwal, Scrinzi & Yakovlev,
General approximator for strong-field ionization rates,
Physical Review A 113, L021101 (2026).

It addresses the long-standing problem of determining accurate, time-resolved ionization rates for atoms in strong laser fields — a quantity fundamental to attosecond science. By fitting a five-parameter kernel to ionization probabilities of few-cycle laser pulses it retrieves sub-optical-cycle ionization dynamics within and beyond the strong-field approximation.

Three calculation methods are provided:

  1. GASFIR — semi-analytic kernel (Eqs. 2, 3, 7 of the paper)

  2. Exact SFA — full numerical saddle-point integration of the SFA kernel

  3. QS — quasistatic (tunneling) limit; connects to ADK/PPT theory

GASFIR Logo

Quick Start

Install GASFIR:

pip install gasfir

Browse available materials

Call get_parameters() with no argument to print all built-in parameter sets:

from gasfir import get_parameters

get_parameters()  # prints a table of all available material keys

Ionization probability (GASFIR kernel)

from gasfir import create_pulse, get_diabatic_ionization_probability, get_parameters

# Create a laser pulse: 800 nm, 1e14 W/cm², CEP=0, 6 optical cycles
laser = create_pulse(800, 1e14, 0, 6)

# Load GASFIR parameters for hydrogen
params = get_parameters("Hydrogen_SFA")

# Compute total ionization probability
prob = get_diabatic_ionization_probability(pulse=laser, param_dict=params)
print(f"Ionization probability: {prob:.6f}")

Time-resolved ionization rate

from gasfir import create_pulse, get_diabatic_ionization_rate, get_parameters
import numpy as np

laser = create_pulse(800, 1e14, 0, 6)
params = get_parameters("Hydrogen_SFA")

t_grid = laser.get_tgrid(dt=1.0)        # time grid in atomic units
rate = get_diabatic_ionization_rate(t_grid, laser, params)
# rate[i] is the instantaneous ionization rate at t_grid[i]

Quasi-static (tunneling) rate

from gasfir import create_pulse, get_quasi_static_rate_for_field, get_parameters

laser = create_pulse(800, 1e14, 0, 6)
params = get_parameters("Hydrogen_SFA")

t_grid = laser.get_tgrid(dt=1.0)
field_values = laser.get_electric_field(t_grid)   # E(t) array in atomic units
qs_rate = get_quasi_static_rate_for_field(field_values, params)

Exact SFA calculation

The exact SFA kernel is selected with kernel_type="exact_SFA" on the same probability function — it is not chosen by swapping the parameter set. The exact SFA kernel is currently implemented for the hydrogen ground state only; use it together with the H_SFA parameters. For hydrogen the GASFIR kernel with H_SFA reproduces the exact SFA result (GASFIR was fit to it), so the two should agree:

from gasfir import create_pulse, get_diabatic_ionization_probability, get_parameters

laser  = create_pulse(800, 1e14, 0, 6)
params = get_parameters("H_SFA")

# exact SFA (hydrogen only) — note the kernel_type argument
prob_sfa = get_diabatic_ionization_probability(
    pulse=laser, param_dict=params, kernel_type="exact_SFA"
)
# GASFIR kernel (default) with the same H_SFA parameters ≈ exact SFA
prob_gasfir = get_diabatic_ionization_probability(pulse=laser, param_dict=params)

print(f"exact SFA = {prob_sfa:.6f},  GASFIR = {prob_gasfir:.6f}")

Batch processing over intensities

from gasfir import create_pulse, get_diabatic_ionization_probability, get_parameters

params = get_parameters("Hydrogen_SFA")
intensities = [5e13, 1e14, 2e14, 5e14]  # W/cm²

results = []
for I in intensities:
    laser = create_pulse(800, I, 0, 6)
    prob = get_diabatic_ionization_probability(pulse=laser, param_dict=params)
    results.append((I, prob))
    print(f"I = {I:.0e} W/cm²  ->  P = {prob:.6f}")

Features

  • Multiple Ionization Methods: Support for GASFIR, exact SFA, and quasi-static calculations

  • Flexible Laser Field Modeling: Create single or multiple laser pulses with various parameters

  • Pre-built Parameter Sets: Ready-to-use parameters for common atoms and molecules

  • High Performance: Optimized with Numba for fast calculations

  • Comprehensive Documentation: Detailed API documentation and examples

  • Extensive Testing: Thorough test suite with high coverage

Indices and tables