Python scripting API
OSS-DBSv2 can be driven entirely from Python without the CLI. This is useful for parameter sweeps, custom post-processing, embedding OSS-DBSv2 in larger pipelines, or interactive exploration in notebooks.
Quick start
The simplest way to run a simulation from Python is to load a JSON settings dictionary and call the same pipeline that the CLI uses:
import json
import ossdbs
ossdbs.set_logger()
with open("input.json") as f:
settings = json.load(f)
from ossdbs.main import main_run
main_run(settings)
This is equivalent to ossdbs input.json on the command line.
Step-by-step API
For more control, individual pipeline steps can be called separately. The typical sequence mirrors the CLI pipeline:
import json
import ngsolve
import ossdbs
from ossdbs.api import (
build_brain_model,
generate_electrodes,
load_images,
prepare_dielectric_properties,
prepare_solver,
prepare_stimulation_signal,
prepare_volume_conductor_model,
run_volume_conductor_model,
set_contact_and_encapsulation_layer_properties,
validate_solver_settings,
)
from ossdbs.fem import ConductivityCF
from ossdbs.model_geometry import ModelGeometry
from ossdbs.utils.settings import Settings
ossdbs.set_logger()
# 1. Load and complete settings
with open("input.json") as f:
input_settings = json.load(f)
settings = Settings(input_settings).complete_settings()
# 2. Load imaging data
mri_image, dti_image = load_images(settings)
# 3. Build geometry
electrodes = generate_electrodes(settings)
brain_model = build_brain_model(settings, mri_image)
geometry = ModelGeometry(brain_model, electrodes)
set_contact_and_encapsulation_layer_properties(settings, geometry)
# 4. Prepare material model
validate_solver_settings(settings, geometry)
dielectric_properties = prepare_dielectric_properties(settings)
materials = settings["MaterialDistribution"]["MRIMapping"]
conductivity = ConductivityCF(
mri_image,
brain_model.brain_region,
dielectric_properties,
materials,
geometry.encapsulation_layers,
complex_data=settings["EQSMode"],
dti_image=dti_image,
)
# 5. Create solver and volume conductor
with ngsolve.TaskManager():
solver = prepare_solver(settings)
volume_conductor = prepare_volume_conductor_model(
settings, geometry, conductivity, solver
)
frequency_domain_signal = prepare_stimulation_signal(settings)
# 6. Refine mesh and solve
volume_conductor.prepare_mesh_refinements(
settings["Mesh"]["MaterialRefinementSteps"]
)
run_volume_conductor_model(
settings, volume_conductor, frequency_domain_signal
)
After the solve, the volume conductor object provides access to the solution:
# Access the potential at the last solved frequency
potential = volume_conductor.potential
# Evaluate at arbitrary points (Nx3 array in mm)
import numpy as np
points = np.array([[0.0, 0.0, 5.0], [0.0, 0.0, 6.0]])
values = volume_conductor.evaluate_potential_at_points(points)
# Impedance data (if ComputeImpedance was enabled)
impedances = volume_conductor.impedances
Key functions
The ossdbs.api module provides the building blocks:
load_images(settings)— load MRI and optional DTI datagenerate_electrodes(settings)— create electrode CAD modelsbuild_brain_model(settings, mri_image)— create the brain geometryset_contact_and_encapsulation_layer_properties(settings, geometry)— apply contact activation, voltages, currents, and surface impedanceprepare_dielectric_properties(settings)— set up the dielectric modelprepare_solver(settings)— instantiate solver and preconditionerprepare_volume_conductor_model(settings, geometry, conductivity, solver)— create the appropriate volume conductor variantprepare_stimulation_signal(settings)— build the frequency-domain signalrun_volume_conductor_model(settings, vc, signal)— solve at all frequencies, export resultsgenerate_point_models(settings)— create lattice or pathway evaluators
Settings handling
The Settings class in ossdbs.utils.settings merges user-provided
values with defaults. Call Settings(input_dict).complete_settings() to
get a fully populated dictionary before passing it to any API function.
This is the same step the CLI performs internally.
API reference
- ossdbs.api.build_brain_model(settings, mri_image: MagneticResonanceImage | None = None, rotate_initial_geo: bool = False) BrainGeometry[source]
Build geometry model of brain.
- ossdbs.api.create_bounding_box(box_parameters: dict) BoundingBox[source]
Create a bounding box around a domain in space.
Notes
box_parameters need to contain the center (three coordinates, all following the style x[mm] (and similar for y and z-component)) and the outer dimensions.
- ossdbs.api.generate_brain_model(settings, rotate_initial_geo: bool = False)[source]
Generate OCC brain model.
- ossdbs.api.generate_electrodes(settings: dict)[source]
Generate an OCC electrode model from the settings dict.
- ossdbs.api.generate_mesh(settings)[source]
Generate a mesh from settings.
Notes
Attention! This mesh is not yet curved!
- ossdbs.api.generate_model_geometry(settings)[source]
Generate a full geometry comprising brain and electrodes.
- ossdbs.api.generate_signal(settings) TimeDomainSignal[source]
Generate a time-domain signal (waveform).
- ossdbs.api.prepare_dielectric_properties(settings: dict) dict[source]
Return dictionary with dielectric properties for each tissue.
- ossdbs.api.prepare_stimulation_signal(settings) FrequencyDomainSignal[source]
Prepare the frequency-domain representation of stimulation signal.
- ossdbs.api.prepare_volume_conductor_model(settings, model_geometry, conductivity, solver) VolumeConductor[source]
Prepare the volume conductor model.
- ossdbs.api.run_stim_sets(settings, geometry, conductivity, solver, frequency_domain_signal)[source]
Run StimSets batch workflow: compute unit solutions per contact.
For each non-ground contact, sets up a unit-current solve (1 A on that contact, all others floating, ground at -1 A), loads the pre-saved mesh, applies HP refinement, and runs the full FEM analysis. Results are stored in per-contact output directories.
- Parameters:
settings (dict) – Complete simulation settings dictionary.
geometry (ModelGeometry) – Geometry with electrode and contact definitions.
conductivity (ConductivityCF) – Conductivity coefficient function.
solver (Solver) – Configured FEM solver.
frequency_domain_signal (FrequencyDomainSignal) – Signal defining the frequencies and amplitudes to solve.
- ossdbs.api.run_volume_conductor_model(settings, volume_conductor, frequency_domain_signal, truncation_time=None)[source]
Run the volume conductor model at all required frequencies.
Solves the FEM system at each frequency in the signal spectrum, optionally computes impedance and currents, exports results, and reconstructs the time-domain solution.
- Parameters:
settings (dict) – Complete simulation settings dictionary.
volume_conductor (VolumeConductor) – Prepared volume conductor model with mesh and FEM space.
frequency_domain_signal (FrequencyDomainSignal) – Signal defining the frequencies and amplitudes to solve.
truncation_time (float, optional) – If set, truncate the time-domain reconstruction to this duration.
- ossdbs.api.set_contact_and_encapsulation_layer_properties(settings, model_geometry)[source]
Update boundary and material values on contacts and encapsulation layers.
- ossdbs.api.validate_solver_settings(settings: dict, model_geometry: ModelGeometry) None[source]
Validate and adjust solver settings based on model configuration.
Notes
BDDC preconditioner does not work well with the FloatingImpedance formulation. This function enforces the use of ‘local’ preconditioner in such cases.