Utilities

This section groups smaller helper functionality used across the package. These utilities are mostly relevant for developers, advanced scripting workflows, and debugging.

Type checking

OSS-DBSv2 validates settings before running larger simulations. This reduces the chance that expensive solves fail only after geometry construction or meshing.

Other utility modules support tasks such as:

  • field post-processing

  • material lookup and settings handling

  • NIfTI image interaction

  • VTK export

For new users, these helpers usually stay in the background. They become more important when extending the package or building custom automation around it.

API reference

Settings and defaults

class ossdbs.utils.settings.Settings(partial_settings: dict)[source]

Bases: object

Default settings of OSS-DBS.

CONTACT_SETTING: ClassVar[dict] = {'Active': False, 'Contact_ID': 0, 'Current[A]': 0.0, 'Floating': False, 'MaxMeshSize': 1000000.0, 'MaxMeshSizeEdge': 1000000.0, 'SurfaceImpedance': {'Model': None, 'Parameters': {}}, 'Voltage[V]': 0.0}
CUSTOM_SETTING: ClassVar[dict] = {'BrainRegion': {'Center': {'x[mm]': 0, 'y[mm]': 0, 'z[mm]': 0}, 'Dimension': {'x[mm]': 40, 'y[mm]': 40, 'z[mm]': 40}, 'Shape': 'Ellipsoid'}}
ELECTRODE_SETTING: ClassVar[dict] = {'Contacts': [], 'Direction': {'x[mm]': 0.0, 'y[mm]': 0.0, 'z[mm]': 1.0}, 'EncapsulationLayer': {}, 'Name': 'BostonScientificVerciseDirected', 'Rotation[Degrees]': 0.0, 'TipPosition': {'x[mm]': 0.0, 'y[mm]': 0.0, 'z[mm]': 0.0}}
ENCAPSULATION_SETTING: ClassVar[dict] = {'DielectricModel': 'ColeCole4', 'DielectricParameters': None, 'Material': 'Gray matter', 'MaxMeshSize': 1000000.0, 'Thickness[mm]': 0.0}
SETTING: ClassVar[dict] = {'ActivationThresholdVTA[V-per-m]': None, 'CalcAxonActivation': False, 'ComputeCurrents': False, 'ComputeImpedance': False, 'DielectricAccuracy': 0.01, 'DielectricModel': {'CustomParameters': None, 'Type': 'ColeCole4'}, 'EQSMode': False, 'Electrodes': [], 'ExportElectrode': False, 'ExportFrequency': None, 'ExportVTK': False, 'ExportVTKSubdivision': 0, 'FEMOrder': 2, 'FailFlag': 'oss', 'ImpedanceAnalysis': {'Enabled': False, 'Frequencies': None, 'IncludeFloating': True}, 'MaterialDistribution': {'DTIPath': '', 'DiffusionTensorActive': False, 'MRIMapping': {'Blood': 4, 'CSF': 3, 'Gray matter': 1, 'Unknown': 0, 'White matter': 2}, 'MRIPath': '', 'WMMasking': True}, 'Mesh': {'AdaptiveMeshRefinement': {'Active': False, 'ErrorTolerance': 0.1, 'MaxIterations': 10}, 'HPRefinement': {'Active': False, 'Factor': 0.125, 'Levels': 2}, 'LoadMesh': False, 'LoadPath': '', 'MaterialRefinementSteps': 0, 'MeshSize': {'Edges': {}, 'Faces': {}, 'Volumes': {}}, 'MeshingHypothesis': {'MaxMeshSize': 1000000.0, 'MeshSizeFilename': '', 'Type': 'Default'}, 'SaveMesh': False, 'SavePath': 'mesh'}, 'ModelSide': 0, 'OutOfCore': False, 'OutputPath': 'Results', 'PathwayFile': None, 'PointModel': {'Lattice': {'Active': False, 'Center': {'x[mm]': 0, 'y[mm]': 0, 'z[mm]': 0}, 'CollapseVTA': False, 'Direction': {'x[mm]': 0, 'y[mm]': 0, 'z[mm]': 1}, 'ExportField': True, 'PointDistance[mm]': 0.1, 'Shape': {'x': 10, 'y': 10, 'z': 10}}, 'Pathway': {'Active': False, 'ExportField': False, 'FileName': ''}, 'VoxelLattice': {'Active': False, 'ExportField': True, 'Shape': {'x': 10, 'y': 10, 'z': 10}, 'TimeDomain': False}}, 'Solver': {'MaximumSteps': 10000, 'Precision': 1e-12, 'Preconditioner': 'bddc', 'PreconditionerKwargs': {}, 'Type': 'CG'}, 'StimSets': {'Active': False, 'StimSetsFile': None}, 'StimulationFolder': '/home/docs/checkouts/readthedocs.org/user_builds/oss-dbsv2/checkouts/latest/docs', 'StimulationSignal': {'CounterAmplitude': 1.0, 'CounterPulseWidth[us]': 0.0, 'CurrentControlled': False, 'CutoffFrequency': 1000000.0, 'Frequency[Hz]': 130.0, 'InterPulseWidth[us]': 0.0, 'ListOfFrequencies': [130.0], 'PulseTopWidth[us]': 0.0, 'PulseWidth[us]': 60.0, 'SpectrumMode': 'FullSpectrum', 'Type': 'Rectangle'}, 'Surfaces': [], 'TruncateAfterActivePartRatio': None}
complete_settings() dict[source]

