API documentation¶
pylsewave.mesh¶
Pylsewave.mesh is a module that contains classes for vessel/tube and vessel network definition.
-
class
pylsewave.mesh.
Vessel
(name, L, R_proximal, R_distal, Wall_thickness, WindKessel=None, Id=0, model_type='Linear_elastic')¶ Bases:
object
Base class for singe-vessel/tube definition
Vessel-Class constructor
Parameters: - name (str) -- name of the vessel/segment
- L (float) -- Length of the segment
- R_proximal (float) -- Proximal radius
- R_distal (float) -- Distal radius
- Wall_thickness (float) -- (float) Wall thickness
- WindKessel (dict) -- Windkessel params (R, L, C)
- Id (int) -- Unique id of the segment
- model_type (str) -- model of the vessel (linear-elastic, visco-elastic)
- priority -- "Linear Elastic" or "Visco-elastic"
Example: >>> myVessel = Vessel(name="Barchial artery", L=100, R_proximal=4.0, R_disWindKessel=None, Id=0, model_type='Linear_elastic')
-
RLC
¶ Getter: returns the RLC params Setter: sets the RLC params Type: dict Raises: ValueError -- if the input is not a dictionary
-
calculate_R0
(x)¶ Method to calculate the reference radius \(R_0(x)\). The default model is
\[R_0(x) = R_{prox} \exp(\log(R_{distal} / R_{prox})(x/L))\]Parameters: x -- (array or float) spatial point to calculate the reference diameter Returns: (array or float) \(R_0(x)\)
-
static
dfdr
(R0, k)¶ This property calculates the df(r0, k)/dr function
Parameters: - R0 -- Reference diameter
- k (ndarray[ndim=1, type=float]) -- the empirical k params
Returns: float or ndarray[ndim=1, type=float]
-
dx
¶ The spatial descretisation of the Vessel object.
When we set this property/attribute of the object, several other Vessel characteristics are calculated, such as
mesh.Vessel.r0()
,mesh.Vessel.f_r0()
, etc.Getter: returns the spatial dx of the vessel Setter: sets the dx of vessel/segment Type: (float)
-
static
f
(R0, k)¶ This property calculates the f(r0, k) function
Parameters: - R0 -- Reference diameter
- k (ndarray[ndim=1, type=float]) -- the empirical k params
Returns: float or ndarray[ndim=1, type=float]
-
f_r0
¶ This property contains the f(r0, k) function
Getter: returns the f value calculated in every node of the 1D vessel
-
id
¶ Getter: returns the unique vessel/segment Id Setter: sets the unique vessel/segment Id Type: (int)
-
interpolate_R0
(value)¶ Method to interpolate initial radius values in points offset from nodes
Parameters: value (float) -- the node offset interpolation points (e.g. 0.5, -0.5) Returns: float or ndarray[ndim=1, type=float]
-
length
¶ Getter: returns the Length of the 1D segment Type: float
-
name
¶ Getter: returns the name of the vessel Type: (string)
-
r0
¶ Getter: returns the reference radius along the 1D segment Type: ndarray[ndim=1, type=float]
-
r_dist
¶ Getter: returns the: \(R_{distal}\) Type: float
-
r_prox
¶ Getter: returns the: \(R_{proximal}\) Type: float
-
set_k_vector
(k)¶ Sets the k vector (k1, k2, k3) params :param k: the k vector :type k: ndarray[ndim=1, type=float] :return: None :raise ValueError: if the k input variable is not ndarray
-
w_th
¶ Getter: returns the \(W_{thickness}\) Type: float
-
class
pylsewave.mesh.
VesselNetwork
(vessels, dx=None, Nx=None)¶ Bases:
object
Class to assemble the single vessels/segments that comprise the Vessel network
Parameters: - vessels -- (list) object containing the Vessel segments
- dx -- (float) global spatial descretisation
- Nx -- (int) global number of nodes per element
Raises: ValueError -- if the vessels are not pylsewave.mesh.Vessel type
-
Nx
¶ Getter: returns the global Nx (number of nodes) Type: int
-
bif
¶ Getter: returns the vessel connectivity (bifurcations) Setter: sets the vessel connectivity (bifurcations) Type: ndarray[ndim=2, type=float]
-
conj
¶ Getter: returns the vessel connectivity (conjunctions) Setter: sets the vessel connectivity (conjunctions) Type: ndarray[ndim=2, type=float]
-
dx
¶ Getter: returns the global dx Type: float
pylsewave.pdes¶
-
class
pylsewave.pdes.
PDEm
(mesh, rho, mu, Re=None)¶ Bases:
object
Base class for PDE system declaration.
Parameters: - mesh (pylsewave.mesh.VesselNetwork) -- the vessel network
- rho (float) -- the density of fluid medium (e.g. blood, water, etc.)
- Re (float) -- the Reynolds number
Raises: TypeError -- if the mesh is not pylsewave.mesh.VesselNetwork
-
compute_c
(A, index=None, vessel_index=None)¶ Method to compute local wave speed
Parameters: - A -- cross-sectional area of vessel
- index (int) -- the number or the node to be calculated
- vessel_index -- the vessel unique Id
Returns: float or ndarray[ndim=1, type=float]
-
flux
(u, x, index=None, vessel_index=None)¶ Method to calculate the Flux term.
Parameters: - u -- the solution vector u
- x -- the spatial points discretising the 1D segments
- index (int) -- the index Id of the node
- vessel_index (int) -- the unique vessel Id
Returns: ndarray[ndim=1, dtype=float] or ndarray[ndim=2, dtype=float]
Raises: NotImplementedError -- This is a base class, define an inherited class to override this method
-
pressure
(a, index=None, vessel_index=None)¶ Method to compute pressure from pressure-area relationship
Parameters: - a -- cross-sectional area of the vessel
- index (int) -- index of the node in vessel
- vessel_index (int) -- the vessel unique Id
Returns: float or ndarray[ndim=1, type=float]
-
set_boundary_layer_th
(T, no_cycles)¶ method to calculate the boundary layer
\[\delta = \sqrt{\frac{\nu T_{cycle}}{2 \pi}}\]Parameters: - T -- (float) Total period of the simulation
- no_cycles -- No of cycles that the simulation will run
Returns: (float) Boundary layer delta
-
source
(u, x, index=None, indexVessel=None)¶ Method to calculate the Source term.
Parameters: - u -- the solution vector u
- x -- the spatial points discretising the 1D segments
- index (int) -- the index Id of the node
- vessel_index (int) -- the unique vessel Id
Returns: ndarray[ndim=1, dtype=float] or ndarray[ndim=2, dtype=float]
Raises: NotImplementedError -- This is a base class, define an inherited class to override this method
-
class
pylsewave.pdes.
PDEsVisco
(mesh, rho, mu, Re=None)¶ Bases:
pylsewave.pdes.PDEsWat
PDEs class inherited from PDEsWat. This PDE system contains the viscous part of the model, as well.
Parameters: - mesh (pylsewave.mesh.VesselNetwork) -- the vessel network
- rho (float) -- the density of fluid medium (e.g. blood, water, etc.)
- Re (float) -- the Reynolds number
Raises: TypeError -- if the mesh is not pylsewave.mesh.VesselNetwork
-
compute_c
(A, index=None, vessel_index=None)¶ Method to compute local wave speed
\[c(A) = \sqrt{\frac{1}{2 \rho} f(R_0, k) \sqrt{\frac{A}{A_0}}}\]Parameters: - A -- cross-sectional area of vessel
- index (int) -- the number or the node to be calculated
- vessel_index -- the vessel unique Id
Returns: float or ndarray[ndim=1, type=float]
-
flux
(u, x, index=None, vessel_index=None)¶ Method to calculate the Flux term.
Parameters: - u -- the solution vector u
- x -- the spatial points discretising the 1D segments
- index (int) -- the index Id of the node
- vessel_index (int) -- the unique vessel Id
Returns: ndarray[ndim=1, dtype=float] or ndarray[ndim=2, dtype=float]
-
pressure
(a, q, index=None, vessel_index=None)¶ Method to compute pressure from pressure-area relationship
\[p(A) = f(R_0, k) \left( \sqrt{\frac{A}{A_0}} - 1 \right)\]Parameters: - a -- cross-sectional area of the vessel
- index (int) -- index of the node in vessel
- vessel_index (int) -- the vessel unique Id
Returns: float or ndarray[ndim=1, type=float]
-
set_boundary_layer_th
(T, no_cycles)¶ method to calculate the boundary layer
\[\delta = \sqrt{\frac{\nu T_{cycle}}{2 \pi}}\]Parameters: - T -- (float) Total period of the simulation
- no_cycles -- No of cycles that the simulation will run
Returns: (float) Boundary layer delta
-
source
(u, x, index=None, indexVessel=None)¶ Method to calculate the Source term.
Parameters: - u -- the solution vector u
- x -- the spatial points discretising the 1D segments
- index (int) -- the index Id of the node
- vessel_index (int) -- the unique vessel Id
Returns: ndarray[ndim=1, dtype=float] or ndarray[ndim=2, dtype=float]
-
class
pylsewave.pdes.
PDEsWat
(mesh, rho, mu, Re=None)¶ Bases:
pylsewave.pdes.PDEm
Class for PDE system as defined in Watanabe et al. 2013, Mathematical model of blood flow in an anatomically detailed arterial network of the arm, In: Mathematical modelling and numerical analysis.
Parameters: - mesh (pylsewave.mesh.VesselNetwork) -- the vessel network
- rho (float) -- the density of fluid medium (e.g. blood, water, etc.)
- Re (float) -- the Reynolds number
Raises: TypeError -- if the mesh is not pylsewave.mesh.VesselNetwork
-
compute_c
(A, index=None, vessel_index=None)¶ Method to compute local wave speed
\[c(A) = \sqrt{\frac{1}{2 \rho} f(R_0, k) \sqrt{\frac{A}{A_0}}}\]Parameters: - A -- cross-sectional area of vessel
- index (int) -- the number or the node to be calculated
- vessel_index -- the vessel unique Id
Returns: float or ndarray[ndim=1, type=float]
-
flux
(u, x, index=None, vessel_index=None)¶ Method to calculate the Flux term.
Parameters: - u -- the solution vector u
- x -- the spatial points discretising the 1D segments
- index (int) -- the index Id of the node
- vessel_index (int) -- the unique vessel Id
Returns: ndarray[ndim=1, dtype=float] or ndarray[ndim=2, dtype=float]
-
pressure
(a, index=None, vessel_index=None)¶ Method to compute pressure from pressure-area relationship
\[p(A) = f(R_0, k) \left( \sqrt{\frac{A}{A_0}} - 1 \right)\]Parameters: - a -- cross-sectional area of the vessel
- index (int) -- index of the node in vessel
- vessel_index (int) -- the vessel unique Id
Returns: float or ndarray[ndim=1, type=float]
-
set_boundary_layer_th
(T, no_cycles)¶ method to calculate the boundary layer
\[\delta = \sqrt{\frac{\nu T_{cycle}}{2 \pi}}\]Parameters: - T -- (float) Total period of the simulation
- no_cycles -- No of cycles that the simulation will run
Returns: (float) Boundary layer delta
-
source
(u, x, index=None, indexVessel=None)¶ Method to calculate the Source term.
Parameters: - u -- the solution vector u
- x -- the spatial points discretising the 1D segments
- index (int) -- the index Id of the node
- vessel_index (int) -- the unique vessel Id
Returns: ndarray[ndim=1, dtype=float] or ndarray[ndim=2, dtype=float]
pylsewave.bcs¶
-
class
pylsewave.bcs.
BCs
(inpdes, inletfun=None)¶ Bases:
object
Base class for BCs definition.
Parameters: - inpdes (pylsewave.pdes.PDEm) -- the PDE class of the problem
- inletfun (pylsewave.interpolate.CubicSpline.eval_spline) -- interpolation fun to calculate the inlet prescribed field
Raises: TypeError -- if the inpdes is not pylsewave.pdes.PDEm type
-
I
(x, vessel_index=None)¶ Sets the initial condition of the state variables
Parameters: - x -- the spatial coordinate of the node
- vessel_index (int) -- the unique Id of the vessel
Returns: ndarray[ndim=1, type=float]
-
U_0
(u, t, dx, dt, vessel_index, out)¶ Method to compute the inlet BCs.
Example: Let assume that we have a MacCormack scheme. Moreover, we solve the PDE system with respect to A, Q. We prescribe also flow waveform at the inlet \(Q_{x=0}(t)\). Then from continuity equation, the cross sectional area \(A_{x=0}^{t=n+1}\) can be calculated as
\[A_{x=0}^{t=n+1} = A_{x=0}^{t=n} - \left(2 \frac{dt}{dx} \left( Q_{x=1}^{t=n+1} - Q_{x=0}^{t=n+1}\right) \right)\]Parameters: - u -- the solution vector u_n
- t -- the time step of the simulation
- dx -- the spatial discretisation of the vessel
- dt -- the space discretisation of the simulation
- vessel_index -- the unique vessel index
- out -- the solution array u
Returns: None
-
U_L
(u, t, dx, dt, vessel_index, out)¶ Class-method to compute the outlet BCs. This method implements the three-element Windkessel model (RCR) along with the iterative method proposed in Kolachalama et al 2007, Predictive Haemodynamics in a One-Dimensional Carotid Artery Bifurcation. Part I Application to Stent Design. In: IEEE Transactions on Biomedical Engineering. This method has also been implemented in VamPy
Parameters: - u -- Solution vector u_n
- t -- time step
- dx -- vessel/segment dx
- dt -- time discretisation dt
- vessel_index -- unique vessel index
- out -- solution vector u
Returns: None
-
bifurcation_Jr
(x, u, dt, vessel_index_list)¶ Method to calculate the Jacobian of the non-linear system formed at the bifurcations. Normally the output will be a 6x6 array/matrix.
Parameters: - x -- the solution vector at iteration i
- u -- the solution vector u_n of internal nodes
- dt -- time step dt
- vessel_index_list -- the unique vessel Id
Returns: ndarray[ndim=2, dtype=float]
Raises: NotImpementedError -- this is abstract class, define an inherited class to override this method.
-
bifurcation_R
(x, u, dt, vessel_index_list)¶ Method to calculate the non-linear system of equations in the bifurcation. Usually, 6 equations are given.
Parameters: - x -- the solution vector x at an iteration i
- u -- the solution vector u_n of internal nodes
- dt -- time step dt
- vessel_index_list -- unique vessel Id
Returns: ndarray[ndim=1, type=float]
Raises: NotImpementedError -- this is abstract class, define an inherited class to override this method.
-
conjuction_Jr
(x, u, dt, vessel_index_list)¶ Method to calculate the Jacobian of the non-linear system formed at the conjunctions. Normally the output will be a 4x4 array/matrix.
Parameters: - x -- the solution vector at iteration i
- u -- the solution vector u_n of internal nodes
- dt -- time step dt
- vessel_index_list -- the unique vessel Id
Returns: ndarray[ndim=2, dtype=float]
Raises: NotImpementedError -- this is abstract class, define an inherited class to override this method.
-
conjuction_R
(x, u, dt, vessel_index_list)¶ Method to calculate the non-linear system of equations in a conjunction. Usually, 4 equations are given.
Parameters: - x -- the solution vector x at an iteration i
- u -- the solution vector u_n of internal nodes
- dt -- time step dt
- vessel_index_list -- unique vessel Id
Returns: ndarray[ndim=1, type=float]
Raises: NotImpementedError -- this is abstract class, define an inherited class to override this method.
-
pdes
¶ Getter: returns the pde system Type: pylsewave.pdes.PDEm
-
class
pylsewave.bcs.
BCsWat
(inpdes, inletfun=None)¶ Bases:
pylsewave.bcs.BCs
Class for BCs definition. pylsewave.bcs.BCs inherited class for the pylsewave.pdes.PDEsWat class.
Parameters: - inpdes (pylsewave.pdes.PDEm) -- the PDE class of the problem
- inletfun (pylsewave.interpolate.CubicSpline.eval_spline) -- interpolation fun to calculate the inlet prescribed field
Raises: TypeError -- if the inpdes is not pylsewave.pdes.PDEm type
-
I
(x, vessel_index=None)¶ Sets the initial condition of the state variables
Parameters: - x -- the spatial coordinate of the node
- vessel_index (int) -- the unique Id of the vessel
Returns: ndarray[ndim=1, type=float]
-
U_0
(u, t, dx, dt, vessel_index, out)¶ Method to compute the inlet BCs. This method assumes that a pressure is prescribed at the inlet.
Note
override the method in case you want to prescribe other type of inlet BCs.
Parameters: - u -- the solution vector u_n
- t -- the time step of the simulation
- dx -- the spatial discretisation of the vessel
- dt -- the space discretisation of the simulation
- vessel_index -- the unique vessel index
- out -- the solution array u
Returns: None
-
U_L
(u, t, dx, dt, vessel_index, out)¶ Class-method to compute the outlet BCs. This method implements the three-element Windkessel model (RCR) along with the iterative method proposed in Kolachalama et al 2007, Predictive Haemodynamics in a One-Dimensional Carotid Artery Bifurcation. Part I Application to Stent Design. In: IEEE Transactions on Biomedical Engineering. This method has also been implemented in VamPy
Parameters: - u -- Solution vector u_n
- t -- time step
- dx -- vessel/segment dx
- dt -- time discretisation dt
- vessel_index -- unique vessel index
- out -- solution vector u
Returns: None
-
bifurcation_Jr
(x, u, dt, vessel_index_list)¶ Method to calculate the Jacobian of the non-linear system formed at the bifurcations. Normally the output will be a 6x6 array/matrix.
Parameters: - x -- the solution vector at iteration i
- u -- the solution vector u_n of internal nodes
- dt -- time step dt
- vessel_index_list -- the unique vessel Id
Returns: ndarray[ndim=2, dtype=float]
-
bifurcation_R
(x, u, dt, vessel_index_list)¶ Method to calculate the non-linear system of equations in the bifurcation. Usually, 6 equations are given.
Parameters: - x -- the solution vector x at an iteration i
- u -- the solution vector u_n of internal nodes
- dt -- time step dt
- vessel_index_list -- unique vessel Id
Returns: ndarray[ndim=1, type=float]
-
conjuction_Jr
(x, u, dt, vessel_index_list)¶ Method to calculate the Jacobian of the non-linear system formed at the conjunctions. Normally the output will be a 4x4 array/matrix.
Parameters: - x -- the solution vector at iteration i
- u -- the solution vector u_n of internal nodes
- dt -- time step dt
- vessel_index_list -- the unique vessel Id
Returns: ndarray[ndim=2, dtype=float]
-
conjuction_R
(x, u, dt, vessel_index_list)¶ Method to calculate the non-linear system of equations in a conjunction. Usually, 4 equations are given.
Parameters: - x -- the solution vector x at an iteration i
- u -- the solution vector u_n of internal nodes
- dt -- time step dt
- vessel_index_list -- unique vessel Id
Returns: ndarray[ndim=1, type=float]
-
static
linear_extrapolation
(x, x1, x2, y1, y2)¶ Static method for linear extrapolation..
Parameters: - x -- location of extrapolation
- x1 -- x value left
- x2 -- y value left
- y1 -- x value right
- y2 -- y value right
Returns: extrapolation value
-
pdes
¶ Getter: returns the pde system Type: pylsewave.pdes.PDEm
pylsewave.fdm¶
pylsewave module containing classes for finite difference numerical computation. E.g. MacCormack, Lax-Wedroff, etc.
Parts of this module have been recycled from Finite difference computing with partial differential equations
-
class
pylsewave.fdm.
BloodWaveLaxWendroff
(inbcs)¶ Bases:
pylsewave.fdm.FDMSolver
Class with 2 step Lax-Wendroff scheme.
\[U_i^{n+1} = U_i^{n} - \frac{\Delta t}{\Delta x} \left( F_{i + \frac{1}{2}}^{n + \frac{1}{2}} - F_{i - \frac{1}{2}}^{n + \frac{1}{2}} \right) + \frac{\Delta t}{2}\left( S_{i + \frac{1}{2}}^{n + \frac{1}{2}} + S_{i - \frac{1}{2}}^{n + \frac{1}{2}} \right)\]where the intermediate values calculated as
\[U_j^{n+\frac{1}{2}} = \frac{U_{j+1/2}^n + U_{j-1/2}^n}{2} - \frac{\Delta t}{2 \Delta x}\left( F_{j+1/2} - F_{j-1/2} \right) + \frac{\Delta t}{4}\left( S_{j+1/2} + S_{j-1/2} \right)\]Parameters: inbcs (pylsewave.pdes.PDEm) -- class with the PDE definition Raises: TypeError -- if the in pdes is not type pylsewave.pdes.PDEm -
cfl_condition
(u, dx, dt, vessel_index)¶ Method to check the pre-defined CFL condition.
Parameters: - u (numpy.ndarray) -- array with the solution u
- dx -- spatial discretisation of the vessel
- dt -- time step
- vessel_index -- unique vessel Id
Returns: bool
-
dt
¶ Getter: returns the dt of the solver Type: float
-
set_BC
(U0_array, UL_array, UBif_array=None, UConj_array=None)¶ Method to set the BCs/connectivity of the network/mesh
Parameters: - U0_array (numpy.ndarray) -- array with vessel Ids that inlet bcs will be prescribed
- UL_array (numpy.ndarray) -- array with vessel Ids that outlet bcs will be prescribed
- UBif_array (numpy.ndarray) -- array with vessel Ids forming bifurcations
- UConj_array (numpy.ndarray) -- array with vessel Ids forming conjunctions
Returns: None
-
set_T
(dt, T, no_cycles)¶ Method to set the total period along with the time stepping of the solver. :param float dt: time step :param float T: total period :param no_cycles: number of the cycles (if periodic) :return: None
-
solve
(casename, user_action=None, version='vectorized')¶ Method to solve the hyperbolic part of the PDE system:
\[U_t + F(U)_x = S(U)\]Parameters: - casename (str) -- the name of the simulation case
- user_action -- a callback function/class to store the solution
- version (str) -- scalar or vectorised computations
Returns: int
-
-
class
pylsewave.fdm.
BloodWaveMacCormack
(inbcs)¶ Bases:
pylsewave.fdm.FDMSolver
Class with the MacCormack scheme implementation.
\[U_i^{\star} = U_i^n - \frac{\Delta t}{\Delta x}\left( F_{i+1}^n - F_i^n \right) + \Delta t S_i^n\]and
\[U_i^{n+1} = \frac{1}{2}\left( U_i^n + U_i^{\star} \right) - \frac{\Delta t}{2 \Delta x}\left( F_i^{\star} - F_{i-1}^{\star} \right) + \frac{\Delta t}{2}S_i^{\star}\]Parameters: inbcs (pylsewave.pdes.PDEm) -- class with the PDE definition Raises: TypeError -- if the in pdes is not type pylsewave.pdes.PDEm -
cfl_condition
(u, dx, dt, vessel_index)¶ Method to check the pre-defined CFL condition.
Parameters: - u (numpy.ndarray) -- array with the solution u
- dx -- spatial discretisation of the vessel
- dt -- time step
- vessel_index -- unique vessel Id
Returns: bool
-
dt
¶ Getter: returns the dt of the solver Type: float
-
set_BC
(U0_array, UL_array, UBif_array=None, UConj_array=None)¶ Method to set the BCs/connectivity of the network/mesh
Parameters: - U0_array (numpy.ndarray) -- array with vessel Ids that inlet bcs will be prescribed
- UL_array (numpy.ndarray) -- array with vessel Ids that outlet bcs will be prescribed
- UBif_array (numpy.ndarray) -- array with vessel Ids forming bifurcations
- UConj_array (numpy.ndarray) -- array with vessel Ids forming conjunctions
Returns: None
-
set_T
(dt, T, no_cycles)¶ Method to set the total period along with the time stepping of the solver. :param float dt: time step :param float T: total period :param no_cycles: number of the cycles (if periodic) :return: None
-
solve
(casename, user_action=None, version='vectorized')¶ Method to solve the hyperbolic part of the PDE system:
\[U_t + F(U)_x = S(U)\]Parameters: - casename (str) -- the name of the simulation case
- user_action -- a callback function/class to store the solution
- version (str) -- scalar or vectorised computations
Returns: int
-
-
class
pylsewave.fdm.
BloodWaveMacCormackGodunov
(inbcs)¶ Bases:
pylsewave.fdm.BloodWaveMacCormack
Parameters: inbcs (pylsewave.pdes.PDEm) -- class with the PDE definition Raises: TypeError -- if the in pdes is not type pylsewave.pdes.PDEm -
cfl_condition
(u, dx, dt, vessel_index)¶ Method to check the pre-defined CFL condition.
Parameters: - u (numpy.ndarray) -- array with the solution u
- dx -- spatial discretisation of the vessel
- dt -- time step
- vessel_index -- unique vessel Id
Returns: bool
-
dt
¶ Getter: returns the dt of the solver Type: float
-
set_BC
(U0_array, UL_array, UBif_array=None, UConj_array=None)¶ Method to set the BCs/connectivity of the network/mesh
Parameters: - U0_array (numpy.ndarray) -- array with vessel Ids that inlet bcs will be prescribed
- UL_array (numpy.ndarray) -- array with vessel Ids that outlet bcs will be prescribed
- UBif_array (numpy.ndarray) -- array with vessel Ids forming bifurcations
- UConj_array (numpy.ndarray) -- array with vessel Ids forming conjunctions
Returns: None
-
set_T
(dt, T, no_cycles)¶ Method to set the total period along with the time stepping of the solver. :param float dt: time step :param float T: total period :param no_cycles: number of the cycles (if periodic) :return: None
-
solve
(casename, user_action=None, version='vectorized')¶ Solve U_t + F_x = S
-
-
class
pylsewave.fdm.
FDMSolver
(inbcs)¶ Bases:
object
Base class for finite-difference computing schemes.
Parameters: inbcs (pylsewave.pdes.PDEm) -- class with the PDE definition Raises: TypeError -- if the in pdes is not type pylsewave.pdes.PDEm -
cfl_condition
(u, dx, dt, vessel_index)¶ Method to check the pre-defined CFL condition.
Parameters: - u (numpy.ndarray) -- array with the solution u
- dx -- spatial discretisation of the vessel
- dt -- time step
- vessel_index -- unique vessel Id
Returns: bool
-
dt
¶ Getter: returns the dt of the solver Type: float
-
set_BC
(U0_array, UL_array, UBif_array=None, UConj_array=None)¶ Method to set the BCs/connectivity of the network/mesh
Parameters: - U0_array (numpy.ndarray) -- array with vessel Ids that inlet bcs will be prescribed
- UL_array (numpy.ndarray) -- array with vessel Ids that outlet bcs will be prescribed
- UBif_array (numpy.ndarray) -- array with vessel Ids forming bifurcations
- UConj_array (numpy.ndarray) -- array with vessel Ids forming conjunctions
Returns: None
-
set_T
(dt, T, no_cycles)¶ Method to set the total period along with the time stepping of the solver. :param float dt: time step :param float T: total period :param no_cycles: number of the cycles (if periodic) :return: None
-
solve
(casename, user_action=None, version='vectorized')¶ Method to solve the hyperbolic part of the PDE system:
\[U_t + F(U)_x = S(U)\]Parameters: - casename (str) -- the name of the simulation case
- user_action -- a callback function/class to store the solution
- version (str) -- scalar or vectorised computations
Returns: int
-
pylsewave.cynum¶
This is a Cython file with classes for pylsewave toolkit. Cython translates everything to C++. This module contains solvers that support OPENMP (via Cython OpenMP) for parallel CPU calcs.
-
class
pylsewave.cynum.
PDEmCs
¶ Bases:
pylsewave.cynum.cPDEm
-
set_boundary_layer_th
()¶ Method to calculate the boundary layer.
\[\delta = \sqrt{\frac{\nu T_{cycle}}{2 \pi}}\]Parameters: - T -- (float) Total period of the simulation
- no_cycles -- No of cycles that the simulation will run
Returns: (float) Boundary layer delta
-
-
class
pylsewave.cynum.
cPDEm
¶ Bases:
object
Parameters: - mesh (pylsewave.mesh.VesselNetwork) -- the vessel network
- rho (float) -- the density of fluid medium (e.g. blood, water, etc.)
- Re (float) -- the Reynolds number
Param: float mu: blood viscosity
-
set_boundary_layer_th
()¶ Method to calculate the boundary layer.
\[\delta = \sqrt{\frac{\nu T_{cycle}}{2 \pi}}\]Parameters: - T -- (float) Total period of the simulation
- no_cycles -- No of cycles that the simulation will run
Returns: (float) Boundary layer delta
-
class
pylsewave.cynum.
cPDEsWat
¶ Bases:
pylsewave.cynum.cPDEm
-
set_boundary_layer_th
()¶ Method to calculate the boundary layer.
\[\delta = \sqrt{\frac{\nu T_{cycle}}{2 \pi}}\]Parameters: - T -- (float) Total period of the simulation
- no_cycles -- No of cycles that the simulation will run
Returns: (float) Boundary layer delta
-
-
class
pylsewave.cynum.
cPDEsWatVisco
¶ Bases:
pylsewave.cynum.cPDEsWat
-
set_boundary_layer_th
()¶ Method to calculate the boundary layer.
\[\delta = \sqrt{\frac{\nu T_{cycle}}{2 \pi}}\]Parameters: - T -- (float) Total period of the simulation
- no_cycles -- No of cycles that the simulation will run
Returns: (float) Boundary layer delta
-