pulse package

Submodules

pulse.dolfin_utils module

class pulse.dolfin_utils.MixedParameter(fun, n, name='material_parameters')

Bases: dolfin.function.function.Function

assign_sub(f, i)

Assign subfunction

Parameters
  • f – The function you want to assign

  • i (int) – The subspace number

pulse.dolfin_utils.QuadratureSpace(mesh, degree, dim=3)

From FEniCS version 1.6 to 2016.1 there was a change in how FunctionSpace is defined for quadrature spaces. This functions checks your dolfin version and returns the correct quadrature space

Parameters
  • mesh (dolfin.Mesh) – The mesh

  • degree (int) – The degree of the element

  • dim (int) – For a mesh of topological dimension 3, dim = 1 would be a scalar function, and dim = 3 would be a vector function.

Returns

The quadrature space

Return type

dolfin.FunctionSpace

class pulse.dolfin_utils.RegionalParameter(meshfunction)

Bases: dolfin.function.function.Function

A regional paramerter defined in terms of a dolfin.MeshFunction

Suppose you have a MeshFunction defined different regions in your mesh, and you want to define different parameters on different regions, then this is what you want.

property function

Return linear combination of coefficents and basis functions

Returns

A function with parameter values at each segment specified by the meshfunction

Return type

:py:class`dolfin.Function

get_values()
property proj_space

Space for projecting the scalars. This is a DG 0 space.

pulse.dolfin_utils.compute_meshvolume(domain=None, dx=Measure('cell', subdomain_id='everywhere'), subdomain_id=None)
pulse.dolfin_utils.get_cavity_volume(geometry, unload=False, chamber='lv', u=None, xshift=0.0)
pulse.dolfin_utils.get_cavity_volume_form(mesh, u=None, xshift=0.0)
pulse.dolfin_utils.get_constant(val, value_size=None, value_rank=0, constant=<class 'dolfin.function.constant.Constant'>)
pulse.dolfin_utils.get_dimesion(u)
pulse.dolfin_utils.get_pressure(problem)

Returns p_lv (and p_rv if BiV mesh)

pulse.dolfin_utils.heaviside(x)

Heaviside function

\[\frac{\mathrm{d}}{\mathrm{d}x} \max\{x,0\}\]
pulse.dolfin_utils.list_sum(lst)

Return the sum of a list, when the convetiional method (like sum) it not working. For example if you have a list of dolfin functions.

Parameters

l (list) – a list of objects

Returns

The sum of the list. The type depends on the type of elemets in the list

pulse.dolfin_utils.map_vector_field(f0, new_mesh, u=None, name='fiber', normalize=True)

Map a vector field (f0) onto a new mesh (new_mesh) where the new mesh can be a moved version of the original one according to some displacement (u). In that case we will just a Piola transform to map the vector field.

pulse.dolfin_utils.normalize_vector_field(u)

Given a vector field, return a vector field with an L2 norm equal to 1.0

pulse.dolfin_utils.subplus(x)

Ramp function

\[\max\{x,0\}\]
pulse.dolfin_utils.update_function(mesh, f)

Given a function \(f\) defined on some domain, update the function so that it now is defined on the domain given in the mesh

pulse.dolfin_utils.vectorfield_to_components(u, S, dim)

pulse.geometry module

class pulse.geometry.CRLBasis(c0: Optional[dolfin.function.function.Function] = None, r0: Optional[dolfin.function.function.Function] = None, l0: Optional[dolfin.function.function.Function] = None)

Bases: object

l0 - longitudinal, c0 - circumferential, r0 - radial

c0: Optional[dolfin.function.function.Function] = None
l0: Optional[dolfin.function.function.Function] = None
r0: Optional[dolfin.function.function.Function] = None
class pulse.geometry.Geometry(mesh, markers=None, marker_functions=None, microstructure=None, crl_basis=None)

Bases: object

Base class for geometry

Parameters
  • mesh (dolfin.mesh) – The mesh

  • markers (dict) – A dictionary with markers for the mesh (optional)

  • marker_functions (pulse.geometry.MarkerFunctions) – A Markerfunction object with Meshfunctions for the mesh

  • microstructure (pulse.geometry.Microstructure) – A Markerfunction object with functions for the fiber, sheet and sheet normal (optional)

  • crl_basis (pulse.geometry.CRLBasis) – A CRLBasis objedt with funcions for the circumferential, radial and longitudinal vectors. (optional)

Geometry can be intanciated directly

Example

import dolfin
mesh = dolfin.UnitCubeMesh(3,3,3)
geometry = Geometry(mesh)

You can load create and instace be loading a geometry from a file

Example

# Geometry is stored in a file "geometry.h5"
geo = Geometry.from_file("geometry.h5")
property c0

Circumferential

property cfun

Cell Function

copy(u=None, factor=1.0)

Return as copy of the geometry with possiblity to move the new geometry according to some factor times some displacement (u)

crl_basis_list()

Return a list of the CRL basis in the order c0, r0, l0. Basis elements that are none will not be included

property ds

Return the surface measure of exterior facets using self.mesh as domain and self.ffun as subdomain_data

property dx

Return the volume measure using self.mesh as domain and self.cfun as subdomain_data

property efun

Edge Function

property f0

Fibers

property facet_normal
property ffun

Facet Function

classmethod from_file(h5name, h5group='', comm=None)
property geo_dim

Geometric Dimension

property l0

Longitudinal

static load_from_file(h5name, h5group, comm)
meshfunction_list()
property meshvol

Return the volume of the whole mesh

property meshvols

Return a list of the volume of each subdomain defined by the cell function (cfun)

microstructure_list()

Fibers, sheet and sheet-normals in a list

property n0

Cross-Sheets

property nbasis
property nregions
property r0

Radial

property regions

Return a set of the unique values of the cell function (cfun)

property s0

Sheets

save(h5name, h5group='', other_functions=None, other_attributes=None, overwrite_file=False, overwrite_group=True)
property top_dim

Topological Dimension

property vfun

Vertex Function

class pulse.geometry.HeartGeometry(mesh, markers=None, marker_functions=None, microstructure=None, crl_basis=None)

Bases: pulse.geometry.Geometry

property base_mean_position

Return mean coordinates of the base (serial only?)

cavity_volume(chamber='lv', u=None)
property is_biv
update_xshift()
class pulse.geometry.MarkerFunctions(vfun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None, efun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None, ffun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None, cfun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None)

Bases: object

# vfun - vertex function # efun - edge function # ffun - facet function # cfun - cell function

cfun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None
efun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None
ffun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None
vfun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None
class pulse.geometry.Microstructure(f0: Optional[dolfin.function.function.Function] = None, s0: Optional[dolfin.function.function.Function] = None, n0: Optional[dolfin.function.function.Function] = None)

Bases: object

f0 - fibers, s0 - sheets, n0 - cross-sheets

f0: Optional[dolfin.function.function.Function] = None
n0: Optional[dolfin.function.function.Function] = None
s0: Optional[dolfin.function.function.Function] = None
pulse.geometry.get_attribute(obj, key1, key2, default=None)

pulse.geometry_utils module

pulse.geometry_utils.get_geometric_dimension(mesh)
pulse.geometry_utils.ggroup(h5group)
pulse.geometry_utils.load_geometry_from_h5(h5name, h5group='', include_sheets=True, comm=<mpi4py.MPI.Intracomm object>)

Load geometry and other mesh data from a h5file to an object. If the file contains muliple fiber fields you can spefify the angles, and if the file contais sheets and cross-sheets this can also be included

Parameters
  • h5name (str) – Name of the h5file

  • h5group (str) – The group within the file

  • include_sheets (bool) – Include sheets and cross-sheets

Returns

An object with geometry data

Return type

object

pulse.geometry_utils.load_local_basis(h5file, lgroup, mesh, geo)
pulse.geometry_utils.load_markers(h5file, mesh, ggroup, dgroup)
pulse.geometry_utils.load_microstructure(h5file, fgroup, mesh, geo, include_sheets=True)
pulse.geometry_utils.mgroup(h5group)
pulse.geometry_utils.move(mesh, u, factor=1.0)

Move mesh according to some displacement times some factor

pulse.geometry_utils.save_fields(h5file, h5group, fields)

Save mictrostructure to h5file

Parameters
  • h5file (dolfin.HDF5File) – The file to write to

  • h5group (str) – The folder inside the file to write to

  • fields (list) – List of dolfin functions to write

pulse.geometry_utils.save_geometry_to_h5(mesh, h5name, h5group='', markers=None, fields=None, local_basis=None, meshfunctions=None, comm=None, other_functions=None, other_attributes=None, overwrite_file=False, overwrite_group=True)

Save geometry and other geometrical functions to a HDF file.

Parameters
  • mesh (dolfin.mesh) – The mesh

  • h5name (str) – Path to the file

  • h5group (str) – Folder within the file. Default is “” which means in the top folder.

  • markers (dict) – A dictionary with markers. See get_markers.

  • fields (list) – A list of functions for the microstructure

  • local_basis (list) – A list of functions for the crl basis

  • meshfunctions (dict) – A dictionary with keys being the dimensions the the values beeing the meshfunctions.

  • comm (dolfin.MPI) – MPI communicator

  • other_functions (dict) – Dictionary with other functions you want to save

  • other_attributes (dict) – Dictionary with other attributes you want to save

  • overwrite_file (bool) – If true, and the file exists, the file will be overwritten (default: False)

  • overwrite_group (bool) – If true and h5group exist, the group will be overwritten.

pulse.geometry_utils.save_local_basis(h5file, h5group, local_basis)

Save local basis functions

Parameters
  • h5file (dolfin.HDF5File) – The file to write to

  • h5group (str) – Folder inside file to store the local basis

  • local_basis (list) – List of dolfin functions with local basis functions

pulse.geometry_utils.save_markers(h5file, h5group, markers)
pulse.geometry_utils.save_meshfunctions(h5file, h5group, mesh, meshfunctions, comm)
pulse.geometry_utils.save_other_attributes(h5file, h5group, other_attributes)
pulse.geometry_utils.save_other_functions(h5file, h5group, other_functions)

pulse.io_utils module

pulse.io_utils.check_h5group(h5name, h5group, delete=False, comm=<mpi4py.MPI.Intracomm object>)
pulse.io_utils.open_h5py(h5name, file_mode='a', comm=<mpi4py.MPI.Intracomm object>)
pulse.io_utils.read_h5file(h5file, obj, group, *args, **kwargs)

pulse.iterate module

class pulse.iterate.Iterator(problem, control, target, old_states=None, old_controls=None, **params)

Bases: object

assign_control(new_control)

Assign a new value to the control

change_step_for_final_iteration(prev_control)

Change step size so that target is reached in the next iteration

change_step_size(factor)
continuation_step()
static default_parameters()
increment_control()
property ncontrols

Number of controls

print_control()
solve()
property step
target_reached()

Check if control and target are the same

pulse.iterate.constant2float(const)

Convert a dolfin.Constant to float

pulse.iterate.copy(f, deepcopy=True, name='copied_function')

Copy a function. This is to ease the integration with dolfin adjoint where copied fuctions are annotated.

pulse.iterate.get_delta(new_control, c0, c1)

Get extrapolation parameter used in the continuation step.

pulse.iterate.get_diff(current, target)

Get difference between current and target value

pulse.iterate.get_initial_step(current, target, nsteps=5)

Estimate the step size needed to step from current to target in nsteps.

pulse.iterate.iterate(problem, control, target, continuation=True, max_adapt_iter=8, adapt_step=True, old_states=None, old_controls=None, max_nr_crash=20, max_iters=40, initial_number_of_steps=5, reinit_each_step=False)

Using the given problem, iterate control to given target.

Parameters
  • problem (pulse.MechanicsProblem) – The problem

  • control (dolfin.Function or dolfin.Constant) – The control

  • target (dolfin.Function, dolfin.Constant, tuple or float) – The target value. Typically a float if target is LVP, a tuple if target is (LVP, RVP) and a function if target is gamma.

  • continuation (bool) – Apply continuation for better guess for newton problem Note: Replay test seems to fail when continuation is True, but taylor test passes

  • max_adapt_iter (int) – If number of iterations is less than this number and adapt_step=True, then adapt control step. Default: 8

  • adapt_step (bool) – Adapt / increase step size when sucessful iterations are achevied.

  • old_states (list) – List of old controls to help speed in the continuation

  • reinit_each_step (bool) – If True reinitialize form at each step.

pulse.iterate.np2str(c_, fmt='{:.3f}')
pulse.iterate.print_control(cs, msg)
pulse.iterate.squeeze(x)
pulse.iterate.step_too_large(current, target, step)

Check if current + step exceeds target

pulse.kinematics module

pulse.kinematics.DeformationGradient(u, isochoric=False)

Return deformation gradient from displacement.

pulse.kinematics.EngineeringStrain(F, isochoric=False)

Equivalent of engineering strain

pulse.kinematics.EulerAlmansiStrain(F, isochoric=False)

Euler-Almansi strain tensor

pulse.kinematics.GreenLagrangeStrain(F, isochoric=False)

Green-Lagrange strain tensor

pulse.kinematics.I1(F, isochoric=False)
pulse.kinematics.I2(F, isochoric=False)
pulse.kinematics.I3(F, isochoric=False)
pulse.kinematics.I4(F, a0=None, isochoric=False)
pulse.kinematics.I5(F, a0, isochoric=False)
pulse.kinematics.I6(F, b0)
pulse.kinematics.I7(F, b0)
pulse.kinematics.I8(F, a0, b0)
pulse.kinematics.InversePiolaTransform(A, F)
pulse.kinematics.IsochoricDeformationGradient(F)

Return the isochoric deformation gradient

pulse.kinematics.Jacobian(F)

Determinant of the deformation gradient

pulse.kinematics.LeftCauchyGreen(F, isochoric=False)

Left Cauchy-Green tensor

pulse.kinematics.PiolaTransform(A, F)
pulse.kinematics.RightCauchyGreen(F, isochoric=False)

Left Cauchy-Green tensor

pulse.kinematics.SecondOrderIdentity(F)

Return identity with same dimension as input

pulse.mechanicsproblem module

class pulse.mechanicsproblem.BoundaryConditions(neumann: Sequence[pulse.mechanicsproblem.NeumannBC] = (), dirichlet: Sequence[Union[Callable[[dolfin.function.functionspace.FunctionSpace], dolfin.fem.dirichletbc.DirichletBC], dolfin.fem.dirichletbc.DirichletBC]] = (), robin: Sequence[pulse.mechanicsproblem.RobinBC] = (), body_force: Sequence[ufl.coefficient.Coefficient] = ())

Bases: object

body_force: Sequence[ufl.coefficient.Coefficient] = ()
dirichlet: Sequence[Union[Callable[[dolfin.function.functionspace.FunctionSpace], dolfin.fem.dirichletbc.DirichletBC], dolfin.fem.dirichletbc.DirichletBC]] = ()
neumann: Sequence[pulse.mechanicsproblem.NeumannBC] = ()
robin: Sequence[pulse.mechanicsproblem.RobinBC] = ()
class pulse.mechanicsproblem.MechanicsProblem(geometry: pulse.geometry.Geometry, material: pulse.material.material_model.Material, bcs=None, bcs_parameters=None, solver_parameters=None)

Bases: object

Base class for mechanics problem

ChachyStress()
FirstPiolaStress()
GreenLagrangeStrain()
SecondPiolaStress()
static default_bcs_parameters()
static default_solver_parameters()
get_displacement(annotate=False)
reinit(state, annotate=False)

Reinitialze state

solve()

Solve the variational problem

\[\delta W = 0\]
class pulse.mechanicsproblem.NeumannBC(traction: Union[float, ufl.coefficient.Coefficient], marker: int, name: str = '')

Bases: object

marker: int
name: str = ''
traction: Union[float, ufl.coefficient.Coefficient]
class pulse.mechanicsproblem.RobinBC(value: Union[float, ufl.coefficient.Coefficient], marker: int)

Bases: object

marker: int
value: Union[float, ufl.coefficient.Coefficient]
exception pulse.mechanicsproblem.SolverDidNotConverge

Bases: Exception

pulse.mechanicsproblem.cardiac_boundary_conditions(geometry, pericardium_spring=0.0, base_spring=0.0, base_bc='fix_x', **kwargs)
pulse.mechanicsproblem.dirichlet_fix_base(W, ffun, marker)

Fix the basal plane.

pulse.mechanicsproblem.dirichlet_fix_base_directional(W, ffun, marker, direction=0)

pulse.numpy_mpi module

These functions are copied from cbcpost https://bitbucket.org/simula_cbc/cbcpost

pulse.numpy_mpi.assign_to_vector(v, a)

Assign the value of the array a to the dolfin vector v

pulse.numpy_mpi.broadcast(array, from_process)

Broadcast array to all processes

pulse.numpy_mpi.compile_extension_module(cpp_code, **kwargs)
pulse.numpy_mpi.distribution(number)

Get distribution of number on all processes

pulse.numpy_mpi.gather(array, on_process=0, flatten=False)

Gather array from all processes on a single process

pulse.numpy_mpi.gather_broadcast(arr)
pulse.numpy_mpi.gather_vector(u, size=None)
pulse.numpy_mpi.mpi_print(mess)

pulse.solver module

class pulse.solver.NonlinearProblem(J, F, bcs, output_matrix=False, output_matrix_path='output', **kwargs)

Bases: dolfin.cpp.nls.NonlinearProblem

F(self: dolfin.cpp.nls.NonlinearProblem, arg0: dolfin.cpp.la.GenericVector, arg1: dolfin.cpp.la.GenericVector) None
J(self: dolfin.cpp.nls.NonlinearProblem, arg0: dolfin.cpp.la.GenericMatrix, arg1: dolfin.cpp.la.GenericVector) None
class pulse.solver.NonlinearSolver(problem: pulse.solver.NonlinearProblem, state, parameters=None)

Bases: object

static default_solver_parameters()
solve()

Solve the problem.

Returns

residual – A measure of the accuracy (convergence and error) of the performed computation.

Return type

_solver.snes (???)

update_parameters(parameters)
pulse.solver.dump_matrix_to_mtx(A: dolfin.cpp.la.PETScMatrix, filename: pathlib.Path) None

pulse.unloader module

This script contains the code for unloading geometry

class pulse.unloader.FixedPointUnloader(problem, pressure, h5name='test.h5', options=None, h5group='', overwrite=False, merge_control='')

Bases: pulse.unloader.MeshUnloader

Assumes that the given geometry is loaded with the given pressure.

Find the reference configuration corresponding to zero pressure using a backward displacement method, described in [2]_. This method assumes that material properties are given. Make sure to change this in utils.py

This also runs in parallel.

References

2

Bols, Joris, et al. “A computational method to assess the in vivo stresses and unloaded configuration of patient-specific blood vessels.” Journal of computational and Applied mathematics 246 (2013): 10-17.

unload_step(u, residual, save=False, it=0)

Unload step

class pulse.unloader.MeshUnloader(problem, pressure, h5name='test.h5', options=None, h5group='', overwrite=False, merge_control='')

Bases: object

property backward_displacement

Return the current backward displacement as a function on the original geometry.

default_parameters()

Default parameters.

initial_solve(save=False)

Inflate the original geometry to the target pressure and return the displacement field.

save(obj, name, h5group='')

Save object to and HDF file.

Parameters
  • obj (dolfin.Mesh or dolfin.Function) – The object you want to save

  • name (str) – Name of the object

  • h5group (str) – The folder you want to save the object withing the HDF file. Default: ‘’

unload(save=False)

Unload the geometry

property unloaded_geometry
class pulse.unloader.RaghavanUnloader(problem, pressure, h5name='test.h5', options=None, h5group='', overwrite=False, merge_control='')

Bases: pulse.unloader.MeshUnloader

Assumes that the given geometry is loaded with the given pressure.

Find the reference configuration corresponding to zero pressure by first inflating the original geometry to the target pressure and subtacting the displacement times some factor k. Use a 1D minimization algorithm from scipy to find k. The method is described in [1]_ This method assumes that material properties are given. Make sure to change this in utils.py

This also runs in parallel.

Parameters
  • geometry (object) – An obeject with attributes mesh, fiber and markers.

  • pressure (float or tuple) – The target pressure. If LV only provide a float, if BiV proble a tuple of the form (p_LV, p_RV).

  • material_parameters (dict, optional) – A dictionary with material parameters

  • h5name (str, optional) – Path to where you want to save the output Default: test.h5

  • options (dict, optional, optional) –

    A dictionary of solver options. All methods accept the following generic options:

    maxiterint

    Maximum number of iterations to perform.

    tolfloat

    Tolerance of termination

    solve_triesint

    Number of attemtps the solver should use to increase the pressure (after pressure reduction)

  • h5group (str) – Subfolder within the HDF file where you save the results

  • remove_old (bool) – If the file with the same h5name exist, delete it before starting

  • params (dict) – These are parameters that is parsed to pulse_adjoint that make the parameters for the solver. See pulse_adjoint.setup_parametere.setup_adjoint_contraction_parameters

References

1

Raghavan, M. L., Baoshun Ma, and Mark F. Fillinger. “Non-invasive determination of zero-pressure geometry of arterial aneurysms.” Annals of biomedical engineering 34.9 (2006): 1414-1419.

unload_step(u, residual, save=True)
pulse.unloader.setup_unloading_parameters()

Parameters for coupled unloading/material parameter estimation.

For info about the different parameters, see the unloading module.

pulse.unloading_utils module

class pulse.unloading_utils.ResidualCalculator(mesh)

Bases: object

calculate_residual(mesh2)
pulse.unloading_utils.continuation_step(params, it_, paramvec)
pulse.unloading_utils.inflate_to_pressure(target_pressure, problem, ntries=5, n=2, annotate=False)
pulse.unloading_utils.load_material_parameter(h5name, h5group, paramvec)
pulse.unloading_utils.load_opt_target(h5name, h5group, key='volume', data='simulated')
pulse.unloading_utils.print_volumes(geometry, logger=<KeywordArgumentAdapter pulse.unloading_utils (INFO)>, txt='original', u=None)
pulse.unloading_utils.solve(target_pressure, problem, pressure, ntries=5, n=2, annotate=False)

pulse.utils module

class pulse.utils.Annotation

Bases: object

Object holding global annotation for dolfin-adjoint

property annotate
class pulse.utils.Enlisted(iterable=(), /)

Bases: tuple

class pulse.utils.MPIFilt(name='')

Bases: logging.Filter

filter(record)

Determine if the specified record is to be logged.

Returns True if the record should be logged, or False otherwise. If deemed appropriate, the record may be modified in-place.

class pulse.utils.Object

Bases: object

exception pulse.utils.UnableToChangePressureExeption

Bases: Exception

pulse.utils.delist(x)
pulse.utils.enlist(x, force_enlist=False)
pulse.utils.getLogger(name)
pulse.utils.get_dolfin_version()
pulse.utils.get_lv_marker(geometry)
pulse.utils.mpi_comm_world()
pulse.utils.value_size(obj)

Module contents

class pulse.ActiveModel(activation=None, f0=None, s0=None, n0=None, model: pulse.material.active_model.ActiveModels = ActiveModels.active_strain, active_isotropy: pulse.material.active_model.ActiveStressModels = ActiveStressModels.transversally, eta: Optional[float] = None, T_ref: Optional[float] = None, isochoric=True)

Bases: object

property Fa
Fe(F)
Wactive(F=None, diff=0)

Active stress energy

property activation

Return the activation paramter. If regional, this will return one parameter for each segment.

property activation_field

Return the activation field. If regional, this will return a piecewise constant function (DG_0)

property dim
class pulse.ActiveModels(value)

Bases: str, enum.Enum

An enumeration.

active_strain = 'active_strain'
active_stress = 'active_stress'
class pulse.BoundaryConditions(neumann: Sequence[pulse.mechanicsproblem.NeumannBC] = (), dirichlet: Sequence[Union[Callable[[dolfin.function.functionspace.FunctionSpace], dolfin.fem.dirichletbc.DirichletBC], dolfin.fem.dirichletbc.DirichletBC]] = (), robin: Sequence[pulse.mechanicsproblem.RobinBC] = (), body_force: Sequence[ufl.coefficient.Coefficient] = ())

Bases: object

body_force: Sequence[ufl.coefficient.Coefficient] = ()
dirichlet: Sequence[Union[Callable[[dolfin.function.functionspace.FunctionSpace], dolfin.fem.dirichletbc.DirichletBC], dolfin.fem.dirichletbc.DirichletBC]] = ()
neumann: Sequence[pulse.mechanicsproblem.NeumannBC] = ()
robin: Sequence[pulse.mechanicsproblem.RobinBC] = ()
class pulse.CRLBasis(c0: Optional[dolfin.function.function.Function] = None, r0: Optional[dolfin.function.function.Function] = None, l0: Optional[dolfin.function.function.Function] = None)

Bases: object

l0 - longitudinal, c0 - circumferential, r0 - radial

c0: Optional[dolfin.function.function.Function] = None
l0: Optional[dolfin.function.function.Function] = None
r0: Optional[dolfin.function.function.Function] = None
class pulse.Constant(value, cell=None, name=None)

Bases: ufl.coefficient.Coefficient

assign(x)
cell()
compute_vertex_values(mesh)
cpp_object()
eval(*args, **kwargs)
id()
name()
rename(name, s)
str(verbose)
value_size()
values()
pulse.DeformationGradient(u, isochoric=False)

Return deformation gradient from displacement.

class pulse.DirichletBC(*args, **kwargs)

Bases: dolfin.cpp.fem.DirichletBC

pulse.EulerAlmansiStrain(F, isochoric=False)

Euler-Almansi strain tensor

class pulse.FixedPointUnloader(problem, pressure, h5name='test.h5', options=None, h5group='', overwrite=False, merge_control='')

Bases: pulse.unloader.MeshUnloader

Assumes that the given geometry is loaded with the given pressure.

Find the reference configuration corresponding to zero pressure using a backward displacement method, described in [2]_. This method assumes that material properties are given. Make sure to change this in utils.py

This also runs in parallel.

References

2

Bols, Joris, et al. “A computational method to assess the in vivo stresses and unloaded configuration of patient-specific blood vessels.” Journal of computational and Applied mathematics 246 (2013): 10-17.

unload_step(u, residual, save=False, it=0)

Unload step

class pulse.Function(*args, **kwargs)

Bases: ufl.coefficient.Coefficient

assign(rhs)

Assign either a Function or linear combination of Functions.

Arguments
rhs (_Function_)

A Function or a linear combination of Functions. If a linear combination is passed all Functions need to be in the same FunctionSpaces.

compute_vertex_values(mesh=None)
copy(deepcopy=False)
cpp_object()
eval(u, x)
eval_cell(u, x, cell)
extrapolate(u)
function_space()

Return the FunctionSpace

get_allow_extrapolation()
id()
interpolate(u)
leaf_node()
name()
num_sub_spaces()
part()
rename(name, s)
restrict(element, cell)

Returns expansion coefficients of function restricted to local cell.

Arguments
elementcpp.fem.FiniteElement

The element.

cellCell

The cell.

root_node()
set_allow_extrapolation(value)
split(deepcopy=False)

Extract any sub functions.

A sub function can be extracted from a discrete function that is in a mixed, vector, or tensor FunctionSpace. The sub function resides in the subspace of the mixed space.

Arguments
deepcopy

Copy sub function vector instead of sharing

sub(i, deepcopy=False)

Return a sub function.

The sub functions are numbered from i = 0..N-1, where N is the total number of sub spaces.

Arguments
iint

The number of the sub function

ufl_evaluate(x, component, derivatives)

Function used by ufl to evaluate the Expression

value_dimension(i)
value_rank()
value_shape()
vector()
class pulse.FunctionAssigner

Bases: pybind11_builtins.pybind11_object

assign(*args, **kwargs)

Overloaded function.

  1. assign(self: dolfin.cpp.function.FunctionAssigner, arg0: dolfin.cpp.function.Function, arg1: dolfin.cpp.function.Function) -> None

  2. assign(self: dolfin.cpp.function.FunctionAssigner, arg0: object, arg1: object) -> None

class pulse.Geometry(mesh, markers=None, marker_functions=None, microstructure=None, crl_basis=None)

Bases: object

Base class for geometry

Parameters
  • mesh (dolfin.mesh) – The mesh

  • markers (dict) – A dictionary with markers for the mesh (optional)

  • marker_functions (pulse.geometry.MarkerFunctions) – A Markerfunction object with Meshfunctions for the mesh

  • microstructure (pulse.geometry.Microstructure) – A Markerfunction object with functions for the fiber, sheet and sheet normal (optional)

  • crl_basis (pulse.geometry.CRLBasis) – A CRLBasis objedt with funcions for the circumferential, radial and longitudinal vectors. (optional)

Geometry can be intanciated directly

Example

import dolfin
mesh = dolfin.UnitCubeMesh(3,3,3)
geometry = Geometry(mesh)

You can load create and instace be loading a geometry from a file

Example

# Geometry is stored in a file "geometry.h5"
geo = Geometry.from_file("geometry.h5")
property c0

Circumferential

property cfun

Cell Function

copy(u=None, factor=1.0)

Return as copy of the geometry with possiblity to move the new geometry according to some factor times some displacement (u)

crl_basis_list()

Return a list of the CRL basis in the order c0, r0, l0. Basis elements that are none will not be included

property ds

Return the surface measure of exterior facets using self.mesh as domain and self.ffun as subdomain_data

property dx

Return the volume measure using self.mesh as domain and self.cfun as subdomain_data

property efun

Edge Function

property f0

Fibers

property facet_normal
property ffun

Facet Function

classmethod from_file(h5name, h5group='', comm=None)
property geo_dim

Geometric Dimension

property l0

Longitudinal

static load_from_file(h5name, h5group, comm)
meshfunction_list()
property meshvol

Return the volume of the whole mesh

property meshvols

Return a list of the volume of each subdomain defined by the cell function (cfun)

microstructure_list()

Fibers, sheet and sheet-normals in a list

property n0

Cross-Sheets

property nbasis
property nregions
property r0

Radial

property regions

Return a set of the unique values of the cell function (cfun)

property s0

Sheets

save(h5name, h5group='', other_functions=None, other_attributes=None, overwrite_file=False, overwrite_group=True)
property top_dim

Topological Dimension

property vfun

Vertex Function

pulse.GreenLagrangeStrain(F, isochoric=False)

Green-Lagrange strain tensor

class pulse.Guccione(activation: Optional[ufl.coefficient.Coefficient] = None, parameters=None, active_model: Union[pulse.material.active_model.ActiveModel, pulse.material.active_model.ActiveModels] = ActiveModels.active_strain, T_ref: Optional[float] = None, eta: float = 0.0, isochoric: bool = True, compressible_model='incompressible', geometry=None, f0=None, s0=None, n0=None, active_isotropy: pulse.material.active_model.ActiveStressModels = ActiveStressModels.transversally, *args, **kwargs)

Bases: pulse.material.material_model.Material

Guccione material model.

SecondPiolaStress(F, p, *args, **kwargs)
static default_parameters()
is_isotropic()

Return True if the material is isotropic.

name = 'guccione'
strain_energy(F_)

UFL form of the strain energy.

class pulse.HeartGeometry(mesh, markers=None, marker_functions=None, microstructure=None, crl_basis=None)

Bases: pulse.geometry.Geometry

property base_mean_position

Return mean coordinates of the base (serial only?)

cavity_volume(chamber='lv', u=None)
property is_biv
update_xshift()
class pulse.HolzapfelOgden(activation: Optional[ufl.coefficient.Coefficient] = None, parameters=None, active_model: Union[pulse.material.active_model.ActiveModel, pulse.material.active_model.ActiveModels] = ActiveModels.active_strain, T_ref: Optional[float] = None, eta: float = 0.0, isochoric: bool = True, compressible_model='incompressible', geometry=None, f0=None, s0=None, n0=None, active_isotropy: pulse.material.active_model.ActiveStressModels = ActiveStressModels.transversally, *args, **kwargs)

Bases: pulse.material.material_model.Material

Orthotropic model by Holzapfel and Ogden

\[\mathcal{W}(I_1, I_{4f_0}) = \frac{a}{2 b} \left( e^{ b (I_1 - 3)} -1 \right) + \frac{a_f}{2 b_f} \left( e^{ b_f (I_{4f_0} - 1)_+^2} -1 \right) + \frac{a_s}{2 b_s} \left( e^{ b_s (I_{4s_0} - 1)_+^2} -1 \right) + \frac{a_fs}{2 b_fs} \left( e^{ b_fs I_{8fs}^2} -1 \right)\]

where

\[(\cdot)_+ = \max\{x,0\}\]

Reference

[1] Holzapfel, Gerhard A., and Ray W. Ogden. “Constitutive modelling of passive myocardium: a structurally based framework for material characterization. “Philosophical Transactions of the Royal Society of London A: Mathematical, Physical and Engineering Sciences 367.1902 (2009): 3445-3475.

W_1(I1, diff=0, *args, **kwargs)

Isotropic contribution.

If diff = 0, return

\[\frac{a}{2 b} \left( e^{ b (I_1 - 3)} -1 \right)\]

If diff = 1, return

\[\frac{a}{b} e^{ b (I_1 - 3)}\]

If diff = 2, return

\[\frac{a b}{2} e^{ b (I_1 - 3)}\]
W_4(I4, direction, diff=0, use_heaviside=False, *args, **kwargs)

Anisotropic contribution.

If diff = 0, return

\[\frac{a_f}{2 b_f} \left( e^{ b_f (I_{4f_0} - 1)_+^2} -1 \right)\]

If diff = 1, return

\[a_f (I_{4f_0} - 1)_+ e^{ b_f (I_{4f_0} - 1)^2}\]

If diff = 2, return

\[a_f h(I_{4f_0} - 1) (1 + 2b(I_{4f_0} - 1)) e^{ b_f (I_{4f_0} - 1)_+^2}\]

where

\[h(x) = \frac{\mathrm{d}}{\mathrm{d}x} \max\{x,0\}\]

is the Heaviside function.

W_8(I8, *args, **kwargs)

Cross fiber-sheet contribution.

static default_parameters()

Default matereial parameter for the Holzapfel Ogden model

Taken from Table 1 row 3 of [1]

name = 'holzapfel_ogden'
strain_energy(F)

Strain-energy density function.

\[\mathcal{W} = \mathcal{W}_1 + \mathcal{W}_{4f} + \mathcal{W}_{\mathrm{active}}\]

where

\[\begin{split}\mathcal{W}_{\mathrm{active}} = \begin{cases} 0 & \text{if acitve strain} \\ \gamma I_{4f} & \text{if active stress} \end{cases}\end{split}\]
Parameters

F (dolfin.Function) – Deformation gradient

pulse.InversePiolaTransform(A, F)
pulse.Jacobian(F)

Determinant of the deformation gradient

pulse.LeftCauchyGreen(F, isochoric=False)

Left Cauchy-Green tensor

class pulse.LinearElastic(activation: Optional[ufl.coefficient.Coefficient] = None, parameters=None, active_model: Union[pulse.material.active_model.ActiveModel, pulse.material.active_model.ActiveModels] = ActiveModels.active_strain, T_ref: Optional[float] = None, eta: float = 0.0, isochoric: bool = True, compressible_model='incompressible', geometry=None, f0=None, s0=None, n0=None, active_isotropy: pulse.material.active_model.ActiveStressModels = ActiveStressModels.transversally, *args, **kwargs)

Bases: pulse.material.material_model.Material

Class for linear elastic material

static default_parameters()
name = 'linear_elastic'
strain_energy(F_)
class pulse.MarkerFunctions(vfun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None, efun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None, ffun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None, cfun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None)

Bases: object

# vfun - vertex function # efun - edge function # ffun - facet function # cfun - cell function

cfun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None
efun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None
ffun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None
vfun: Optional[dolfin.mesh.meshfunction.MeshFunction] = None
class pulse.Material(activation: Optional[ufl.coefficient.Coefficient] = None, parameters=None, active_model: Union[pulse.material.active_model.ActiveModel, pulse.material.active_model.ActiveModels] = ActiveModels.active_strain, T_ref: Optional[float] = None, eta: float = 0.0, isochoric: bool = True, compressible_model='incompressible', geometry=None, f0=None, s0=None, n0=None, active_isotropy: pulse.material.active_model.ActiveStressModels = ActiveStressModels.transversally, *args, **kwargs)

Bases: abc.ABC

Initialize material model

Parameters
  • f0 (:py:class`dolfin.Function`) – Fiber field

  • gamma (:py:class`dolfin.Function`) – Activation parameter

  • params (dict) – Material parameters

  • active_model (str) – Active model - active strain or active stress

  • s0 (:py:class`dolfin.Function`) – Sheets

  • n0 (:py:class`dolfin.Function`) – Sheet - normal

  • T_ref (float) – Scale factor for activation parameter (default = 1.0)

  • dev_iso_split (bool) – Decouple deformation into deviatoric and isochoric deformations

  • eta (float) – Fraction of transverse active tesion for active stress formulation. 0 = active only along fiber, 1 = equal forces in all directions (default=0.0).

CauchyStress(F, p=None, deviatoric=False)
property Fa
Fe(F)
FirstPiolaStress(F, p=None, *args, **kwargs)
SecondPiolaStress(F, p=None, deviatoric=False, *args, **kwargs)
property T_ref
Wactive(F=None, diff=0)

Active stress energy

property activation

Return the activation paramter. If regional, this will return one parameter for each segment.

property activation_field

Return the activation field. If regional, this will return a piecewise constant function (DG_0)

property active_model
compressibility(p, J)
copy(geometry=None)
abstract static default_parameters()
property eta
property f0
property isochoric
property material_model
property n0
property name
property s0
abstract strain_energy(F)
update_geometry(geometry)
class pulse.MechanicsProblem(geometry: pulse.geometry.Geometry, material: pulse.material.material_model.Material, bcs=None, bcs_parameters=None, solver_parameters=None)

Bases: object

Base class for mechanics problem

ChachyStress()
FirstPiolaStress()
GreenLagrangeStrain()
SecondPiolaStress()
static default_bcs_parameters()
static default_solver_parameters()
get_displacement(annotate=False)
reinit(state, annotate=False)

Reinitialze state

solve()

Solve the variational problem

\[\delta W = 0\]
class pulse.Mesh

Bases: dolfin.cpp.common.Variable

DOLFIN Mesh object

bounding_box_tree(self: dolfin.cpp.mesh.Mesh) dolfin::BoundingBoxTree
build_mapping(self: dolfin.cpp.mesh.Mesh, arg0: dolfin.cpp.mesh.Mesh) None
cell_name(self: dolfin.cpp.mesh.Mesh) str
cell_orientations(self: dolfin.cpp.mesh.Mesh) List[int]
cells(self: dolfin.cpp.mesh.Mesh) numpy.ndarray
color(*args, **kwargs)

Overloaded function.

  1. color(self: dolfin.cpp.mesh.Mesh, arg0: str) -> List[int]

  2. color(self: dolfin.cpp.mesh.Mesh, arg0: List[int]) -> List[int]

coordinates(self: dolfin.cpp.mesh.Mesh) numpy.ndarray[numpy.float64[m, n], flags.writeable, flags.c_contiguous]
data(self: dolfin.cpp.mesh.Mesh) dolfin::MeshData

Data associated with a mesh

domains(self: dolfin.cpp.mesh.Mesh) dolfin::MeshDomains
geometric_dimension()

Returns geometric dimension for ufl interface

geometry(self: dolfin.cpp.mesh.Mesh) dolfin.cpp.mesh.MeshGeometry

Mesh geometry

hash(self: dolfin.cpp.mesh.Mesh) int
hmax(self: dolfin.cpp.mesh.Mesh) float
hmin(self: dolfin.cpp.mesh.Mesh) float
id(self: dolfin.cpp.mesh.Mesh) int
init(*args, **kwargs)

Overloaded function.

  1. init(self: dolfin.cpp.mesh.Mesh) -> None

  2. init(self: dolfin.cpp.mesh.Mesh, arg0: int) -> int

  3. init(self: dolfin.cpp.mesh.Mesh, arg0: int, arg1: int) -> None

init_cell_orientations(*args, **kwargs)

Overloaded function.

  1. init_cell_orientations(self: dolfin.cpp.mesh.Mesh, arg0: dolfin.cpp.function.Expression) -> None

  2. init_cell_orientations(self: dolfin.cpp.mesh.Mesh, arg0: object) -> None

init_global(self: dolfin.cpp.mesh.Mesh, arg0: int) None
mpi_comm(self: dolfin.cpp.mesh.Mesh) MPICommWrapper
num_cells(self: dolfin.cpp.mesh.Mesh) int

Number of cells

num_edges(self: dolfin.cpp.mesh.Mesh) int

Number of edges

num_entities(self: dolfin.cpp.mesh.Mesh, arg0: int) int

Number of mesh entities

num_entities_global(self: dolfin.cpp.mesh.Mesh, arg0: int) int
num_faces(self: dolfin.cpp.mesh.Mesh) int

Number of faces

num_facets(self: dolfin.cpp.mesh.Mesh) int

Number of facets

num_vertices(self: dolfin.cpp.mesh.Mesh) int

Number of vertices

order(self: dolfin.cpp.mesh.Mesh) None
ordered(self: dolfin.cpp.mesh.Mesh) bool
rmax(self: dolfin.cpp.mesh.Mesh) float
rmin(self: dolfin.cpp.mesh.Mesh) float
rotate(*args, **kwargs)

Overloaded function.

  1. rotate(self: dolfin.cpp.mesh.Mesh, arg0: float, arg1: int, arg2: dolfin::Point) -> None

  2. rotate(self: dolfin.cpp.mesh.Mesh, angle: float, axis: int = 2) -> None

scale(self: dolfin.cpp.mesh.Mesh, arg0: float) None
smooth(self: dolfin.cpp.mesh.Mesh, num_iterations: int = 1) None
smooth_boundary(self: dolfin.cpp.mesh.Mesh, arg0: int, arg1: bool) None
snap_boundary(self: dolfin.cpp.mesh.Mesh, subdomain: dolfin::SubDomain, harmonic_smoothing: bool = True) None
topology(self: dolfin.cpp.mesh.Mesh) dolfin.cpp.mesh.MeshTopology

Mesh topology

translate(self: dolfin.cpp.mesh.Mesh, arg0: dolfin::Point) None
type(self: dolfin.cpp.mesh.Mesh) dolfin.cpp.mesh.CellType
ufl_cell()
ufl_coordinate_element()

Return the finite element of the coordinate vector field of this domain.

ufl_domain()

Returns the ufl domain corresponding to the mesh.

ufl_id(self: dolfin.cpp.mesh.Mesh) int
class pulse.MeshUnloader(problem, pressure, h5name='test.h5', options=None, h5group='', overwrite=False, merge_control='')

Bases: object

property backward_displacement

Return the current backward displacement as a function on the original geometry.

default_parameters()

Default parameters.

initial_solve(save=False)

Inflate the original geometry to the target pressure and return the displacement field.

save(obj, name, h5group='')

Save object to and HDF file.

Parameters
  • obj (dolfin.Mesh or dolfin.Function) – The object you want to save

  • name (str) – Name of the object

  • h5group (str) – The folder you want to save the object withing the HDF file. Default: ‘’

unload(save=False)

Unload the geometry

property unloaded_geometry
class pulse.Microstructure(f0: Optional[dolfin.function.function.Function] = None, s0: Optional[dolfin.function.function.Function] = None, n0: Optional[dolfin.function.function.Function] = None)

Bases: object

f0 - fibers, s0 - sheets, n0 - cross-sheets

f0: Optional[dolfin.function.function.Function] = None
n0: Optional[dolfin.function.function.Function] = None
s0: Optional[dolfin.function.function.Function] = None
class pulse.MixedParameter(fun, n, name='material_parameters')

Bases: dolfin.function.function.Function

assign_sub(f, i)

Assign subfunction

Parameters
  • f – The function you want to assign

  • i (int) – The subspace number

class pulse.NeoHookean(activation: Optional[ufl.coefficient.Coefficient] = None, parameters=None, active_model: Union[pulse.material.active_model.ActiveModel, pulse.material.active_model.ActiveModels] = ActiveModels.active_strain, T_ref: Optional[float] = None, eta: float = 0.0, isochoric: bool = True, compressible_model='incompressible', geometry=None, f0=None, s0=None, n0=None, active_isotropy: pulse.material.active_model.ActiveStressModels = ActiveStressModels.transversally, *args, **kwargs)

Bases: pulse.material.material_model.Material

Class for Neo Hookean material

W_1(I_1, diff=0, dim=3, *args, **kwargs)
static default_parameters()
name = 'neo_hookean'
strain_energy(F)

Strain-energy density function.

\[\mathcal{W} = \mathcal{W}_1 + \mathcal{W}_{4f} + \mathcal{W}_{\mathrm{active}}\]

where

\[\begin{split}\mathcal{W}_{\mathrm{active}} = \begin{cases} 0 & \text{if acitve strain} \\ \gamma I_{4f} & \text{if active stress} \end{cases}\end{split}\]
Parameters

F (dolfin.Function) – Deformation gradient

class pulse.NeumannBC(traction: Union[float, ufl.coefficient.Coefficient], marker: int, name: str = '')

Bases: object

marker: int
name: str = ''
traction: Union[float, ufl.coefficient.Coefficient]
class pulse.NonlinearProblem(J, F, bcs, output_matrix=False, output_matrix_path='output', **kwargs)

Bases: dolfin.cpp.nls.NonlinearProblem

F(self: dolfin.cpp.nls.NonlinearProblem, arg0: dolfin.cpp.la.GenericVector, arg1: dolfin.cpp.la.GenericVector) None
J(self: dolfin.cpp.nls.NonlinearProblem, arg0: dolfin.cpp.la.GenericMatrix, arg1: dolfin.cpp.la.GenericVector) None
class pulse.NonlinearSolver(problem: pulse.solver.NonlinearProblem, state, parameters=None)

Bases: object

static default_solver_parameters()
solve()

Solve the problem.

Returns

residual – A measure of the accuracy (convergence and error) of the performed computation.

Return type

_solver.snes (???)

update_parameters(parameters)
pulse.PiolaTransform(A, F)
pulse.QuadratureSpace(mesh, degree, dim=3)

From FEniCS version 1.6 to 2016.1 there was a change in how FunctionSpace is defined for quadrature spaces. This functions checks your dolfin version and returns the correct quadrature space

Parameters
  • mesh (dolfin.Mesh) – The mesh

  • degree (int) – The degree of the element

  • dim (int) – For a mesh of topological dimension 3, dim = 1 would be a scalar function, and dim = 3 would be a vector function.

Returns

The quadrature space

Return type

dolfin.FunctionSpace

class pulse.RaghavanUnloader(problem, pressure, h5name='test.h5', options=None, h5group='', overwrite=False, merge_control='')

Bases: pulse.unloader.MeshUnloader

Assumes that the given geometry is loaded with the given pressure.

Find the reference configuration corresponding to zero pressure by first inflating the original geometry to the target pressure and subtacting the displacement times some factor k. Use a 1D minimization algorithm from scipy to find k. The method is described in [1]_ This method assumes that material properties are given. Make sure to change this in utils.py

This also runs in parallel.

Parameters
  • geometry (object) – An obeject with attributes mesh, fiber and markers.

  • pressure (float or tuple) – The target pressure. If LV only provide a float, if BiV proble a tuple of the form (p_LV, p_RV).

  • material_parameters (dict, optional) – A dictionary with material parameters

  • h5name (str, optional) – Path to where you want to save the output Default: test.h5

  • options (dict, optional, optional) –

    A dictionary of solver options. All methods accept the following generic options:

    maxiterint

    Maximum number of iterations to perform.

    tolfloat

    Tolerance of termination

    solve_triesint

    Number of attemtps the solver should use to increase the pressure (after pressure reduction)

  • h5group (str) – Subfolder within the HDF file where you save the results

  • remove_old (bool) – If the file with the same h5name exist, delete it before starting

  • params (dict) – These are parameters that is parsed to pulse_adjoint that make the parameters for the solver. See pulse_adjoint.setup_parametere.setup_adjoint_contraction_parameters

References

1

Raghavan, M. L., Baoshun Ma, and Mark F. Fillinger. “Non-invasive determination of zero-pressure geometry of arterial aneurysms.” Annals of biomedical engineering 34.9 (2006): 1414-1419.

unload_step(u, residual, save=True)
class pulse.RegionalParameter(meshfunction)

Bases: dolfin.function.function.Function

A regional paramerter defined in terms of a dolfin.MeshFunction

Suppose you have a MeshFunction defined different regions in your mesh, and you want to define different parameters on different regions, then this is what you want.

property function

Return linear combination of coefficents and basis functions

Returns

A function with parameter values at each segment specified by the meshfunction

Return type

:py:class`dolfin.Function