Complete dictionary provided by user with default settings.

Type checking

class ossdbs.utils.type_check.TypeChecker[source]

Bases: object

Check types of input dictionary.

CONTACT_SETTING: ClassVar[dict] = {'Active': <class 'bool'>, 'Contact_ID': <class 'int'>, 'Current[A]': (<class 'int'>, <class 'float'>), 'Floating': <class 'bool'>, 'SurfaceImpedance': {'Model': (<class 'NoneType'>, <class 'str'>), 'Parameters': <class 'dict'>}, 'Voltage[V]': (<class 'int'>, <class 'float'>)}
ELECTRODE_SETTING: ClassVar[dict] = {'Contacts': <class 'list'>, 'Direction': {'x[mm]': (<class 'int'>, <class 'float'>), 'y[mm]': (<class 'int'>, <class 'float'>), 'z[mm]': (<class 'int'>, <class 'float'>)}, 'Name': <class 'str'>, 'Rotation[Degrees]': (<class 'int'>, <class 'float'>), 'TipPosition': {'x[mm]': (<class 'int'>, <class 'float'>), 'y[mm]': (<class 'int'>, <class 'float'>), 'z[mm]': (<class 'int'>, <class 'float'>)}}
TYPES: ClassVar[dict] = {'ActivationThresholdVTA[V-per-m]': (<class 'NoneType'>, <class 'float'>), 'AdaptiveMeshRefinement': <class 'bool'>, 'BrainRegion': {'Center': {'x[mm]': (<class 'int'>, <class 'float'>), 'y[mm]': (<class 'int'>, <class 'float'>), 'z[mm]': (<class 'int'>, <class 'float'>)}, 'Dimension': {'x[mm]': (<class 'int'>, <class 'float'>), 'y[mm]': (<class 'int'>, <class 'float'>), 'z[mm]': (<class 'int'>, <class 'float'>)}, 'Shape': <class 'str'>}, 'CalcAxonActivation': <class 'bool'>, 'Contacts': {'MaxMeshSize': <class 'float'>, 'MaxMeshSizeEdge': <class 'float'>}, 'DielectricAccuracy': <class 'float'>, 'DielectricModel': {'Type': <class 'str'>}, 'EQSMode': <class 'bool'>, 'Electrodes': <class 'list'>, 'EncapsulationLayer': {'Material': <class 'str'>, 'MaxMeshSize': (<class 'int'>, <class 'float'>), 'Thickness[mm]': (<class 'int'>, <class 'float'>)}, 'FEMOrder': <class 'int'>, 'FailFlag': <class 'str'>, 'MaterialDistribution': {'DTIPath': <class 'str'>, 'DiffusionTensorActive': <class 'bool'>, 'MRIPath': <class 'str'>, 'WMMasking': <class 'bool'>}, 'Mesh': {'LoadMesh': <class 'bool'>, 'LoadPath': <class 'str'>, 'MeshElementOrder': <class 'int'>, 'MeshingHypothesis': {'MaxMeshSize': (<class 'int'>, <class 'float'>), 'MeshSizeFilename': <class 'str'>, 'Type': <class 'str'>}, 'SaveMesh': <class 'bool'>, 'SavePath': <class 'str'>}, 'ModelSide': <class 'int'>, 'OutOfCore': <class 'bool'>, 'OutputPath': <class 'str'>, 'PathwayFile': (<class 'NoneType'>, <class 'str'>), 'PointModel': {'Lattice': {'Center': {'x[mm]': (<class 'int'>, <class 'float'>), 'y[mm]': (<class 'int'>, <class 'float'>), 'z[mm]': (<class 'int'>, <class 'float'>)}, 'CollapseVTA': <class 'bool'>, 'Direction': {'x[mm]': (<class 'int'>, <class 'float'>), 'y[mm]': (<class 'int'>, <class 'float'>), 'z[mm]': (<class 'int'>, <class 'float'>)}, 'ExportField': <class 'bool'>, 'PointDistance[mm]': (<class 'int'>, <class 'float'>), 'Shape': {'x': <class 'int'>, 'y': <class 'int'>, 'z': <class 'int'>}}, 'Pathway': {'Active': <class 'bool'>, 'ExportField': <class 'bool'>, 'FileName': <class 'str'>}, 'VoxelLattice': {'Active': <class 'bool'>, 'ExportField': <class 'bool'>, 'Shape': {'x': <class 'int'>, 'y': <class 'int'>, 'z': <class 'int'>}, 'TimeDomain': <class 'bool'>}}, 'Solver': {'MaximumSteps': <class 'int'>, 'Precision': (<class 'int'>, <class 'float'>), 'Preconditioner': <class 'str'>, 'PreconditionerKwargs': <class 'dict'>, 'Type': <class 'str'>}, 'StimSets': {'Active': <class 'bool'>, 'StimSetsFile': (<class 'NoneType'>, <class 'str'>)}, 'StimulationFolder': <class 'str'>, 'StimulationSignal': {'CounterAmplitude': <class 'float'>, 'CounterPulseWidth[us]': (<class 'int'>, <class 'float'>), 'CutoffFrequency': <class 'float'>, 'Frequency[Hz]': (<class 'int'>, <class 'float'>), 'InterPulseWidth[us]': (<class 'int'>, <class 'float'>), 'PulseTopWidth[us]': (<class 'int'>, <class 'float'>), 'PulseWidth[us]': (<class 'int'>, <class 'float'>), 'SpectrumMode': <class 'str'>, 'Type': <class 'str'>}, 'TruncateAfterActivePartRatio': (<class 'NoneType'>, <class 'float'>)}
classmethod check(settings: dict) None[source]

