Volume conductor model
The volume conductor model is the numerical core of OSS-DBSv2. It combines the geometry, tissue properties, stimulation signal, and solver configuration to compute electric potentials and fields in the simulated domain.
Role in the workflow
Most standalone and Lead-DBS-driven runs eventually reach the same core steps:
build the geometry and mesh
assign conductivity information from MRI and optional DTI data
define boundary conditions at electrode contacts and outer surfaces
solve the resulting finite-element system
export fields, impedance estimates, and pointwise analyses
Floating and non-floating contacts
OSS-DBSv2 distinguishes between several conductor variants:
VolumeConductorNonFloatingfor standard active and passive contact setupsVolumeConductorFloatingwhen floating contacts are presentVolumeConductorFloatingImpedancewhen floating contacts also carry surface impedance information
These variants share the same general workflow but differ in how contact boundary conditions are enforced.
Stimulation modes
OSS-DBSv2 supports both voltage-controlled and current-controlled stimulation.
The CurrentControlled flag in StimulationSignal selects the mode.
Depending on the mode, contacts can be active (carrying a Dirichlet
boundary condition), floating (equipotential via Lagrange multiplier), or
floating with surface impedance (coupled via a Robin boundary condition).
Not every combination is valid. Below, the supported cases are summarised.
Voltage-controlled stimulation
In voltage-controlled mode (CurrentControlled: false), the Voltage[V]
value on each active contact is imposed directly as a Dirichlet boundary
condition.
Two or more active contacts: each contact is set to its prescribed voltage. At least one should be at 0 V to serve as ground. This is the standard bipolar or multipolar voltage-controlled configuration.
One active contact: a single Dirichlet condition is imposed. Floating contacts may carry the return current. Without a second active contact or floating contacts, the solution is trivially zero everywhere.
Floating contacts (plain): floating contacts are constrained to be equipotential surfaces via a Lagrange multiplier in the
VolumeConductorFloatingformulation. The potential on each floating contact is determined by current conservation.Floating contacts with surface impedance: the
VolumeConductorFloatingImpedanceformulation replaces the equipotential constraint with a Robin boundary condition that models the electrode–tissue interface impedance. The contact potential is a free variable coupled to the tissue potential through the admittance \(1/Z_s\).
Current-controlled stimulation
In current-controlled mode (CurrentControlled: true), each contact
specifies a Current[A] value. The sum of all currents (active and floating)
must be zero (Kirchhoff’s current law). Internally, OSS-DBSv2 uses an
admittance-matrix superposition to translate prescribed currents into the
potentials that the FEM solver needs.
Two active contacts (bipolar): one contact must be grounded (
Voltage[V]: 0). The other is assigned a pseudo-voltage for an intermediate Dirichlet solve; the final solution is scaled so that the prescribed current flows. Floating contacts may be present in either the plain or surface-impedance variant.One active contact (monopolar with floating): the single active contact must be grounded (
Voltage[V]: 0). All other contacts are floating with prescribed currents. In the plain floating case, each floating contact is equipotential. In the surface-impedance case, Robin BCs are used instead.No active contacts (all floating with surface impedance): this is only valid in the
VolumeConductorFloatingImpedanceformulation. Because there is no Dirichlet boundary to fix the voltage gauge, a Lagrange multiplier constrains the sum of floating potentials to zero (\(\sum_k u_k = 0\)). Currents are imposed on the floating contacts directly.More than two active contacts: currently not supported in current-controlled mode. Use at most one grounded active contact plus floating contacts to model multipolar configurations with more than two electrodes.
Constraints and validation
The code enforces several consistency checks before solving:
In current-controlled mode, the sum of all prescribed currents must be zero.
In multipolar current-controlled mode with plain floating contacts, exactly one active contact must be grounded.
The
FloatingImpedanceformulation requires every floating contact to carry a surface impedance model. Mixed setups (some floating contacts with impedance, some without) are rejected at construction time.Active contacts must not carry a surface impedance model in the
FloatingImpedanceformulation.
Impedance computation
When ComputeImpedance is enabled, OSS-DBSv2 computes the complex impedance
at each frequency. For current-controlled setups with multiple contacts, a full
admittance matrix is assembled by superposition and then inverted. For setups
with surface impedance, the dissipation across the electrode–tissue interface is
included in the power balance so that the impedance estimate accounts for the
interface losses.
Surface impedance
Real DBS electrodes do not make perfect electrical contact with surrounding tissue. A thin electrochemical interface forms at the metal–tissue boundary, adding a frequency-dependent impedance that affects current flow and measured impedance values. OSS-DBSv2 models this interface as a surface impedance applied as a Robin boundary condition on the contact surface.
Physical background
The interface is often described by an equivalent circuit (Randles-type model). Common elements include a resistor for charge-transfer resistance and a constant-phase element (CPE) for the double-layer capacitance. The surface impedance \(Z_s\) relates the potential jump across the interface to the normal current density:
where \(u\) is the tissue-side potential and \(u_{\text{contact}}\) is the contact potential. The surface impedance has units of \(\Omega \cdot \text{mm}^2\) (impedance times contact area) because the FEM geometry uses millimetres.
JSON configuration
Surface impedance is set per contact in the Electrodes section:
{
"Contact_ID": 1,
"Active": true,
"Voltage[V]": 1.0,
"Floating": false,
"SurfaceImpedance": {
"Model": "R",
"Parameters": {"R": 500.0}
}
}
The Model string and Parameters dictionary are passed to the
impedancefitter library, which
evaluates the equivalent circuit at each frequency. This means any model
supported by impedancefitter can be used.
Available models
Commonly used models include:
R— pure resistance. Parameters:{"R": <value in Ohm>}.RC— resistor in parallel with a capacitor. Parameters:{"R": <Ohm>, "C": <Farad>}.CPE_dl— constant-phase element for the double layer. Parameters:{"dl_k": <value>, "dl_alpha": <exponent 0–1>}. This is the model used by Lempka et al. (2009) for DBS electrode interfaces.
At zero frequency (DC), capacitive and CPE elements have infinite impedance. Only the resistive path contributes, so a pure CPE model yields no interface effect at DC.
Where surface impedance can be applied
On active (nonfloating) contacts: the
VolumeConductorNonFloatingformulation adds a Robin term for each active contact that carries aSurfaceImpedanceentry. The Dirichlet voltage is replaced by a Robin condition coupling the tissue potential to the prescribed contact voltage.On floating contacts: the
VolumeConductorFloatingImpedanceformulation couples each floating contact to the tissue via the Robin condition. The contact potential becomes a free variable solved for by the FEM system. Every floating contact must carry a surface impedance model in this formulation.
See Stimulation modes for the full list of valid contact configurations.
Effect on impedance estimates
When ComputeImpedance is enabled, the power balance includes both the
volume dissipation (tissue) and the interface dissipation (surface impedance
layer). Without accounting for the interface term, the computed impedance would
overestimate the true value because part of the dissipated power is missing.
Mesh and linear algebra
The volume conductor stack also includes:
conductivity coefficient functions derived from the material model
iterative and direct solvers
several preconditioner strategies for larger problems
Solver selection
OSS-DBSv2 provides three solver types:
CG (Conjugate Gradient) — the default. Efficient for symmetric positive definite systems, which covers most real-valued and EQS volume conductor problems. Use this as the starting point.
GMRES — required for non-symmetric or complex-valued systems. When
EQSModeis enabled, the system matrix becomes complex. GMRES also handles theFloatingImpedanceformulation where the Lagrange multiplier and Robin terms could break strict symmetry.Direct — factorises the system matrix with a direct solver like PARDISO. Robust and deterministic but memory-intensive. Practical for small to medium problems (up to ~200k DOFs) or as a reference for verifying iterative solver results.
Preconditioner selection
The preconditioner accelerates convergence of the iterative solvers (CG and GMRES). Available options:
bddc (default) — Balancing Domain Decomposition by Constraints. The most efficient option for standard (nonfloating) problems, especially at higher FEM orders. Not compatible with the
FloatingImpedanceformulation — the code automatically falls back tolocalin that case.local — NGSolve’s built-in Jacobi (diagonal) preconditioner. Works with all formulations. Sufficient for moderate problem sizes. If the solver produces unexpected results with
local, trycustomized_localinstead.customized_local — a custom Jacobi implementation that explicitly constructs the inverse diagonal. More robust than the native
localpreconditioner for problems with mixed spaces (e.g. FloatingImpedance with Lagrange multipliers). Use this iflocalgives warnings about non-finite values.h1amg — Algebraic multigrid. Can be faster than
localfor very large problems but is less tested in this codebase.multigrid — Geometric multigrid. Requires a mesh hierarchy; not commonly used in practice.
direct — Uses a direct solve as preconditioner. Expensive but can be useful for debugging convergence problems.
Recommended configurations
Scenario |
Solver |
Preconditioner |
|---|---|---|
Standard (real-valued) |
CG |
bddc |
EQS mode (complex) |
CG or GMRES |
bddc |
FloatingImpedance |
CG or GMRES |
local |
FloatingImpedance + EQS |
CG or GMRES |
customized_local |
Small problem / debugging |
Direct |
(not applicable) |
Convergence parameters
MaximumSteps(default: 10000) — upper limit on Krylov iterations. If the solver hits this limit, it raises aRuntimeError. Increase this value for very large problems or tight precision targets.Precision(default: 1e-12) — relative residual tolerance. For most applications, 1e-8 to 1e-10 is sufficient. Very tight values (1e-12 and below) increase iteration count without visible improvement in the solution.PrintRates— whentrue, the solver prints the residual norm at each iteration. Useful for monitoring convergence.
Troubleshooting
Solver does not converge: try increasing
MaximumStepsor relaxingPrecision. If that does not help, switch fromlocaltocustomized_localor from CG to GMRES.Warning about non-finite preconditioner values: switch from
localtocustomized_local. This is most common with theFloatingImpedanceformulation.BDDC fails with FloatingImpedance: this is expected — the code automatically switches to
local. No action needed.Very slow convergence: check that the mesh is not excessively refined in regions far from the electrode. Consider HP refinement instead of uniform mesh refinement for better cost-accuracy tradeoff.
Mesh refinement
Mesh quality strongly affects the accuracy of the FEM solution, especially near electrode contacts where the field gradient is steepest. OSS-DBSv2 provides several levels of mesh control.
Global mesh hypothesis
The MeshingHypothesis block sets the baseline element size:
"MeshingHypothesis": {
"Type": "Moderate"
}
Supported types are VeryCoarse, Coarse, Moderate, Fine,
VeryFine, and Default. A Custom type is also available, in which
case a CustomParameters dictionary is passed directly to Netgen.
Optional parameters can be added regardless of type:
MaxMeshSize— global upper bound on element size (in mm)CurvatureSafety— controls refinement near curved surfacesGrading— mesh grading factor (higher = faster size transition)MeshSizeFilename— path to a Netgen mesh-size file that specifies target element sizes at individual points. This is especially useful for refining the mesh along neuron trajectories in pathway activation studies. The file can be generated from aPathwayPointModelwithwrite_netgen_meshsize_file():from ossdbs.api import generate_point_models, load_images mri_image, _ = load_images(settings) point_models = generate_point_models(settings) pw = point_models[0] # PathwayPointModel pw.write_netgen_meshsize_file( meshsize=min(mri_image.voxel_sizes), filename="meshsizes.txt", )
Then reference the file in the JSON input:
"MeshingHypothesis": { "Type": "Fine", "MeshSizeFilename": "meshsizes.txt" }
Local mesh sizes
The MeshSize block allows targeted refinement on named geometry entities:
"MeshSize": {
"Edges": {"E1C1": 0.05},
"Faces": {"E1C1": 0.1, "E1C2": 0.1},
"Volumes": {"EncapsulationLayer_1": 0.2}
}
Keys are boundary or region names as they appear in the geometry (e.g.
E1C1 for electrode 1, contact 1). Contact-level MaxMeshSize and
MaxMeshSizeEdge in the Contacts JSON block provide a shorthand for
the same mechanism.
Mesh element order
MeshElementOrder (default: 2) sets the polynomial order of the curved mesh
geometry. Higher orders better approximate the true electrode surface but
increase mesh generation cost. This is independent of the FEM polynomial order
set by FEMOrder.
HP refinement
HP refinement combines geometric (h) and polynomial (p) refinement near boundaries, concentrating degrees of freedom where the solution varies most:
"HPRefinement": {
"Active": true,
"Levels": 2,
"Factor": 0.125
}
Levels— number of refinement layers added near boundariesFactor— geometric grading factor; smaller values concentrate elements more tightly at the boundary
HP refinement provides the best cost-accuracy tradeoff for DBS simulations, typically reaching near-best accuracy at 60–120k degrees of freedom. It is applied after any material-based bisection refinement so that the specialised HP element types are not disrupted by subsequent bisection steps.
Material-based mesh refinement
When tissue boundaries are poorly resolved by the initial mesh, bisection-based
refinement can be applied. This is controlled by the
MaterialRefinementSteps parameter (default: 0) and uniformly bisects
elements that straddle material interfaces.
Mesh reuse
For parameter studies where only the stimulation changes but the geometry stays fixed, mesh generation can be skipped:
"Mesh": {
"SaveMesh": true,
"LoadMesh": true,
"LoadPath": "path/to/saved_mesh.vol.gz"
}
When loading a mesh, HP refinement is re-applied automatically if configured.
Practical guidance
Start with
Moderatehypothesis andHPRefinementwith 2 levels.Check convergence by comparing results across refinement levels before drawing scientific conclusions.
Use local face/edge mesh sizes on contacts when the global hypothesis does not provide enough resolution near the electrode.
For large studies, save the base mesh and reload it to avoid repeated mesh generation.
API reference
Volume conductor model
- class ossdbs.fem.volume_conductor.volume_conductor_model.VolumeConductor(geometry: ModelGeometry, conductivity: ConductivityCF, solver: Solver, order: int, meshing_parameters: dict, output_path: str = 'Results')[source]
Bases:
ABCTemplate class of a volume conductor.
- Parameters:
geometry (ModelGeometry) – Model geometry, brain with implanted electrodes
conductivity (ConductivityCF) – Material information
solver (Solver) – Solver (linear algebra part)
order (int) – Order of solver and mesh (curved elements)
meshing_parameters (dict) – Dictionary with setting for meshing
output_path (str) – Path to store solution
- apply_h_refinements(material_mesh_refinement_steps: int = 0)[source]
Apply only h-refinements (material bisection).
Use this when the mesh will be saved for later reuse (e.g. StimSets). HP refinement is deferred so it can be applied on each loaded mesh instance.
- apply_hp_and_update_space()[source]
Apply deferred HP refinement and rebuild the FEM space.
Call after loading an h-refined mesh to complete the refinement pipeline.
- compute_impedance() complex[source]
Compute scalar impedance at most recent solution.
For two active contacts, the scalar impedance is computed by volume integration. This approach is superior to integration of the normal current density. It has been described in [Zimmermann2021a].
Multicontact configurations are not handled here. The full admittance-matrix analysis is a separate analysis tool (see
docs/impedance_analyzer_plan.md).References
[Zimmermann2021a]Zimmermann, J., et al. (2021). Frontiers in Bioengineering and Biotechnology, 9, 765516. https://doi.org/10.3389/fbioe.2021.765516
- Returns:
Scalar impedance between the two active contacts.
- Return type:
complex
- abstract compute_solution(frequency: float) None[source]
Compute solution at frequency.
- Parameters:
frequency (float) – Frequency at which solution is computed
- property conductivity: ngsolve.CoefficientFunction
Return conductivity of latest solution.
- property conductivity_cf: ConductivityCF
Returns the coefficient function of the conductivity.
- property current_controlled: bool
Return if stimulation is current-controlled.
- property current_density: ngsolve.GridFunction
Return current density in A/mm^2.
- property electric_field: ngsolve.GridFunction
Compute electric field from potential.
- estimate_currents() dict[source]
Estimate currents by integration of normal component.
Notes
Meant for debugging purposes. If singularities are present, this method will not be accurate.
- evaluate_field_at_points(lattice: ndarray) ndarray[source]
Return electric field components at specifed 3-D coordinates.
- Parameters:
lattice (np.ndarray) – Nx3 numpy.ndarray of lattice points
Notes
Requires that points outside of the computational domain have been filtered!
- evaluate_potential_at_points(lattice: ndarray) ndarray[source]
Return electric potential at specifed 3-D coordinates.
- Parameters:
lattice (np.ndarray) – Nx3 numpy.ndarray of lattice points
Notes
Requires that points outside of the computational domain have been filtered!
- export_material_distribution_to_vtk(subdivision: int = 0) None[source]
Write material distribution to VTK file.
- export_solution_to_vtk(freq_idx: int, multisine_mode: bool = False, subdivision: int = 0) None[source]
Export potential and field at frequency to VTK.
- Parameters:
freq_idx (int) – Index of frequency
multisine_mode (bool) – If rectangular pulse is used (multisine_mode = False)
subdivision (int) – Element subdivision count forwarded to ngsolve.VTKOutput.
- flux_space() ngsolve.comp.HDiv[source]
Return a flux space on the mesh.
- Return type:
ngsolve.HDiv
Notes
The HDiv space is returned with a minimum order of 1. It is needed for the a-posteriori error estimator needed for adaptive mesh refinement.
- property frequency: float
Most recent frequency, not equal to the frequency of the signal!.
- get_scale_factor(freq_idx: int) float[source]
Scale solution by signal amplitude at a frequency given by index.
Notes
In voltage-controlled mode, only the amplitude of the Fourier coefficient is used. In current-controlled mode without using floating conductors, the impedance is also considered.
- h1_space(boundaries: list[str], is_complex: bool) ngsolve.H1[source]
Return a h1 space on the mesh.
- Parameters:
boundaries (list of str) – list of boundary names.
is_complex (bool) – Whether to use complex arithmetic
- Return type:
ngsolve.H1
- property impedances: ndarray
Scalar impedance per frequency (1-D array, length
n_freq).Populated by
run_full_analysiswhencompute_impedanceis True. Multicontact admittance / impedance matrices are produced byossdbs.fem.analysis.ImpedanceAnalyzer, not here.
- property is_complex: bool
Return the state of the data type for spaces. True if complex, False otherwise.
- Return type:
bool
- property model_geometry: ModelGeometry
The underlying model geometry used for mesh generation.
- number_space() ngsolve.comp.NumberSpace[source]
Return a number space on the mesh.
- Returns:
ngsolve.NumberSpace – Space with only one single (global) DOF.
TODO check if needed
- property output_path: str
Returns the path to output.
- property potential: ngsolve.GridFunction
Return solution at most recent frequency.
- prepare_mesh_refinements(material_mesh_refinement_steps: int = 0)[source]
Apply material and HP mesh refinements.
- refine_mesh_by_material(material_mesh_refinement_steps: int) None[source]
Check materials and refine mesh if more than one material per element.
- run_full_analysis(frequency_domain_signal: FrequencyDomainSignal, compute_impedance: bool = False, export_vtk: bool = False, point_models: list[ossdbs.point_analysis.point_model.PointModel] | None = None, activation_threshold: float | None = None, dielectric_threshold: float = 0.01, out_of_core: bool = False, export_frequency: float | None = None, adaptive_mesh_refinement_settings: dict | None = None, truncation_time: float | None = None, estimate_currents: bool | None = False, vtk_subdivision: int = 0) dict[source]
Run volume conductor model at all frequencies.
- Parameters:
compute_impedance (bool) – If True, the impedance will be computed at each frequency.
export_vtk (bool) – VTK export for visualization in ParaView
point_models (list[PointModel]) – list of PointModel to extract solution for VTA / PAM
activation_threshold (float) – If VTA is estimated by threshold, provide it here. Its unit must be V/m!
dielectric_threshold (float) – Threshold for accuracy of dielectric properties
out_of_core (bool) – Indicate whether point model shall be done out-of-core
export_frequency (float) – Frequency at which the VTK file should be exported. Otherwise, median frequency is used.
frequency_domain_signal (FrequencyDomainSignal) – Frequency-domain representation of stimulation signal
adaptive_mesh_refinement_settings (dict) – Perform adaptive mesh refinement (only at first frequency)
truncation_time (float) – Time until which result will be written to hard drive
estimate_currents (bool) – Get current estimate per contact by integration of normal component
vtk_subdivision (int) – Element subdivision count forwarded to
ngsolve.VTKOutput.
Notes
The volume conductor model is run at all frequencies and the time-domain signal is computed (if relevant).
- setup_timings_dict(export_vtk: bool, point_models: list[ossdbs.point_analysis.point_model.PointModel]) dict[source]
Setup dictionary to save execution times estimate.
- property signal: FrequencyDomainSignal
Returns the frequency-domain representation of stimulation signal.
- threshold_frequency_domain_Efield(scale_factor: float, activation_threshold: float) float[source]
Determine volume of E-field above threshold at current frequency.
- vtk_export(freq_idx: int, multisine_mode: bool = False, subdivision: int = 0) None[source]
Export all relevant properties to VTK.
- Parameters:
freq_idx (int) – Index of frequency
multisine_mode (bool) – If rectangular pulse is used (multisine_mode = False)
subdivision (int) – Element subdivision count forwarded to ngsolve.VTKOutput.
- class ossdbs.fem.volume_conductor.nonfloating.VolumeConductorNonFloating(geometry: ModelGeometry, conductivity: ConductivityCF, solver: Solver, order: int, meshing_parameters: dict, output_path: str = 'Results')[source]
Bases:
VolumeConductorModel for representing a volume conductor which evaluates the potential.
- compute_solution(frequency: float) ngsolve.comp.GridFunction[source]
Evaluate electrical potential of volume conductor.
- Parameters:
frequency (float) – Frequency [Hz] of the input signal.
- Returns:
Data object representing the potential of volume conductor and floating values of floating contacts.
- Return type:
Potential
- class ossdbs.fem.volume_conductor.floating.VolumeConductorFloating(geometry: ModelGeometry, conductivity: ConductivityCF, solver: Solver, order: int, meshing_parameters: dict, output_path: str = 'Results')[source]
Bases:
VolumeConductorVolume conductor with floating conductors.
- compute_solution(frequency: float) ngsolve.comp.GridFunction[source]
Evaluate electrical potential of volume conductor.
- Parameters:
frequency (float) – Frequency [Hz] of the input signal.
- Returns:
Data object representing the potential of volume conductor and floating values of floating contacts.
- Return type:
Potential
- class ossdbs.fem.volume_conductor.floating_impedance.VolumeConductorFloatingImpedance(geometry: ModelGeometry, conductivity: ConductivityCF, solver: Solver, order: int, meshing_parameters: dict, output_path: str = 'Results')[source]
Bases:
VolumeConductorVolume conductor with floating conductors and surface impedances.
- compute_solution(frequency: float) ngsolve.comp.GridFunction[source]
Evaluate electrical potential of volume conductor.
- Parameters:
frequency (float) – Frequency [Hz] of the input signal.
- Returns:
Data object representing the potential of volume conductor and floating values of floating contacts.
- Return type:
Potential
Material parameters
- class ossdbs.fem.volume_conductor.conductivity.ConductivityCF(mri_image: MagneticResonanceImage, brain_bounding_box: BoundingBox, dielectric_properties: dict, materials: dict, encapsulation_layers: list | None = None, complex_data: bool = False, dti_image: DiffusionTensorImage | None = None, wm_masking: bool = True)[source]
Bases:
objectConductivity wrapper.
- create_dti_voxel_cf(dti_data, dti_voxel_bounding_box, dti_image)[source]
Map DTI image on NGSolve VoxelCoefficient function.
Notes
The DTI images assume symmetry of the tensor. This is not (yet) the case in NGSolve.
- property dielectric_properties: dict
Return dielectric properties.
- property dti_voxel_distribution: ngsolve.VoxelCoefficient
Return DTI data as VoxelCoefficient before scaling by conductivity.
- property is_complex: bool
Return if complex arithmetic is used.
- property is_tensor: bool
Return if conductivity is a tensor.
- material_distribution(mesh: Mesh) ngsolve.VoxelCoefficient[source]
Return MRI image projected onto mesh.
Notes
The encapsulation layer material is set to a fixed value.
- property materials: dict
Return dictionary of material names.
Mesh
- class ossdbs.fem.mesh.Mesh(geometry: netgen.occ.OCCGeometry, order: int)[source]
Bases:
objectClass for interacting with the mesh for FEM.
- Parameters:
geometry (netgen.occ.OCCGeometry) –
order (int) – Order of mesh elements.
complex_datatype (bool) – True for complex data type, False otherwise.
- apply_hp_refinement() None[source]
Apply deferred HP refinement.
Must be called after all bisection-based refinement steps (e.g. material refinement) are complete.
- property boundaries: list
Get list of boundary names.
- boundary_coefficients(boundaries: dict) ngsolve.fem.CoefficientFunction[source]
Return a boundary coefficient function.
- Return type:
ngsolve.fem.CoefficientFunction
- curve(order: int) None[source]
Curve mesh and overwrite mesh order.
- Parameters:
order (int) – Order of curved mesh
- property geometry: netgen.occ.OCCGeometry
Underlying CAD geometry of mesh.
- get_meshing_parameters(mesh_parameters: dict)[source]
Prepare NGSolve meshing parameters deviating from default.
- property hp_refinement_applied: bool
Whether HP refinement has been applied to this mesh.
- material_coefficients(materials: dict) ngsolve.fem.CoefficientFunction[source]
Return a boundary coefficient function.
- Return type:
ngsolve.fem.CoefficientFunction
- property materials: list
Get list of material names.
- property n_elements: int
Number of elements.
- property ngsolvemesh: ngsolve.Mesh
Return mesh as a ngsolve object.
- Return type:
ngsolve.comp.Mesh
- not_included(points: ndarray) ndarray[source]
Check each point in collection for collision with geometry. True if point is included in geometry, false otherwise.
- Parameters:
points (np.ndarray) – Array of point coordinates (x, y, z).
- Returns:
Array representing the state of collision for each point. True if point is included in geometry, False otherwise.
- Return type:
np.ndarray
- property order: int
Order of curved mesh.
- refine(at_surface: bool = False) None[source]
Refine the mesh.
- Parameters:
at_surface (bool) – Whether to mark surface elements.
- refine_at_boundaries(boundaries: list) None[source]
Refine the mesh by the boundaries.
- Parameters:
boundaries (list of str) – Collection of boundary names.
- refine_by_error_cf(error_cf: ngsolve.GridFunction) list[source]
Refine the mesh by the error at each mesh element.
- Parameters:
error_cf (ngsolve.GridFunction) – Estimated error from
- refine_by_material_cf(material_cf: ngsolve.GridFunction) list[source]
Refine the mesh by checking the material cf. Each element-wise integral divided by area should yield an integer value. If not, it needs to be refined.
- Parameters:
material_cf (ngsolve.GridFunction) – CF with material indices (need to be integers)
Solver
- class ossdbs.fem.solver.CGSolver(precond_par: ~ossdbs.fem.preconditioner.Preconditioner = <ossdbs.fem.preconditioner.BDDCPreconditioner object>, maxsteps: int = 10000, precision: float = 1e-12)[source]
Bases:
SolverConjugate gradient solver.
- bvp(bilinear_form: ngsolve.BilinearForm, linear_form: ngsolve.LinearForm, grid_function: ngsolve.GridFunction) None[source]
Solve boundary-value problem.
- Parameters:
bilinear_form (ngsolve.BilinearForm) – Bilinear form (left-hand side)
linear_form (ngsolve.LinearForm) – Linear form (right-hand side)
grid_function (ngsolve.GridFunction) – Solution vector
- class ossdbs.fem.solver.DirectSolver(precond_par: ~ossdbs.fem.preconditioner.Preconditioner = <ossdbs.fem.preconditioner.BDDCPreconditioner object>, maxsteps: int = 10000, precision: float = 1e-12)[source]
Bases:
SolverDirect solver.
- bvp(bilinear_form: ngsolve.BilinearForm, linear_form: ngsolve.LinearForm, grid_function: ngsolve.GridFunction) None[source]
Solve boundary-value problem.
- Parameters:
bilinear_form (ngsolve.BilinearForm) – Bilinear form (left-hand side)
linear_form (ngsolve.LinearForm) – Linear form (right-hand side)
grid_function (ngsolve.GridFunction) – Solution vector
Notes
TODO have a second look
- class ossdbs.fem.solver.GMRESSolver(precond_par: ~ossdbs.fem.preconditioner.Preconditioner = <ossdbs.fem.preconditioner.BDDCPreconditioner object>, maxsteps: int = 10000, precision: float = 1e-12)[source]
Bases:
SolverGMRes solver.
- bvp(bilinear_form: ngsolve.BilinearForm, linear_form: ngsolve.LinearForm, grid_function: ngsolve.GridFunction) None[source]
Solve boundary-value problem.
- Parameters:
bilinear_form (ngsolve.BilinearForm) – Bilinear form (left-hand side)
linear_form (ngsolve.LinearForm) – Linear form (right-hand side)
grid_function (ngsolve.GridFunction) – Solution vector
- class ossdbs.fem.solver.Solver(precond_par: ~ossdbs.fem.preconditioner.Preconditioner = <ossdbs.fem.preconditioner.BDDCPreconditioner object>, maxsteps: int = 10000, precision: float = 1e-12)[source]
Bases:
ABCComputes the boundary volume problem.
- DEFAULT_PRECONDITIONER = <ossdbs.fem.preconditioner.BDDCPreconditioner object>
- abstract bvp(bilinear_form: ngsolve.BilinearForm, linear_form: ngsolve.LinearForm, grid_function: ngsolve.GridFunction) None[source]
Solve boundary-value problem.
- Parameters:
bilinear_form (ngsolve.BilinearForm) – Bilinear form (left-hand side)
linear_form (ngsolve.LinearForm) – Linear form (right-hand side)
grid_function (ngsolve.GridFunction) – Solution vector
Preconditioner
- class ossdbs.fem.preconditioner.AMGPreconditioner[source]
Bases:
PreconditionerAlgebraic multigrid preconditioner.
- class ossdbs.fem.preconditioner.BDDCPreconditioner(coarsetype: str = 'h1amg')[source]
Bases:
PreconditionerBDDC preconditioner, highly efficient for higher-order FEM.
- class ossdbs.fem.preconditioner.CustomizedLocalPreconditioner[source]
Bases:
PreconditionerCustom diagonal Jacobi preconditioner.
- class ossdbs.fem.preconditioner.DirectPreconditioner(inverse: str = '')[source]
Bases:
PreconditionerDirect solver.
- class ossdbs.fem.preconditioner.JacobiPreconditioner(*args: Any, **kwargs: Any)[source]
Bases:
BaseMatrixSimple diagonal Jacobi preconditioner as an NGSolve operator.
- property height: int
Return the operator height.
- property width: int
Return the operator width.
- class ossdbs.fem.preconditioner.LocalPreconditioner[source]
Bases:
PreconditionerJacobi preconditioner.
- class ossdbs.fem.preconditioner.MultigridPreconditioner[source]
Bases:
PreconditionerGeometric multigrid preconditioner.