solvers ¤
Root finding and transient solvers.
Modules:
| Name | Description |
|---|---|
ac_sweep | AC small-signal frequency sweep returning S-parameters. |
assembly | Assembly functions for the transient circuit solver. |
circuit_diffeq | Stripped-down ODE integrator for circuit simulation. |
harmonic_balance | Harmonic Balance solver for periodic steady-state circuit analysis. |
linear | Circuit linear solver strategies. |
transient | Transient solvers to be used with Diffrax. |
Classes:
| Name | Description |
|---|---|
BDF2FactorizedTransientSolver | BDF2 upgrade of :class: |
BDF2RefactoringTransientSolver | BDF2 upgrade of :class: |
BDF2VectorizedTransientSolver | BDF2 upgrade of :class: |
CircuitLinearSolver | Abstract base for all circuit linear solvers. |
DenseSolver | Solves the system using dense matrix factorization (LU). |
FactorizedTransientSolver | Transient solver using a Modified Newton (frozen-Jacobian) scheme. |
KLUSolver | Solves the system using the KLU sparse solver (via |
KLUSplitLinear | KLU split solver paired with Modified Newton (frozen-Jacobian) for linear convergence. |
KLUSplitQuadratic | KLU split solver paired with full Newton for quadratic convergence via |
RefactoringTransientSolver | Transient solver with full Newton (quadratic) convergence using |
SDIRK3FactorizedTransientSolver | 3rd-order A-stable SDIRK3 solver with frozen-Jacobian Newton across all stages. |
SDIRK3RefactoringTransientSolver | 3rd-order A-stable SDIRK3 solver with KLU refactor at each Newton iteration. |
SDIRK3VectorizedTransientSolver | 3rd-order A-stable SDIRK3 solver using full Newton-Raphson at each stage. |
SparseSolver | Solves the system using JAX's Iterative BiCGStab solver. |
VectorizedTransientSolver | Transient solver that works strictly on FLAT (Real) vectors. |
Functions:
| Name | Description |
|---|---|
analyze_circuit | Initializes a linear solver strategy for circuit analysis. |
assemble_system_complex | Assemble the residual vectors and effective Jacobian values for an unrolled complex system. |
assemble_system_real | Assemble the residual vectors and effective Jacobian values for a real system. |
circuit_diffeqsolve | Stripped-down :func: |
setup_ac_sweep | Configure and return a callable for AC small-signal S-parameter sweep. |
setup_harmonic_balance | Configure and return a Harmonic Balance callable. |
setup_transient | Configures and returns a function for executing transient analysis. |
BDF2FactorizedTransientSolver ¤
Bases: FactorizedTransientSolver
BDF2 upgrade of :class:FactorizedTransientSolver (frozen-Jacobian Newton).
Factors J_eff once at the predictor state and reuses it across all Newton iterations, trading quadratic convergence for cheaper per-iteration cost. The BDF2 Jacobian scaling α₀/h is used when factoring.
BDF2RefactoringTransientSolver ¤
Bases: RefactoringTransientSolver
BDF2 upgrade of :class:RefactoringTransientSolver (KLU refactor per iteration).
Full quadratic Newton convergence via klujax.refactor at each iteration, combined with BDF2 time discretisation for 2nd-order accuracy.
BDF2VectorizedTransientSolver ¤
Bases: VectorizedTransientSolver
BDF2 upgrade of :class:VectorizedTransientSolver.
Implements variable-step BDF2 via the companion method. On the first step Backward Euler is used automatically; from step 2 onward BDF2 is activated. The Jacobian scaling changes from 1/h (BE) to α₀/h (BDF2) where α₀ = (1 + 2ω)/(1 + ω) and ω = h_n/h_{n-1}.
solver_state is a 3-tuple (y_nm1, h_nm1, q_nm1). h_nm1 is initialised to +inf so that ω = 0 on the first step, making the BDF2 formula reduce to Backward Euler via IEEE 754 arithmetic (no branching). q_nm1 caches Q(y_{n-1}) to avoid recomputing it each step.
CircuitLinearSolver ¤
Bases: AbstractLinearSolver
Abstract base for all circuit linear solvers.
Attributes:
| Name | Type | Description |
|---|---|---|
ground_indices | Array | Indices of nodes connected to ground (forced to 0V). |
is_complex | bool | If True, the system is 2N×2N (real/imag unrolled); otherwise N×N. |
Methods:
| Name | Description |
|---|---|
assume_full_rank | Indicate if the solver assumes the operator is full rank. |
compute | Satisfies the lineax API; call |
init | Initialize the solver state (no-op for stateless solvers). |
solve_dc | DC operating point via damped Newton-Raphson. |
solve_dc_auto | DC Operating Point with automatic homotopy fallback. |
solve_dc_checked | DC Operating Point with convergence status. |
solve_dc_gmin | DC Operating Point via GMIN stepping (homotopy rescue). |
solve_dc_source | DC Operating Point via source stepping (homotopy rescue). |
compute ¤
Satisfies the lineax API; call _solve_impl directly for internal use.
init ¤
solve_dc ¤
solve_dc(
component_groups: dict[str, Any],
y_guess: Array,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC operating point via damped Newton-Raphson.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess vector (shape | required |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector. |
Source code in circulax/solvers/linear.py
solve_dc_auto ¤
solve_dc_auto(
component_groups: dict[str, Any],
y_guess: Array,
g_start: float = 0.01,
n_gmin: int = 10,
n_source: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point with automatic homotopy fallback.
Attempts a direct Newton solve first. If it fails to converge, falls back to GMIN stepping followed by source stepping — all inside a single JIT-compiled kernel via jax.lax.cond.
Strategy: 1. _run_newton — plain damped Newton from y_guess. 2. On failure: solve_dc_gmin (GMIN stepping) starting from y_guess, then solve_dc_source (source stepping) from the GMIN result.
Because jax.lax.cond evaluates both branches at trace time but only one at runtime, this compiles to a single kernel with no Python- level branching.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
g_start | float | Starting leakage for GMIN stepping (rescue branch). | 0.01 |
n_gmin | int | Number of GMIN steps in the rescue branch. | 10 |
n_source | int | Number of source steps in the rescue branch. | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector. |
Source code in circulax/solvers/linear.py
solve_dc_checked ¤
solve_dc_checked(
component_groups: dict[str, Any],
y_guess: Array,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> tuple[Array, Array]
DC Operating Point with convergence status.
Identical to :meth:solve_dc but additionally returns a boolean JAX scalar indicating whether the Newton-Raphson fixed-point iteration reported success. Because the flag is a JAX array (not a Python bool) it can be consumed inside compiled programs:
.. code-block:: python
y, converged = solver.solve_dc_checked(groups, y0)
# Outside JIT — inspect in Python:
if not converged:
y = solver.solve_dc_gmin(groups, y0)
# Inside JIT — branch without Python-level control flow:
y = jax.lax.cond(converged, lambda: y, lambda: solver.solve_dc_gmin(groups, y0))
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess vector (shape | required |
rtol | float | Relative tolerance for | 1e-06 |
atol | float | Absolute tolerance for | 1e-06 |
max_steps | int | Maximum Newton iterations. | 100 |
Returns:
| Type | Description |
|---|---|
tuple[Array, Array] |
|
Source code in circulax/solvers/linear.py
solve_dc_gmin ¤
solve_dc_gmin(
component_groups: dict[str, Any],
y_guess: Array,
g_start: float = 0.01,
n_steps: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point via GMIN stepping (homotopy rescue).
Steps the diagonal regularisation conductance g_leak logarithmically from g_start down to self.g_leak, using each converged solution as the warm start for the next step. The large initial g_leak linearises highly nonlinear components (diodes above threshold, lasers) that would otherwise cause Newton to diverge from a flat 0V start.
Implemented with jax.lax.scan — fully JIT/grad/vmap-compatible.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
g_start | float | Starting leakage conductance (large value, e.g. | 0.01 |
n_steps | int | Number of log-uniform steps from | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector after the full GMIN schedule. |
Source code in circulax/solvers/linear.py
solve_dc_source ¤
solve_dc_source(
component_groups: dict[str, Any],
y_guess: Array,
n_steps: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point via source stepping (homotopy rescue).
Ramps all source amplitudes (components tagged with amplitude_param) from 10 % to 100 % of their netlist values, using each converged solution as the warm start for the next step. This guides Newton through the nonlinear region without the large initial step from 0V to full excitation.
Implemented with jax.lax.scan — fully JIT/grad/vmap-compatible.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
n_steps | int | Number of uniformly-spaced steps from 0.1 to 1.0. | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector at full source amplitude. |
Source code in circulax/solvers/linear.py
DenseSolver ¤
Bases: CircuitLinearSolver
Solves the system using dense matrix factorization (LU).
Best For
- Small to Medium circuits (N < 2000).
- Wavelength sweeps (AC Analysis) on GPU.
- Systems where VMAP parallelism is critical.
Attributes:
| Name | Type | Description |
|---|---|---|
static_rows | Array | Row indices for placing values into dense matrix. |
static_cols | Array | Column indices. |
g_leak | float | Leakage conductance added to diagonal to prevent singularity. |
Methods:
| Name | Description |
|---|---|
assume_full_rank | Indicate if the solver assumes the operator is full rank. |
compute | Satisfies the lineax API; call |
from_component_groups | Factory method to pre-calculate indices for the dense matrix. |
init | Initialize the solver state (no-op for stateless solvers). |
solve_dc | DC operating point via damped Newton-Raphson. |
solve_dc_auto | DC Operating Point with automatic homotopy fallback. |
solve_dc_checked | DC Operating Point with convergence status. |
solve_dc_gmin | DC Operating Point via GMIN stepping (homotopy rescue). |
solve_dc_source | DC Operating Point via source stepping (homotopy rescue). |
compute ¤
Satisfies the lineax API; call _solve_impl directly for internal use.
from_component_groups classmethod ¤
from_component_groups(
component_groups: dict[str, Any],
num_vars: int,
*,
is_complex: bool = False,
g_leak: float = 1e-09
) -> DenseSolver
Factory method to pre-calculate indices for the dense matrix.
Source code in circulax/solvers/linear.py
init ¤
solve_dc ¤
solve_dc(
component_groups: dict[str, Any],
y_guess: Array,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC operating point via damped Newton-Raphson.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess vector (shape | required |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector. |
Source code in circulax/solvers/linear.py
solve_dc_auto ¤
solve_dc_auto(
component_groups: dict[str, Any],
y_guess: Array,
g_start: float = 0.01,
n_gmin: int = 10,
n_source: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point with automatic homotopy fallback.
Attempts a direct Newton solve first. If it fails to converge, falls back to GMIN stepping followed by source stepping — all inside a single JIT-compiled kernel via jax.lax.cond.
Strategy: 1. _run_newton — plain damped Newton from y_guess. 2. On failure: solve_dc_gmin (GMIN stepping) starting from y_guess, then solve_dc_source (source stepping) from the GMIN result.
Because jax.lax.cond evaluates both branches at trace time but only one at runtime, this compiles to a single kernel with no Python- level branching.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
g_start | float | Starting leakage for GMIN stepping (rescue branch). | 0.01 |
n_gmin | int | Number of GMIN steps in the rescue branch. | 10 |
n_source | int | Number of source steps in the rescue branch. | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector. |
Source code in circulax/solvers/linear.py
solve_dc_checked ¤
solve_dc_checked(
component_groups: dict[str, Any],
y_guess: Array,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> tuple[Array, Array]
DC Operating Point with convergence status.
Identical to :meth:solve_dc but additionally returns a boolean JAX scalar indicating whether the Newton-Raphson fixed-point iteration reported success. Because the flag is a JAX array (not a Python bool) it can be consumed inside compiled programs:
.. code-block:: python
y, converged = solver.solve_dc_checked(groups, y0)
# Outside JIT — inspect in Python:
if not converged:
y = solver.solve_dc_gmin(groups, y0)
# Inside JIT — branch without Python-level control flow:
y = jax.lax.cond(converged, lambda: y, lambda: solver.solve_dc_gmin(groups, y0))
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess vector (shape | required |
rtol | float | Relative tolerance for | 1e-06 |
atol | float | Absolute tolerance for | 1e-06 |
max_steps | int | Maximum Newton iterations. | 100 |
Returns:
| Type | Description |
|---|---|
tuple[Array, Array] |
|
Source code in circulax/solvers/linear.py
solve_dc_gmin ¤
solve_dc_gmin(
component_groups: dict[str, Any],
y_guess: Array,
g_start: float = 0.01,
n_steps: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point via GMIN stepping (homotopy rescue).
Steps the diagonal regularisation conductance g_leak logarithmically from g_start down to self.g_leak, using each converged solution as the warm start for the next step. The large initial g_leak linearises highly nonlinear components (diodes above threshold, lasers) that would otherwise cause Newton to diverge from a flat 0V start.
Implemented with jax.lax.scan — fully JIT/grad/vmap-compatible.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
g_start | float | Starting leakage conductance (large value, e.g. | 0.01 |
n_steps | int | Number of log-uniform steps from | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector after the full GMIN schedule. |
Source code in circulax/solvers/linear.py
solve_dc_source ¤
solve_dc_source(
component_groups: dict[str, Any],
y_guess: Array,
n_steps: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point via source stepping (homotopy rescue).
Ramps all source amplitudes (components tagged with amplitude_param) from 10 % to 100 % of their netlist values, using each converged solution as the warm start for the next step. This guides Newton through the nonlinear region without the large initial step from 0V to full excitation.
Implemented with jax.lax.scan — fully JIT/grad/vmap-compatible.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
n_steps | int | Number of uniformly-spaced steps from 0.1 to 1.0. | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector at full source amplitude. |
Source code in circulax/solvers/linear.py
FactorizedTransientSolver ¤
Bases: VectorizedTransientSolver
Transient solver using a Modified Newton (frozen-Jacobian) scheme.
At each timestep the system Jacobian is assembled and factored once at a predicted state, then reused across all Newton iterations. Compared to a full Newton-Raphson solver this trades quadratic convergence for a much cheaper per-iteration cost — one triangular solve instead of a full factorisation — making it efficient for circuits where the Jacobian varies slowly between steps.
Convergence is linear rather than quadratic, so newton_max_steps is set higher than a standard Newton solver would require. Adaptive damping min(1, 0.5 / max|δy|) is applied at each iteration to stabilise convergence in stiff or strongly nonlinear regions.
Both real and complex assembly paths are supported; the complex path concatenates real and imaginary parts into a single real-valued vector, allowing purely real linear algebra kernels to be reused for frequency-domain-style analyses.
Requires a :class:~circulax.solvers.linear.KLUSplitFactorSolver as the linear_solver — use analyze_circuit(..., backend="klu_split_factor").
KLUSolver ¤
Bases: CircuitLinearSolver
Solves the system using the KLU sparse solver (via klujax).
Best For
- Large circuits (N > 5000) running on CPU.
- DC Operating Points of massive meshes.
- Cases where DenseSolver runs out of memory (OOM).
Note
Does NOT support vmap (batching) automatically.
Methods:
| Name | Description |
|---|---|
assume_full_rank | Indicate if the solver assumes the operator is full rank. |
compute | Satisfies the lineax API; call |
from_component_groups | Factory method to pre-hash indices for sparse coalescence. |
init | Initialize the solver state (no-op for stateless solvers). |
solve_dc | DC operating point via damped Newton-Raphson. |
solve_dc_auto | DC Operating Point with automatic homotopy fallback. |
solve_dc_checked | DC Operating Point with convergence status. |
solve_dc_gmin | DC Operating Point via GMIN stepping (homotopy rescue). |
solve_dc_source | DC Operating Point via source stepping (homotopy rescue). |
compute ¤
Satisfies the lineax API; call _solve_impl directly for internal use.
from_component_groups classmethod ¤
from_component_groups(
component_groups: dict[str, Any],
num_vars: int,
*,
is_complex: bool = False,
g_leak: float = 1e-09
) -> KLUSolver
Factory method to pre-hash indices for sparse coalescence.
Source code in circulax/solvers/linear.py
init ¤
solve_dc ¤
solve_dc(
component_groups: dict[str, Any],
y_guess: Array,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC operating point via damped Newton-Raphson.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess vector (shape | required |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector. |
Source code in circulax/solvers/linear.py
solve_dc_auto ¤
solve_dc_auto(
component_groups: dict[str, Any],
y_guess: Array,
g_start: float = 0.01,
n_gmin: int = 10,
n_source: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point with automatic homotopy fallback.
Attempts a direct Newton solve first. If it fails to converge, falls back to GMIN stepping followed by source stepping — all inside a single JIT-compiled kernel via jax.lax.cond.
Strategy: 1. _run_newton — plain damped Newton from y_guess. 2. On failure: solve_dc_gmin (GMIN stepping) starting from y_guess, then solve_dc_source (source stepping) from the GMIN result.
Because jax.lax.cond evaluates both branches at trace time but only one at runtime, this compiles to a single kernel with no Python- level branching.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
g_start | float | Starting leakage for GMIN stepping (rescue branch). | 0.01 |
n_gmin | int | Number of GMIN steps in the rescue branch. | 10 |
n_source | int | Number of source steps in the rescue branch. | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector. |
Source code in circulax/solvers/linear.py
solve_dc_checked ¤
solve_dc_checked(
component_groups: dict[str, Any],
y_guess: Array,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> tuple[Array, Array]
DC Operating Point with convergence status.
Identical to :meth:solve_dc but additionally returns a boolean JAX scalar indicating whether the Newton-Raphson fixed-point iteration reported success. Because the flag is a JAX array (not a Python bool) it can be consumed inside compiled programs:
.. code-block:: python
y, converged = solver.solve_dc_checked(groups, y0)
# Outside JIT — inspect in Python:
if not converged:
y = solver.solve_dc_gmin(groups, y0)
# Inside JIT — branch without Python-level control flow:
y = jax.lax.cond(converged, lambda: y, lambda: solver.solve_dc_gmin(groups, y0))
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess vector (shape | required |
rtol | float | Relative tolerance for | 1e-06 |
atol | float | Absolute tolerance for | 1e-06 |
max_steps | int | Maximum Newton iterations. | 100 |
Returns:
| Type | Description |
|---|---|
tuple[Array, Array] |
|
Source code in circulax/solvers/linear.py
solve_dc_gmin ¤
solve_dc_gmin(
component_groups: dict[str, Any],
y_guess: Array,
g_start: float = 0.01,
n_steps: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point via GMIN stepping (homotopy rescue).
Steps the diagonal regularisation conductance g_leak logarithmically from g_start down to self.g_leak, using each converged solution as the warm start for the next step. The large initial g_leak linearises highly nonlinear components (diodes above threshold, lasers) that would otherwise cause Newton to diverge from a flat 0V start.
Implemented with jax.lax.scan — fully JIT/grad/vmap-compatible.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
g_start | float | Starting leakage conductance (large value, e.g. | 0.01 |
n_steps | int | Number of log-uniform steps from | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector after the full GMIN schedule. |
Source code in circulax/solvers/linear.py
solve_dc_source ¤
solve_dc_source(
component_groups: dict[str, Any],
y_guess: Array,
n_steps: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point via source stepping (homotopy rescue).
Ramps all source amplitudes (components tagged with amplitude_param) from 10 % to 100 % of their netlist values, using each converged solution as the warm start for the next step. This guides Newton through the nonlinear region without the large initial step from 0V to full excitation.
Implemented with jax.lax.scan — fully JIT/grad/vmap-compatible.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
n_steps | int | Number of uniformly-spaced steps from 0.1 to 1.0. | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector at full source amplitude. |
Source code in circulax/solvers/linear.py
KLUSplitLinear ¤
Bases: KLUSplitSolver
KLU split solver paired with Modified Newton (frozen-Jacobian) for linear convergence.
Extends :class:KLUSplitSolver with an explicit numeric factorization step so the Jacobian can be factored once per time step and reused across all Newton iterations within that step (Modified Newton / frozen-Jacobian scheme). Use together with :class:~circulax.solvers.transient.FactorizedTransientSolver.
Best For
- Large circuits (N > 5000) running on CPU where the Jacobian changes slowly.
- Transient simulations with many Newton iterations per step.
Methods:
| Name | Description |
|---|---|
assume_full_rank | Indicate if the solver assumes the operator is full rank. |
compute | Satisfies the lineax API; call |
factor_jacobian | Factor the Jacobian and return a numeric handle for repeated solves. |
from_component_groups | Factory — delegates to :meth: |
init | Initialize the solver state (no-op for stateless solvers). |
solve_dc | DC operating point via damped Newton-Raphson. |
solve_dc_auto | DC Operating Point with automatic homotopy fallback. |
solve_dc_checked | DC Operating Point with convergence status. |
solve_dc_gmin | DC Operating Point via GMIN stepping (homotopy rescue). |
solve_dc_source | DC Operating Point via source stepping (homotopy rescue). |
solve_with_frozen_jacobian | Solve using a pre-computed numeric factorization. |
compute ¤
Satisfies the lineax API; call _solve_impl directly for internal use.
factor_jacobian ¤
Factor the Jacobian and return a numeric handle for repeated solves.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
all_vals | Array | Flattened non-zero Jacobian values (COO format). | required |
Returns:
| Type | Description |
|---|---|
Array | Opaque numeric handle ( |
Array | meth: |
Array |
|
Source code in circulax/solvers/linear.py
from_component_groups classmethod ¤
from_component_groups(
component_groups: dict[str, Any],
num_vars: int,
*,
is_complex: bool = False,
g_leak: float = 1e-09
) -> KLUSplitLinear
Factory — delegates to :meth:KLUSplitSolver.from_component_groups.
Source code in circulax/solvers/linear.py
init ¤
solve_dc ¤
solve_dc(
component_groups: dict[str, Any],
y_guess: Array,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC operating point via damped Newton-Raphson.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess vector (shape | required |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector. |
Source code in circulax/solvers/linear.py
solve_dc_auto ¤
solve_dc_auto(
component_groups: dict[str, Any],
y_guess: Array,
g_start: float = 0.01,
n_gmin: int = 10,
n_source: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point with automatic homotopy fallback.
Attempts a direct Newton solve first. If it fails to converge, falls back to GMIN stepping followed by source stepping — all inside a single JIT-compiled kernel via jax.lax.cond.
Strategy: 1. _run_newton — plain damped Newton from y_guess. 2. On failure: solve_dc_gmin (GMIN stepping) starting from y_guess, then solve_dc_source (source stepping) from the GMIN result.
Because jax.lax.cond evaluates both branches at trace time but only one at runtime, this compiles to a single kernel with no Python- level branching.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
g_start | float | Starting leakage for GMIN stepping (rescue branch). | 0.01 |
n_gmin | int | Number of GMIN steps in the rescue branch. | 10 |
n_source | int | Number of source steps in the rescue branch. | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector. |
Source code in circulax/solvers/linear.py
solve_dc_checked ¤
solve_dc_checked(
component_groups: dict[str, Any],
y_guess: Array,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> tuple[Array, Array]
DC Operating Point with convergence status.
Identical to :meth:solve_dc but additionally returns a boolean JAX scalar indicating whether the Newton-Raphson fixed-point iteration reported success. Because the flag is a JAX array (not a Python bool) it can be consumed inside compiled programs:
.. code-block:: python
y, converged = solver.solve_dc_checked(groups, y0)
# Outside JIT — inspect in Python:
if not converged:
y = solver.solve_dc_gmin(groups, y0)
# Inside JIT — branch without Python-level control flow:
y = jax.lax.cond(converged, lambda: y, lambda: solver.solve_dc_gmin(groups, y0))
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess vector (shape | required |
rtol | float | Relative tolerance for | 1e-06 |
atol | float | Absolute tolerance for | 1e-06 |
max_steps | int | Maximum Newton iterations. | 100 |
Returns:
| Type | Description |
|---|---|
tuple[Array, Array] |
|
Source code in circulax/solvers/linear.py
solve_dc_gmin ¤
solve_dc_gmin(
component_groups: dict[str, Any],
y_guess: Array,
g_start: float = 0.01,
n_steps: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point via GMIN stepping (homotopy rescue).
Steps the diagonal regularisation conductance g_leak logarithmically from g_start down to self.g_leak, using each converged solution as the warm start for the next step. The large initial g_leak linearises highly nonlinear components (diodes above threshold, lasers) that would otherwise cause Newton to diverge from a flat 0V start.
Implemented with jax.lax.scan — fully JIT/grad/vmap-compatible.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
g_start | float | Starting leakage conductance (large value, e.g. | 0.01 |
n_steps | int | Number of log-uniform steps from | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector after the full GMIN schedule. |
Source code in circulax/solvers/linear.py
solve_dc_source ¤
solve_dc_source(
component_groups: dict[str, Any],
y_guess: Array,
n_steps: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point via source stepping (homotopy rescue).
Ramps all source amplitudes (components tagged with amplitude_param) from 10 % to 100 % of their netlist values, using each converged solution as the warm start for the next step. This guides Newton through the nonlinear region without the large initial step from 0V to full excitation.
Implemented with jax.lax.scan — fully JIT/grad/vmap-compatible.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
n_steps | int | Number of uniformly-spaced steps from 0.1 to 1.0. | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector at full source amplitude. |
Source code in circulax/solvers/linear.py
solve_with_frozen_jacobian ¤
Solve using a pre-computed numeric factorization.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
residual | Array | The right-hand side vector | required |
numeric | Array | Handle returned by :meth: | required |
Returns:
| Type | Description |
|---|---|
Solution | class: |
Source code in circulax/solvers/linear.py
KLUSplitQuadratic ¤
Bases: KLUSplitLinear
KLU split solver paired with full Newton for quadratic convergence via klu_refactor.
Extends :class:KLUSplitLinear with :meth:refactor_jacobian, which updates the numeric LU factorization in-place using klujax.refactor. The sparsity pattern is fixed for a given circuit topology, so KLU reuses the existing memory allocation and fill-reducing permutation — only the L/U values are recomputed. This gives full Newton (quadratic) convergence at a fraction of the cost of re-calling klu_factor at every iteration.
Use together with :class:~circulax.solvers.transient.RefactoringTransientSolver.
Best For
- Large circuits on CPU with nonlinear devices where quadratic convergence is desired.
- Transient simulations where the Jacobian changes significantly between Newton iterates.
Methods:
| Name | Description |
|---|---|
assume_full_rank | Indicate if the solver assumes the operator is full rank. |
compute | Satisfies the lineax API; call |
factor_jacobian | Factor the Jacobian and return a numeric handle for repeated solves. |
from_component_groups | Factory — delegates to :meth: |
init | Initialize the solver state (no-op for stateless solvers). |
refactor_jacobian | Update the numeric factorization in-place with new Jacobian values. |
solve_dc | DC operating point via damped Newton-Raphson. |
solve_dc_auto | DC Operating Point with automatic homotopy fallback. |
solve_dc_checked | DC Operating Point with convergence status. |
solve_dc_gmin | DC Operating Point via GMIN stepping (homotopy rescue). |
solve_dc_source | DC Operating Point via source stepping (homotopy rescue). |
solve_with_frozen_jacobian | Solve using a pre-computed numeric factorization. |
compute ¤
Satisfies the lineax API; call _solve_impl directly for internal use.
factor_jacobian ¤
Factor the Jacobian and return a numeric handle for repeated solves.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
all_vals | Array | Flattened non-zero Jacobian values (COO format). | required |
Returns:
| Type | Description |
|---|---|
Array | Opaque numeric handle ( |
Array | meth: |
Array |
|
Source code in circulax/solvers/linear.py
from_component_groups classmethod ¤
from_component_groups(
component_groups: dict[str, Any],
num_vars: int,
*,
is_complex: bool = False,
g_leak: float = 1e-09
) -> KLUSplitQuadratic
Factory — delegates to :meth:KLUSplitSolver.from_component_groups.
Source code in circulax/solvers/linear.py
init ¤
refactor_jacobian ¤
Update the numeric factorization in-place with new Jacobian values.
Reuses the existing memory allocation and fill-reducing permutation from the symbolic analysis; only the L/U values are recomputed. Faster than calling :meth:~KLUSplitLinear.factor_jacobian from scratch each Newton iteration.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
all_vals | Array | Flattened non-zero Jacobian values (COO format). | required |
numeric | Array | Existing handle returned by :meth: | required |
Returns:
| Type | Description |
|---|---|
Array | Refreshed numeric handle (same underlying C++ object, now connected in the |
Array | XLA computation graph so the refactor cannot be eliminated as dead code). |
Source code in circulax/solvers/linear.py
solve_dc ¤
solve_dc(
component_groups: dict[str, Any],
y_guess: Array,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC operating point via damped Newton-Raphson.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess vector (shape | required |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector. |
Source code in circulax/solvers/linear.py
solve_dc_auto ¤
solve_dc_auto(
component_groups: dict[str, Any],
y_guess: Array,
g_start: float = 0.01,
n_gmin: int = 10,
n_source: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point with automatic homotopy fallback.
Attempts a direct Newton solve first. If it fails to converge, falls back to GMIN stepping followed by source stepping — all inside a single JIT-compiled kernel via jax.lax.cond.
Strategy: 1. _run_newton — plain damped Newton from y_guess. 2. On failure: solve_dc_gmin (GMIN stepping) starting from y_guess, then solve_dc_source (source stepping) from the GMIN result.
Because jax.lax.cond evaluates both branches at trace time but only one at runtime, this compiles to a single kernel with no Python- level branching.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
g_start | float | Starting leakage for GMIN stepping (rescue branch). | 0.01 |
n_gmin | int | Number of GMIN steps in the rescue branch. | 10 |
n_source | int | Number of source steps in the rescue branch. | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector. |
Source code in circulax/solvers/linear.py
solve_dc_checked ¤
solve_dc_checked(
component_groups: dict[str, Any],
y_guess: Array,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> tuple[Array, Array]
DC Operating Point with convergence status.
Identical to :meth:solve_dc but additionally returns a boolean JAX scalar indicating whether the Newton-Raphson fixed-point iteration reported success. Because the flag is a JAX array (not a Python bool) it can be consumed inside compiled programs:
.. code-block:: python
y, converged = solver.solve_dc_checked(groups, y0)
# Outside JIT — inspect in Python:
if not converged:
y = solver.solve_dc_gmin(groups, y0)
# Inside JIT — branch without Python-level control flow:
y = jax.lax.cond(converged, lambda: y, lambda: solver.solve_dc_gmin(groups, y0))
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess vector (shape | required |
rtol | float | Relative tolerance for | 1e-06 |
atol | float | Absolute tolerance for | 1e-06 |
max_steps | int | Maximum Newton iterations. | 100 |
Returns:
| Type | Description |
|---|---|
tuple[Array, Array] |
|
Source code in circulax/solvers/linear.py
solve_dc_gmin ¤
solve_dc_gmin(
component_groups: dict[str, Any],
y_guess: Array,
g_start: float = 0.01,
n_steps: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point via GMIN stepping (homotopy rescue).
Steps the diagonal regularisation conductance g_leak logarithmically from g_start down to self.g_leak, using each converged solution as the warm start for the next step. The large initial g_leak linearises highly nonlinear components (diodes above threshold, lasers) that would otherwise cause Newton to diverge from a flat 0V start.
Implemented with jax.lax.scan — fully JIT/grad/vmap-compatible.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
g_start | float | Starting leakage conductance (large value, e.g. | 0.01 |
n_steps | int | Number of log-uniform steps from | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector after the full GMIN schedule. |
Source code in circulax/solvers/linear.py
solve_dc_source ¤
solve_dc_source(
component_groups: dict[str, Any],
y_guess: Array,
n_steps: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point via source stepping (homotopy rescue).
Ramps all source amplitudes (components tagged with amplitude_param) from 10 % to 100 % of their netlist values, using each converged solution as the warm start for the next step. This guides Newton through the nonlinear region without the large initial step from 0V to full excitation.
Implemented with jax.lax.scan — fully JIT/grad/vmap-compatible.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
n_steps | int | Number of uniformly-spaced steps from 0.1 to 1.0. | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector at full source amplitude. |
Source code in circulax/solvers/linear.py
solve_with_frozen_jacobian ¤
Solve using a pre-computed numeric factorization.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
residual | Array | The right-hand side vector | required |
numeric | Array | Handle returned by :meth: | required |
Returns:
| Type | Description |
|---|---|
Solution | class: |
Source code in circulax/solvers/linear.py
RefactoringTransientSolver ¤
Bases: FactorizedTransientSolver
Transient solver with full Newton (quadratic) convergence using klujax.refactor.
At each timestep the Jacobian is factored once at the predicted state to allocate the numeric handle. Each Newton iteration then calls klujax.refactor — which reuses the existing memory and fill-reducing permutation but recomputes L/U values for the current iterate J(y_k) — followed by a triangular solve. This gives full quadratic Newton convergence at a fraction of the cost of re-factoring from scratch each iteration.
Convergence is quadratic so newton_max_steps is set to 20, matching :class:VectorizedTransientSolver. Adaptive damping min(1, 0.5 / max|δy|) is applied at each iteration to stabilise convergence in stiff or strongly nonlinear regions.
Requires :class:~circulax.solvers.linear.KLUSplitQuadratic as the linear_solver — use analyze_circuit(..., backend="klu_split").
SDIRK3FactorizedTransientSolver ¤
Bases: FactorizedTransientSolver
3rd-order A-stable SDIRK3 solver with frozen-Jacobian Newton across all stages.
Factors J_eff = dF/dy + (1/(γh))·dQ/dy once at the predictor state, then reuses it for all Newton iterations in all three SDIRK stages. This is the recommended backend for large sparse circuits — the single factorisation is shared across all stages because SDIRK's constant diagonal γ gives the same effective Jacobian at every stage.
SDIRK3RefactoringTransientSolver ¤
Bases: RefactoringTransientSolver
3rd-order A-stable SDIRK3 solver with KLU refactor at each Newton iteration.
Provides full quadratic Newton convergence via klujax.refactor within each stage, combined with SDIRK3 time discretisation for 3rd-order accuracy.
SDIRK3VectorizedTransientSolver ¤
Bases: VectorizedTransientSolver
3rd-order A-stable SDIRK3 solver using full Newton-Raphson at each stage.
Uses Alexander's L-stable 3-stage SDIRK tableau with the companion method. Each timestep performs 3 sequential Newton solves (one per stage) with the Jacobian reassembled at every iteration. The same solver_state 2-tuple (y_prev, dt_prev) as Backward Euler is used — SDIRK3 is a one-step method.
SparseSolver ¤
Bases: CircuitLinearSolver
Solves the system using JAX's Iterative BiCGStab solver.
Best For
- Large Transient Simulations on GPU (uses previous step as warm start).
- Systems where N is too large for Dense, but we need VMAP support.
Attributes:
| Name | Type | Description |
|---|---|---|
diag_mask | Array | Mask to extract diagonal elements for preconditioning. |
Methods:
| Name | Description |
|---|---|
assume_full_rank | Indicate if the solver assumes the operator is full rank. |
compute | Satisfies the lineax API; call |
from_component_groups | Factory method to prepare indices and diagonal mask. |
init | Initialize the solver state (no-op for stateless solvers). |
solve_dc | DC operating point via damped Newton-Raphson. |
solve_dc_auto | DC Operating Point with automatic homotopy fallback. |
solve_dc_checked | DC Operating Point with convergence status. |
solve_dc_gmin | DC Operating Point via GMIN stepping (homotopy rescue). |
solve_dc_source | DC Operating Point via source stepping (homotopy rescue). |
compute ¤
Satisfies the lineax API; call _solve_impl directly for internal use.
from_component_groups classmethod ¤
from_component_groups(
component_groups: dict[str, Any],
num_vars: int,
*,
is_complex: bool = False,
g_leak: float = 1e-09
) -> SparseSolver
Factory method to prepare indices and diagonal mask.
Source code in circulax/solvers/linear.py
init ¤
solve_dc ¤
solve_dc(
component_groups: dict[str, Any],
y_guess: Array,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC operating point via damped Newton-Raphson.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess vector (shape | required |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector. |
Source code in circulax/solvers/linear.py
solve_dc_auto ¤
solve_dc_auto(
component_groups: dict[str, Any],
y_guess: Array,
g_start: float = 0.01,
n_gmin: int = 10,
n_source: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point with automatic homotopy fallback.
Attempts a direct Newton solve first. If it fails to converge, falls back to GMIN stepping followed by source stepping — all inside a single JIT-compiled kernel via jax.lax.cond.
Strategy: 1. _run_newton — plain damped Newton from y_guess. 2. On failure: solve_dc_gmin (GMIN stepping) starting from y_guess, then solve_dc_source (source stepping) from the GMIN result.
Because jax.lax.cond evaluates both branches at trace time but only one at runtime, this compiles to a single kernel with no Python- level branching.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
g_start | float | Starting leakage for GMIN stepping (rescue branch). | 0.01 |
n_gmin | int | Number of GMIN steps in the rescue branch. | 10 |
n_source | int | Number of source steps in the rescue branch. | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector. |
Source code in circulax/solvers/linear.py
solve_dc_checked ¤
solve_dc_checked(
component_groups: dict[str, Any],
y_guess: Array,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> tuple[Array, Array]
DC Operating Point with convergence status.
Identical to :meth:solve_dc but additionally returns a boolean JAX scalar indicating whether the Newton-Raphson fixed-point iteration reported success. Because the flag is a JAX array (not a Python bool) it can be consumed inside compiled programs:
.. code-block:: python
y, converged = solver.solve_dc_checked(groups, y0)
# Outside JIT — inspect in Python:
if not converged:
y = solver.solve_dc_gmin(groups, y0)
# Inside JIT — branch without Python-level control flow:
y = jax.lax.cond(converged, lambda: y, lambda: solver.solve_dc_gmin(groups, y0))
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess vector (shape | required |
rtol | float | Relative tolerance for | 1e-06 |
atol | float | Absolute tolerance for | 1e-06 |
max_steps | int | Maximum Newton iterations. | 100 |
Returns:
| Type | Description |
|---|---|
tuple[Array, Array] |
|
Source code in circulax/solvers/linear.py
solve_dc_gmin ¤
solve_dc_gmin(
component_groups: dict[str, Any],
y_guess: Array,
g_start: float = 0.01,
n_steps: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point via GMIN stepping (homotopy rescue).
Steps the diagonal regularisation conductance g_leak logarithmically from g_start down to self.g_leak, using each converged solution as the warm start for the next step. The large initial g_leak linearises highly nonlinear components (diodes above threshold, lasers) that would otherwise cause Newton to diverge from a flat 0V start.
Implemented with jax.lax.scan — fully JIT/grad/vmap-compatible.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
g_start | float | Starting leakage conductance (large value, e.g. | 0.01 |
n_steps | int | Number of log-uniform steps from | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector after the full GMIN schedule. |
Source code in circulax/solvers/linear.py
solve_dc_source ¤
solve_dc_source(
component_groups: dict[str, Any],
y_guess: Array,
n_steps: int = 10,
rtol: float = 1e-06,
atol: float = 1e-06,
max_steps: int = 100,
) -> Array
DC Operating Point via source stepping (homotopy rescue).
Ramps all source amplitudes (components tagged with amplitude_param) from 10 % to 100 % of their netlist values, using each converged solution as the warm start for the next step. This guides Newton through the nonlinear region without the large initial step from 0V to full excitation.
Implemented with jax.lax.scan — fully JIT/grad/vmap-compatible.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component_groups | dict[str, Any] | Compiled circuit components. | required |
y_guess | Array | Initial guess (typically | required |
n_steps | int | Number of uniformly-spaced steps from 0.1 to 1.0. | 10 |
rtol | float | Relative tolerance for each inner Newton solve. | 1e-06 |
atol | float | Absolute tolerance for each inner Newton solve. | 1e-06 |
max_steps | int | Max Newton iterations per step. | 100 |
Returns:
| Type | Description |
|---|---|
Array | Converged solution vector at full source amplitude. |
Source code in circulax/solvers/linear.py
VectorizedTransientSolver ¤
Bases: AbstractSolver
Transient solver that works strictly on FLAT (Real) vectors.
Delegates complexity handling to the 'linear_solver' strategy.
analyze_circuit ¤
analyze_circuit(
groups: list,
num_vars: int,
backend: str = "default",
*,
is_complex: bool = False,
g_leak: float = 1e-09
) -> CircuitLinearSolver
Initializes a linear solver strategy for circuit analysis.
Factory that selects and configures the numerical backend for solving the linear system derived from a circuit's topology.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
groups | list | A list of component groups that define the circuit's structure and properties. | required |
num_vars | int | The total number of variables in the linear system. | required |
backend | str | The name of the solver backend to use. Supported backends are 'klu', 'klu_split', 'dense', and 'sparse'. Defaults to 'default', which uses 'klu_split'. | 'default' |
is_complex | bool | A flag indicating whether the circuit analysis involves complex numbers. Defaults to False. | False |
Returns:
| Name | Type | Description |
|---|---|---|
CircuitLinearSolver | CircuitLinearSolver | An instance of a circuit linear solver strategy |
CircuitLinearSolver | configured for the specified backend and circuit parameters. |
Raises:
| Type | Description |
|---|---|
ValueError | If the specified backend is not supported. |
Source code in circulax/solvers/linear.py
assemble_system_complex ¤
assemble_system_complex(
y_guess: Array,
component_groups: dict,
t1: float,
dt: float,
source_scale: float = 1.0,
alpha: float = 1.0,
) -> tuple[Array, Array, Array]
Assemble the residual vectors and effective Jacobian values for an unrolled complex system.
The complex state vector is stored in unrolled (block) format: the first half of y_guess holds the real parts of all node voltages/states, the second half holds the imaginary parts. This avoids JAX's limited support for complex-valued sparse linear solvers by keeping all arithmetic real.
The Jacobian is split into four real blocks — RR, RI, IR, II — representing the partial derivatives of the real and imaginary residual components with respect to the real and imaginary state components respectively. The blocks are concatenated in RR→RI→IR→II order to match the sparsity index layout produced during compilation.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
y_guess | Array | Unrolled state vector of shape | required |
component_groups | dict | Compiled component groups returned by :func: | required |
t1 | float | Time at which the system is being evaluated. | required |
dt | float | Timestep duration, used to scale the reactive Jacobian blocks. | required |
source_scale | float | Multiplicative scale applied to source amplitudes (components whose | 1.0 |
alpha | float | Jacobian scaling factor for the reactive blocks. Use | 1.0 |
Returns:
| Type | Description |
|---|---|
Array | A three-tuple |
Array |
|
Array |
|
tuple[Array, Array, Array] |
|
Source code in circulax/solvers/assembly.py
258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 | |
assemble_system_real ¤
assemble_system_real(
y_guess: Array,
component_groups: dict,
t1: float,
dt: float,
source_scale: float = 1.0,
alpha: float = 1.0,
) -> tuple[Array, Array, Array]
Assemble the residual vectors and effective Jacobian values for a real system.
For each component group, evaluates the physics at t1 and computes the forward-mode Jacobian via jax.jacfwd. The effective Jacobian combines the resistive and reactive contributions as J_eff = df/dy + (alpha/dt) * dq/dy, where alpha=1 recovers Backward Euler and alpha=3/2 (uniform step) gives BDF2.
Components are processed in sorted key order to ensure a deterministic non-zero layout in the sparse Jacobian, which is required for the factorisation step.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
y_guess | Array | Current state vector of shape | required |
component_groups | dict | Compiled component groups returned by :func: | required |
t1 | float | Time at which the system is being evaluated. | required |
dt | float | Timestep duration, used to scale the reactive Jacobian block. | required |
source_scale | float | Multiplicative scale applied to source amplitudes (components whose | 1.0 |
alpha | float | Jacobian scaling factor for the reactive block. Use | 1.0 |
Returns:
| Type | Description |
|---|---|
Array | A three-tuple |
Array |
|
Array |
|
tuple[Array, Array, Array] |
|
Source code in circulax/solvers/assembly.py
circuit_diffeqsolve ¤
circuit_diffeqsolve(
terms: PyTree[AbstractTerm],
solver,
t0: RealScalarLike,
t1: RealScalarLike,
dt0: RealScalarLike,
y0: PyTree[ArrayLike],
args: PyTree[Any] = None,
*,
saveat: SaveAt = SaveAt(t1=True),
stepsize_controller: AbstractStepSizeController = ConstantStepSize(),
max_steps: int | None = 4096,
throw: bool = True,
checkpoints: int | None = None
) -> Solution
Stripped-down :func:diffrax.diffeqsolve for circuit simulation.
Identical calling convention to diffrax.diffeqsolve except:
adjoint,event,progress_meter,solver_state,controller_state, andmade_jumparguments are absent.t0 < t1is assumed; no backward-time integration.- No SDE/CDE terms.
saveat.denseandsaveat.stepsare not supported.
checkpoints controls the number of binomial checkpoints used by RecursiveCheckpointAdjoint (None = auto from max_steps).
Source code in circulax/solvers/circuit_diffeq.py
310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 | |
setup_ac_sweep ¤
setup_ac_sweep(
groups: dict[str, Any], num_vars: int, port_nodes: list[int], *, z0: float = 50.0
) -> Callable[[Array, Array], Array]
Configure and return a callable for AC small-signal S-parameter sweep.
Linearises the circuit DAE at the DC operating point and sweeps over an array of frequencies, returning the complex S-parameter matrix at each frequency. The returned callable is compatible with :func:jax.jit and :func:jax.vmap.
The analysis solves Y(jω) · V = RHS at each frequency, where::
Y(jω) = G + jωC + Y_fdomain(f) + port_terminations + ground_penalty
G = ∂F/∂y and C = ∂Q/∂y are extracted once at the DC operating point. Y_fdomain(f) is the admittance contribution from frequency-domain components, re-evaluated at each frequency.
S-parameter convention — matched-load verification:
- Matched load (Z_circuit = Z0) → S11 = 0
- Open circuit (Z_circuit → ∞) → S11 = +1
- Short circuit (Z_circuit = 0) → S11 = −1
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
groups | dict[str, Any] | Compiled component groups from :func: | required |
num_vars | int | Total number of scalar unknowns (second return value of :func: | required |
port_nodes | list[int] | Global node indices for each circuit port, in the desired port ordering. Obtain from the port-to-node map returned by :func: | required |
z0 | float | Reference impedance in ohms, applied uniformly to all ports. Defaults to 50.0. | 50.0 |
Returns:
| Type | Description |
|---|---|
Callable[[Array, Array], Array] | A callable |
Callable[[Array, Array], Array] |
|
Callable[[Array, Array], Array] |
|
Callable[[Array, Array], Array] |
|
Callable[[Array, Array], Array] | Compatible with :func: |
Source code in circulax/solvers/ac_sweep.py
53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 | |
setup_harmonic_balance ¤
setup_harmonic_balance(
groups: dict[str, Any],
num_vars: int,
freq: float,
num_harmonics: int = 5,
*,
is_complex: bool = False,
g_leak: float = 1e-09,
osc_node: int | None = None,
amplitude_tries: Array | None = None
) -> Callable[[Array], tuple[Array, Array]]
Configure and return a Harmonic Balance callable.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
groups | dict[str, Any] | Compiled component groups from :func: | required |
num_vars | int | Total number of state variables. | required |
freq | float | Fundamental frequency in Hz. | required |
num_harmonics | int | Number of harmonics; uses K = 2*N+1 time points. | 5 |
is_complex | bool |
| False |
g_leak | float | Diagonal regularisation to prevent singular Jacobians at floating nodes. | 1e-09 |
osc_node | int | None | State-vector index of the oscillator node (from | None |
amplitude_tries | Array | None | Amplitudes (V) to try when | None |
Returns:
| Type | Description |
|---|---|
Callable[[Array], tuple[Array, Array]] |
|
Callable[[Array], tuple[Array, Array]] | compatible with |
Source code in circulax/solvers/harmonic_balance.py
132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 | |
setup_transient ¤
setup_transient(
groups: list,
linear_strategy: CircuitLinearSolver,
transient_solver: AbstractSolver = None,
) -> Callable[..., Solution]
Configures and returns a function for executing transient analysis.
This function acts as a factory, preparing a transient solver that is pre-configured with the circuit's linear strategy. It returns a callable that executes the time-domain simulation using diffrax.diffeqsolve.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
groups | list | A list of component groups that define the circuit. | required |
linear_strategy | CircuitLinearSolver | The configured linear solver strategy, typically obtained from | required |
transient_solver | optional | The transient solver class to use. If None, | None |
Returns:
| Type | Description |
|---|---|
Callable[..., Solution] | Callable[..., Any]: A function that executes the transient analysis. |
Callable[..., Solution] | This returned function accepts the following arguments: t0 (float): The start time of the simulation. t1 (float): The end time of the simulation. dt0 (float): The initial time step for the solver. y0 (ArrayLike): The initial state vector of the system. saveat (diffrax.SaveAt, optional): Specifies time points at which to save the solution. Defaults to None. max_steps (int, optional): The maximum number of steps the solver can take. Defaults to 100000. throw (bool, optional): If True, the solver will raise an error on failure. Defaults to False. term (diffrax.AbstractTerm, optional): The term defining the ODE. Defaults to a zero-value ODETerm. stepsize_controller (diffrax.AbstractStepSizeController, optional): The step size controller. Defaults to |
Source code in circulax/solvers/transient.py
856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 | |