Check types in settings dictionary.

NIfTI image handling

class ossdbs.utils.nifti1image.DiffusionTensorImage(file_path: str)[source]

Bases: Nifti1Image

DTI image.

property components

Component and respective entry in flattened array.

get_component(component: str) ndarray[source]

Get component of DTI tensor.

Notes

Determine which index to extract based on data shape (x,y,z,6) Note the pairs (xy,yx), (yz,zy), and (xz,zx) map to the same indices

class ossdbs.utils.nifti1image.MagneticResonanceImage(file_path: str)[source]

Bases: Nifti1Image

MRI image.

class ossdbs.utils.nifti1image.Nifti1Image(file_path: str)[source]

Bases: object

Interface for NIfTI-1 images via nibabel.

file_path

File path of a NIfTI-1 image (.nii or .nii.gz).

Type:

str

property affine

Get affine transformation as 4x4 matrix.

property bounding_box: BoundingBox

Return the bounding box of the voxel data in voxel space.

Return type:

BoundingBox

property data: ndarray

Return the data of the nifti1 image.

Return type:

np.memmap

get_voxel_bounding_box(bounding_box: BoundingBox)[source]

Get bounding box around voxels in real space.

Notes

The dimensions of the box is different from the number of voxels scaled by the voxel size. The box can be skewed in reality (described by the affine transformation).

property header: nibabel.nifti1.Nifti1Header

Return the header of the nifti1 image.

Return type:

nibabel.nifti1.Nifti1Header

property trafo_cf

Get the transformation matrix as NGSolve function.

property trafo_matrix: ndarray

Get trafo matrix, i.e. a 3x3 matrix.

property translation: ndarray

Returns the lowest cartesian coordinates of the voxel data.

Return type:

tuple

property voxel_sizes: ndarray

Return the voxel size in xyz-direction.

property xyz_shape: tuple

Returns the number of voxels in x-, y- and z-direction.

Return type:

tuple

class ossdbs.utils.nifti1image.VTAImage(file_path: str)[source]

Bases: MagneticResonanceImage

VTA image containing only boolean data.

compute_dice_coefficent(reference_image) float[source]

Compute dice coefficient with image of same shape.

get_vta_volume() float[source]

Compute volume of VTA.

VTK export

class ossdbs.utils.vtk_export.FieldSolution(solution: ngsolve.CoefficientFunction, label: str, mesh: ngsolve.comp.Mesh, is_complex: bool)[source]

Bases: object

Data structure for NGSolve solutions.

is_complex: bool
label: str
mesh: ngsolve.comp.Mesh
save(filename: str, subdivision: int = 0) None[source]

Save solution to VTK file.

subdivision controls how many times each element is subdivided before sampling — higher values resolve the higher-order polynomial basis at the cost of larger files. 0 writes only the corner DOFs and produces visibly faceted plots in ParaView for HP-refined meshes.

solution: ngsolve.CoefficientFunction

Field computation

ossdbs.utils.field_computation.compute_field_magnitude(fields: ndarray)[source]

Compute magnitude of field vector.

ossdbs.utils.field_computation.compute_field_magnitude_from_components(Ex: ndarray, Ey: ndarray, Ez: ndarray)[source]

Compute magnitude of field vector.

VTA collapse

ossdbs.utils.collapse_vta.get_collapsed_VTA(field_on_points: ndarray, implantation_coordinate: ndarray, lead_direction: ndarray, lead_diam: float)[source]

Postprocess probing points to remove the electrode by inward sideways collapse.

Notes

The point coordinates are pulled into the electrode space to ‘counteract’ tissue misplacement due to the implantation. This breaks the grid regularity, so all nifti files will be created externally.

Parameters:
  • field_on_points (numpy.ndarray) – Scalar values on the lattice, Nx7

  • implantation_coordinate (numpy array) – Center of the first contact

  • lead_direction (numpy array) – Direction of the lead (from head to tail)

  • lead_diam (float) – Diameter of the electrode that will be compensated by the inward collapse

Returns:

field_on_points_collided – Array of scalar values on adjusted lattice, Nx7

Return type:

numpy.ndarray

Materials

Material-related utilities.

ossdbs.utils.materials.have_dielectric_properties_changed(dielectric_properties: dict, is_complex: bool, old_freq: float, new_freq: float, threshold: float) bool[source]

Check if dielectric properties have changed significantly between frequencies.

This function compares the dielectric properties at two frequencies and determines whether the change exceeds a threshold. This is useful for optimizing FEM solves by reusing the mesh when properties haven’t changed significantly.

Parameters:
  • dielectric_properties (dict) – Dictionary mapping material names to dielectric models. Each model should have a conductivity(omega) method.

  • is_complex (bool) – Whether to compare complex conductivity values.

  • old_freq (float) – Previous frequency in Hz.

  • new_freq (float) – New frequency in Hz.

  • threshold (float) – Relative error threshold (e.g., 0.01 for 1%).

Returns:

True if properties have changed more than threshold, False otherwise.

Return type:

bool

Example

>>> changed = have_dielectric_properties_changed(
...     dielectric_properties={'Gray matter': gray_model},
...     is_complex=True,
...     old_freq=100.0,
...     new_freq=1000.0,
...     threshold=0.01,
... )