openCARP
Doxygen code documentation for the open cardiac electrophysiology simulator openCARP
electric_integrators.cc
Go to the documentation of this file.
1 // ----------------------------------------------------------------------------
2 // openCARP is an open cardiac electrophysiology simulator.
3 //
4 // Copyright (C) 2020 openCARP project
5 //
6 // This program is licensed under the openCARP Academic Public License (APL)
7 // v1.0: You can use and redistribute it and/or modify it in non-commercial
8 // academic environments under the terms of APL as published by the openCARP
9 // project v1.0, or (at your option) any later version. Commercial use requires
10 // a commercial license (info@opencarp.org).
11 //
12 // This program is distributed without any warranty; see the openCARP APL for
13 // more details.
14 //
15 // You should have received a copy of the openCARP APL along with this program
16 // and can find it online: http://www.opencarp.org/license
17 // ----------------------------------------------------------------------------
18 
27 #include "electric_integrators.h"
28 #include "sf_interface.h"
29 
30 namespace opencarp {
31 
32 #define ELEC_STIFFNESS 1
33 #define ELEC_MASS 2
34 
35 inline int get_preferred_int_order(SF::elem_t & etype, int mat_type)
36 {
37  switch(mat_type) {
38  default:
39  case ELEC_MASS: return 2;
40 
41  case ELEC_STIFFNESS: {
42  switch(etype) {
43  default:
44  case SF::Quad:
45  case SF::Hexa:
46  case SF::Prism:
47  case SF::Pyramid:
48  return 2;
49 
50  case SF::Line:
51  case SF::Tetra:
52  case SF::Tri:
53  return 1;
54  }
55  }
56  }
57 }
58 
60  const int* reg,
61  const double vol,
62  const double* eVals,
63  const char* msg_start)
64 {
65  if(msg_start)
66  printf("\n%s", msg_start);
67 
68  if(reg)
69  printf("gregion[%d], ", *reg);
70 
71  printf("element %ld: ", (long int) elem.global_element_index(SF::NBR_ELEM_REF));
72 
73  for(int i=0; i<elem.num_nodes(); i++)
74  printf("%d ", elem.global_node(i, SF::NBR_REF));
75 
76  printf("\nVolume: %g", vol);
77 
78  if(eVals) {
79  SF::Point fib = elem.fiber(), she = elem.sheet();
80  printf(", fib: (%g, %g, %g), she: (%g, %g, %g), G: (%g, %g, %g)",
81  fib.x, fib.y, fib.z, she.x, she.y, she.z, eVals[0]*1e3, eVals[1]*1e3, eVals[2]*1e3);
82  }
83 
84  printf("\n\n");
85 }
86 
88  const bool is_bath,
89  double* evals)
90 {
91  if(is_bath == false) {
92  switch(emat.g) {
93  case intra_cond:
94  evals[0] = emat.InVal[0], evals[1] = emat.InVal[1], evals[2] = emat.InVal[2];
95  break;
96 
97  case extra_cond:
98  evals[0] = emat.ExVal[0], evals[1] = emat.ExVal[1], evals[2] = emat.ExVal[2];
99  break;
100 
101  case sum_cond:
102  evals[0] = emat.InVal[0] + emat.ExVal[0],
103  evals[1] = emat.InVal[1] + emat.ExVal[1],
104  evals[2] = emat.InVal[2] + emat.ExVal[2];
105  break;
106 
107  case para_cond:
108  evals[0] = (emat.InVal[0] * emat.ExVal[0]) / (emat.InVal[0] + emat.ExVal[0]),
109  evals[1] = (emat.InVal[1] * emat.ExVal[1]) / (emat.InVal[1] + emat.ExVal[1]),
110  evals[2] = (emat.InVal[2] * emat.ExVal[2]) / (emat.InVal[2] + emat.ExVal[2]);
111  break;
112  }
113  }
114  else
115  evals[0] = emat.BathVal[0], evals[1] = emat.BathVal[1], evals[2] = emat.BathVal[2];
116 }
117 
118 void get_conductivity(const double* evals,
119  const SF::Point *f,
120  const SF::Point *s,
121  SF::dmat<double> & cond)
122 {
123  cond.set_size(3,3);
124 
125  // longitudonal, transversal, and normal direction eigenvalues
126  double gl = evals[0], gt = evals[1], gn = evals[2];
127 
128  if(s == nullptr) {
129  // we compute the non-orthotropic conductivity:
130  // sigma = gl f f^T + gt (I - f f^T)
131  // = gt I + (gl - gt) f f^T
132  SF::outer_prod(*f, *f, gl - gt, cond[0]);
133  cond[0][0] += gt, cond[1][1] += gt, cond[2][2] += gt;
134  }
135  else {
136  // we compute the orthotropic conductivity:
137  // sigma = gl f f^T + gt s s^T + gn n n^T
138  SF::Point n = cross(*f, *s);
139 
140  SF::outer_prod(*f, *f, gl, cond[0], false);
141  SF::outer_prod(*s, *s, gt, cond[0], true);
142  SF::outer_prod( n, n, gn, cond[0], true);
143  }
144 }
145 
146 void print_points(SF::Point* pts, int npts, int eidx)
147 {
148  for(int i=0; i<npts; i++){
149  if(pts[i].x!=pts[i].x or pts[i].y!=pts[i].y or pts[i].z!=pts[i].z){
150  printf(" eidx = %d print_points ( %g %g %g ) \n", eidx, pts[i].x, pts[i].y, pts[i].z);
151  printf("\n");
152  }
153  }
154 }
155 
157 {
158  const size_t eidx = elem.element_index();
159  SF::elem_t type = elem.type();
160  mesh_int_t nnodes = elem.num_nodes();
161  const int int_order = get_preferred_int_order(type, ELEC_STIFFNESS);
162 
163  const int reg = material.regionIDs.size() ? material.regionIDs[eidx] : 0;
164  const RegionSpecs & reg_spec = material.regions[reg];
165 
166  assert(reg_spec.material->material_type == ElecMat);
167  const elecMaterial & mat = *static_cast<elecMaterial*>(reg_spec.material);
168  SF::Point fib = elem.fiber(), she = elem.sheet();
169 
170  bool is_bath = SF::inner_prod(fib, fib) < 0.01; // bath elements have no fiber direction
171  bool has_sheet = elem.has_sheet(); // do we have a sheet direciton
172 
173  double eVals[3];
174  get_conductivity_evals(mat, is_bath, eVals);
175 
176  // if we have a region scale
177  if(material.el_scale.size()) {
178  eVals[0] *= material.el_scale[eidx];
179  eVals[1] *= material.el_scale[eidx];
180  eVals[2] *= material.el_scale[eidx];
181  }
182 
183  if(is_bath) {
184  fib = {1, 0, 0};
185  she = {0, 1, 0};
186  }
187 
188  SF::get_transformed_pts(elem, lpts, fib);
189 
190 #ifdef DEBUG_INT
191  print_points(lpts, nnodes);
192 #endif
193 
194  get_conductivity(eVals, &fib, has_sheet ? &she : nullptr, cond);
195 
196  int ndof = nnodes;
197  shape.set_size(4, ndof), gshape.set_size(4, ndof);
198  buff.assign(ndof, ndof, 0.0);
199 
200  int nint = 0;
201  elem.integration_points(int_order, ipts, w, nint);
202 
203  double vol = 0.0;
204  bool have_nan = false;
205 
206  // Integration loop
207  for(int iidx=0; iidx < nint; iidx++)
208  {
209  double detj;
210  // get shape function and its spatial derivatives
211  SF::reference_shape(type, ipts[iidx], shape);
212  SF::jacobian_matrix(shape, nnodes, lpts, J);
213  SF::invert_jacobian_matrix(type, J, detj);
214 
215 #ifdef DEBUG_INT
216  SF::dmat<double> Jbuff(3,3);
217  memcpy(Jbuff.data(), J, 9*sizeof(double));
218  std::string j_name = "J_" + std::to_string(iidx);
219  Jbuff.disp(j_name.c_str());
220 #endif
221 
222  SF::shape_deriv(J, shape, ndof, gshape);
223 
224  double dx = w[iidx] * detj;
225  vol += dx;
226 
227  for (int k = 0; k < ndof; k++) {
228  // derivatives of shape function k
229  double skx = gshape[1][k], sky = gshape[2][k], skz = gshape[3][k];
230  double cx = cond[0][0] * skx + cond[0][1] * sky + cond[0][2] * skz;
231  double cy = cond[1][0] * skx + cond[1][1] * sky + cond[1][2] * skz;
232  double cz = cond[2][0] * skx + cond[2][1] * sky + cond[2][2] * skz;
233 
234  for (int l = 0; l < ndof; l++) {
235  // derivatives of shape function l
236  double slx = gshape[1][l], sly = gshape[2][l], slz = gshape[3][l];
237  double val = - (cx*slx + cy * sly + cz * slz) * dx;
238  buff[k][l] += val;
239 
240  if(val != val) have_nan = true;
241  }
242  }
243  }
244 
245 
246 #ifdef DEBUG_INT
247  printf("\n vol: %g \n\n", vol);
248  buff.disp("Element Matrix:");
249 #endif
250 
251  if(have_nan) {
252  print_element_info(elem, &reg, vol, eVals, "NaN detected: ");
253  exit(EXIT_FAILURE);
254  }
255 }
256 
258 {
259  row_dpn = 1;
260  col_dpn = 1;
261 }
262 
263 
265 {
266  SF::elem_t type = elem.type();
267  mesh_int_t nnodes = elem.num_nodes();
268  const int int_order = get_preferred_int_order(type, ELEC_MASS);
269 
270  int ndof = nnodes;
271  shape.set_size(4, ndof);
272  buff.assign(ndof, ndof, 0.0);
273 
274  SF::Point fib = {1, 0, 0}; // this is a dummy fibre to satisfy get_transformed_pts()
275 
276  double vol = 0.0;
277  bool have_nan = false;
278  int nint;
279 
280  elem.integration_points(int_order, ipts, w, nint);
282 
283  // Integration loop
284  for(int iidx=0; iidx < nint; iidx++)
285  {
286  double detj;
287  // get shape function and its spatial derivatives
288  SF::reference_shape(type, ipts[iidx], shape);
289  SF::jacobian_matrix(shape, nnodes, lpts, J);
290  SF::invert_jacobian_matrix(type, J, detj);
291 
292  double dx = w[iidx] * detj;
293  vol += dx;
294 
295  for (int k = 0; k < ndof; k++) {
296  double sk = shape[0][k];
297  for (int l = 0; l < ndof; l++) {
298  double sl = shape[0][l];
299  double val = sk * sl * dx;
300  buff[k][l] += val;
301 
302  if(val != val) have_nan = true;
303  }
304  }
305  }
306 
307 #ifdef DEBUG_INT
308  printf("\n vol: %g \n\n", vol);
309  buff.disp("Mass Matrix:");
310 #endif
311 
312  if(have_nan) {
313  print_element_info(elem, NULL, vol, NULL, "NaN detected: ");
314  exit(EXIT_FAILURE);
315  }
316 }
317 
318 void mass_integrator::dpn(mesh_int_t & row_dpn, mesh_int_t & col_dpn)
319 {
320  row_dpn = 1;
321  col_dpn = 1;
322 }
323 
324 } // namespace opencarp
325 
int mesh_int_t
Definition: SF_container.h:46
void assign(const dmat< S > &m)
copy a mtrix.
Definition: dense_mat.hpp:135
void set_size(const short irows, const short icols)
set the matrix dimensions
Definition: dense_mat.hpp:73
S * data()
Definition: dense_mat.hpp:264
void disp(const char *name)
Definition: dense_mat.hpp:254
Comfort class. Provides getter functions to access the mesh member variables more comfortably.
Definition: SF_fem_utils.h:704
Point fiber() const
Get element fiber direction.
Definition: SF_fem_utils.h:850
void integration_points(const short order, Point *ip, double *w, int &nint) const
Definition: SF_fem_utils.h:917
elem_t type() const
Getter function for the element type.
Definition: SF_fem_utils.h:771
const T & global_node(short nidx) const
Access the connectivity information.
Definition: SF_fem_utils.h:805
size_t global_element_index() const
Get currently selected element index.
Definition: SF_fem_utils.h:895
size_t element_index() const
Get currently selected element index.
Definition: SF_fem_utils.h:885
bool has_sheet() const
Check if a sheet direction is present.
Definition: SF_fem_utils.h:875
T num_nodes() const
Getter function for the number of nodes.
Definition: SF_fem_utils.h:761
Point sheet() const
Get element sheet direction.
Definition: SF_fem_utils.h:862
size_t size() const
The current size of the vector.
Definition: SF_vector.h:104
void operator()(const SF::element_view< mesh_int_t, mesh_real_t > &elem, SF::dmat< double > &buff)
compute the element matrix for a given element.
void dpn(mesh_int_t &row_dpn, mesh_int_t &col_dpn)
return (by reference) the row and column dimensions
void dpn(mesh_int_t &row_dpn, mesh_int_t &col_dpn)
return (by reference) the row and column dimensions
void operator()(const SF::element_view< mesh_int_t, mesh_real_t > &elem, SF::dmat< double > &buff)
compute the element matrix for a given element.
#define ELEC_STIFFNESS
#define ELEC_MASS
void shape_deriv(const double *iJ, const dmat< double > &rshape, const int ndof, dmat< double > &shape)
Compute shape derivatives for an element, based on the shape derivatives of the associated reference ...
Definition: SF_fem_utils.h:674
double inner_prod(const Point &a, const Point &b)
Definition: SF_container.h:90
void outer_prod(const Point &a, const Point &b, const double s, double *buff, const bool add=false)
Definition: SF_container.h:95
void jacobian_matrix(const dmat< double > &rshape, const int npts, const Point *pts, double *J)
Compute Jacobian matrix from the real element to the reference element.
Definition: SF_fem_utils.h:610
elem_t
element type enum
Definition: SF_container.h:53
@ Line
Definition: SF_container.h:61
@ Tri
Definition: SF_container.h:60
@ Prism
Definition: SF_container.h:58
@ Pyramid
Definition: SF_container.h:57
@ Tetra
Definition: SF_container.h:54
@ Quad
Definition: SF_container.h:59
@ Hexa
Definition: SF_container.h:55
void get_transformed_pts(const element_view< T, S > &view, Point *loc_pts, Point &trsf_fibre, bool orthogonal=true)
void reference_shape(const elem_t type, const Point ip, dmat< double > &rshape)
Compute shape function and its derivatives on a reference element.
Definition: SF_fem_utils.h:383
@ NBR_ELEM_REF
The element numbering of the reference mesh (the one stored on HD).
Definition: SF_container.h:204
@ NBR_REF
The nodal numbering of the reference mesh (the one stored on HD).
Definition: SF_container.h:201
void invert_jacobian_matrix(const elem_t type, double *J, double &detJ)
Definition: SF_fem_utils.h:633
void print_element_info(const SF::element_view< mesh_int_t, mesh_real_t > &elem, const int *reg, const double vol, const double *eVals, const char *msg_start)
void print_points(SF::Point *pts, int npts, int eidx)
@ extra_cond
Definition: fem_types.h:43
@ sum_cond
Definition: fem_types.h:43
@ intra_cond
Definition: fem_types.h:43
@ para_cond
Definition: fem_types.h:43
vec3< V > cross(const vec3< V > &a, const vec3< V > &b)
Definition: vect.h:144
bool phys_defined(int physreg)
function to check if certain physics are defined
void get_conductivity_evals(const elecMaterial &emat, const bool is_bath, double *evals)
int get_preferred_int_order(SF::elem_t &etype, int mat_type)
void get_conductivity(const double *evals, const SF::Point *f, const SF::Point *s, SF::dmat< double > &cond)
@ ElecMat
Definition: fem_types.h:39
Interface to SlimFem.
#define PHYSREG_INTRA_ELEC
Definition: sf_interface.h:84
Point and vector struct.
Definition: SF_container.h:65
double y
Definition: SF_container.h:67
double z
Definition: SF_container.h:68
double x
Definition: SF_container.h:66
SF::vector< RegionSpecs > regions
array with region params
Definition: fem_types.h:126
SF::vector< double > el_scale
optionally provided per-element params scale
Definition: fem_types.h:127
SF::vector< int > regionIDs
elemental region IDs (multiple tags are grouped in one region)
Definition: fem_types.h:125
region based variations of arbitrary material parameters
Definition: fem_types.h:93
physMaterial * material
material parameter description
Definition: fem_types.h:98
double ExVal[3]
extracellular conductivity eigenvalues
Definition: fem_types.h:62
cond_t g
rule to build conductivity tensor
Definition: fem_types.h:64
double InVal[3]
intracellular conductivity eigenvalues
Definition: fem_types.h:61
double BathVal[3]
bath conductivity eigenvalues
Definition: fem_types.h:63
physMat_t material_type
ID of physics material.
Definition: fem_types.h:53