get_values()
property proj_space

Space for projecting the scalars. This is a DG 0 space.

pulse.RightCauchyGreen(F, isochoric=False)

Left Cauchy-Green tensor

class pulse.RobinBC(value: Union[float, ufl.coefficient.Coefficient], marker: int)

Bases: object

marker: int
value: Union[float, ufl.coefficient.Coefficient]
pulse.SecondOrderIdentity(F)

Return identity with same dimension as input

class pulse.StVenantKirchhoff(activation: Optional[ufl.coefficient.Coefficient] = None, parameters=None, active_model: Union[pulse.material.active_model.ActiveModel, pulse.material.active_model.ActiveModels] = ActiveModels.active_strain, T_ref: Optional[float] = None, eta: float = 0.0, isochoric: bool = True, compressible_model='incompressible', geometry=None, f0=None, s0=None, n0=None, active_isotropy: pulse.material.active_model.ActiveStressModels = ActiveStressModels.transversally, *args, **kwargs)

Bases: pulse.material.material_model.Material

Class for linear elastic material

static default_parameters()
name = 'saint_venant_kirchhoff'
strain_energy(F_)
pulse.as_backend_type(x)

Return Matrix and Vector backend instance. Not required for other types as pybind11 automatically downcasts objects to the derived type.

pulse.assemble(form, tensor=None, form_compiler_parameters=None, add_values=False, finalize_tensor=True, keep_diagonal=False, backend=None)

Assemble the given form and return the corresponding tensor.

Arguments

Depending on the input form, which may be a functional, linear form, bilinear form or higher rank form, a scalar value, a vector, a matrix or a higher rank tensor is returned.

In the simplest case, no additional arguments are needed. However, additional arguments may and must in some cases be provided as outlined below.

The form can be either a UFL Form or a DOLFIN Form.

If the form defines integrals over different subdomains, MeshFunctions over the corresponding topological entities defining the subdomains can be provided.

The implementation of the returned tensor is determined by the default linear algebra backend. This can be overridden by specifying a different backend.

Each call to assemble() will create a new tensor. If the tensor argument is provided, this will be used instead. Sparsity pattern of tensor is reset iff tensor.empty() holds.

If the keep_diagonal is set to True, assembler ensures that potential zeros on a matrix diagonal are kept in sparsity pattern so every diagonal entry can be changed in a future (for example by ident() or ident_zeros()).

Specific form compiler parameters can be provided by the form_compiler_parameters argument. Form compiler parameters can also be controlled using the global parameters stored in parameters[“form_compiler”].

Examples of usage

The standard stiffness matrix A and a load vector b can be assembled as follows:

A = assemble(inner(grad(u),grad(v))*dx)
b = assemble(f*v*dx)

To prescribe the domains over which integrals will be evaluated, create a Measure with the MeshFunction passed as subdomain_data. For instance, using a mesh function marking parts of the boundary:

# MeshFunction marking boundary parts
boundary_markers = MeshFunction("size_t", mesh, mesh.topology().dim()-1)
# ... fill values in boundary_markers

# Measures with references to cell and boundary markers
ds = Measure("ds", subdomain_data=boundary_markers)

# Sample variational forms
a = inner(grad(u), grad(v))*dx + p*u*v*ds(0)
L = f*v*dx - g*v*ds(1) + p*q*v*ds(0)

A = assemble(a)
b = assemble(L)

For interior facet integrals (dS), cell markers can be used to specify which cell is ‘+’ and which is ‘-‘. The ‘+’ and ‘-‘ sides are chosen such that the cell marker value in the cell at the ‘+’ side cell is larger than the cell marker value in the cell at the ‘-‘ side cell. If the values are equal or the cell markers are not provided, the sides are chosen arbitrarily.

Note that currently, cell markers must be associated with a cell type integral (dx), and if you don’t have such an integral a workaround is to add the integral of something over an empty domain such as ‘f*dx(99)’ with 99 a number not occuring among the cell markers. A better interface to handle this case will be provided later.

# MeshFunctions marking boundary and cell parts
boundary_markers = MeshFunction("size_t", mesh, mesh.topology().dim()-1)
cell_markers = MeshFunction("size_t", mesh, mesh.topology().dim())
# ... fill values in boundary_markers

# Measures with references to cell and boundary markers
ds = Measure("ds", domain=mesh, subdomain_data=boundary_markers)
dx = Measure("dx", domain=mesh, subdomain_data=cell_markers)

# Sample variational forms
a = inner(grad(u), grad(v))*dx + p*u*v*ds(0)
L = v*dx(99) - g*v*ds(1) + p*q*v*ds(0)

A = assemble(a)
b = assemble(L)

To ensure that the assembled matrix has the right type, one may use the tensor argument:

A = PETScMatrix()
assemble(a, tensor=A)

The form a is now assembled into the PETScMatrix A.

pulse.interpolate(v, V)

Return interpolation of a given function into a given finite element space.

Arguments
v

a Function or a :py:class:`MultiMeshFunction

<dolfin.function.function.MultiMeshFunction>` or

an Expression

V

a FunctionSpace (standard, mixed, etc.) or a MultiMeshFunctionSpace.

Example of usage

v = Expression("sin(pi*x[0])")
V = FunctionSpace(mesh, "Lagrange", 1)
Iv = interpolate(v, V)
pulse.project(v, V=None, bcs=None, mesh=None, function=None, solver_type='lu', preconditioner_type='default', form_compiler_parameters=None)

Return projection of given expression v onto the finite element space V.

Arguments
v

a Function or an Expression

bcs

Optional argument list of DirichletBC

V

Optional argument FunctionSpace

mesh

Optional argument mesh.

solver_type

see solve for options.

preconditioner_type

see solve for options.

form_compiler_parameters

see Parameters for more information.

Example of usage

v = Expression("sin(pi*x[0])")
V = FunctionSpace(mesh, "Lagrange", 1)
Pv = project(v, V)

This is useful for post-processing functions or expressions which are not readily handled by visualization tools (such as for example discontinuous functions).

pulse.set_log_level(level)