Choral
Data Types | Functions/Subroutines
cardiomodel_mod Module Reference

DERIVED TYPE cardioModel More...

Data Types

interface  cardiomodel
 Derived type for tissue properties in cardiac electrophysiology. More...
 
interface  print
 print a short description More...
 

Functions/Subroutines

type(cardiomodel) function cardiomodel_create (fibre, Am, conduc_type)
  constructor for cardioModel type More...
 
subroutine set_conductivities (gil, git, gel, get, type)
 set conductivities More...
 
subroutine, public scale_conductivities (cm, sc)
 rescale the conductivities More...
 
subroutine cardiomodel_print (cm)
 print cardio model More...
 
subroutine, public monodomain_assemble (M, S, cm, X_h, qd_M, qd_S)
 Assemble the monodomain diffusion matrices. More...
 
real(rp) function usigmav (u, v, fd, gl, gt)
 scalar product u^T Sigma v More...
 

Detailed Description

DERIVED TYPE cardioModel

Tissue model in cardiac electrophysiology

Author
Charles Pierre

Function/Subroutine Documentation

◆ cardiomodel_create()

type(cardiomodel) function cardiomodel_mod::cardiomodel_create ( procedure(r3tor3)  fibre,
real(rp), intent(in)  Am,
integer, intent(in)  conduc_type 
)
private

constructor for cardioModel type

  • OUTPUT:
    • cm = cardiac tissue model
  • INPUT:
    • fbre = fibre orientation (procedural)
    • Am = cell membrane surface to volume ratio
    • conduc_type = conductivity denomination

Definition at line 99 of file cardioModel_mod.f90.

Here is the call graph for this function:

◆ cardiomodel_print()

subroutine cardiomodel_mod::cardiomodel_print ( type(cardiomodel), intent(in)  cm)
private

print cardio model

Definition at line 173 of file cardioModel_mod.f90.

◆ monodomain_assemble()

subroutine, public cardiomodel_mod::monodomain_assemble ( type(csr), intent(out)  M,
type(csr), intent(out)  S,
class(cardiomodel), intent(in)  cm,
type(fespace), intent(in)  X_h,
integer, intent(in)  qd_M,
integer, intent(in)  qd_S 
)

Assemble the monodomain diffusion matrices.

OUTPUT :

  • M = mass matrix for the monodomain model
  • S = stiffness matrix for the monodomain model

INPUT :

  • cm = cardio model
  • X = finite element mesh node coordinates
  • fe = finite element method type
  • qd_M, qd_S = quadrature methods to assemble M and S

Definition at line 198 of file cardioModel_mod.f90.

Here is the call graph for this function:

◆ scale_conductivities()

subroutine, public cardiomodel_mod::scale_conductivities ( type(cardiomodel), intent(inout)  cm,
real(rp), intent(in)  sc 
)

rescale the conductivities

Definition at line 159 of file cardioModel_mod.f90.

◆ set_conductivities()

subroutine cardiomodel_mod::set_conductivities ( real(rp), intent(out)  gil,
real(rp), intent(out)  git,
real(rp), intent(out)  gel,
real(rp), intent(out)  get,
integer, intent(in)  type 
)
private

set conductivities

Definition at line 121 of file cardioModel_mod.f90.

Here is the call graph for this function:

◆ usigmav()

real(rp) function cardiomodel_mod::usigmav ( real(rp), dimension(3), intent(in)  u,
real(rp), dimension(3), intent(in)  v,
real(rp), dimension(3), intent(in)  fd,
real(rp), intent(in)  gl,
real(rp), intent(in)  gt 
)
private

scalar product u^T Sigma v

with sigma = Diag( gl, gt, gt) in an orthonormal bass (e1, e2, e3) where e1 is colinear to fd

fd = fibre direction gl, gt = transverse and longitudinal conductivities

Definition at line 256 of file cardioModel_mod.f90.

Here is the call graph for this function: