32 #define ELEC_STIFFNESS 1
63 const char* msg_start)
66 printf(
"\n%s", msg_start);
69 printf(
"gregion[%d], ", *reg);
76 printf(
"\nVolume: %g", vol);
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);
91 if(is_bath ==
false) {
94 evals[0] = emat.
InVal[0], evals[1] = emat.
InVal[1], evals[2] = emat.
InVal[2];
98 evals[0] = emat.
ExVal[0], evals[1] = emat.
ExVal[1], evals[2] = emat.
ExVal[2];
126 double gl = evals[0], gt = evals[1], gn = evals[2];
133 cond[0][0] += gt, cond[1][1] += gt, cond[2][2] += gt;
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);
178 eVals[0] *= material.
el_scale[eidx];
179 eVals[1] *= material.
el_scale[eidx];
180 eVals[2] *= material.
el_scale[eidx];
198 buff.
assign(ndof, ndof, 0.0);
204 bool have_nan =
false;
207 for(
int iidx=0; iidx < nint; iidx++)
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());
224 double dx = w[iidx] * detj;
227 for (
int k = 0; k < ndof; 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;
234 for (
int l = 0; l < ndof; 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;
240 if(val != val) have_nan =
true;
247 printf(
"\n vol: %g \n\n", vol);
248 buff.
disp(
"Element Matrix:");
272 buff.
assign(ndof, ndof, 0.0);
277 bool have_nan =
false;
284 for(
int iidx=0; iidx < nint; iidx++)
292 double dx = w[iidx] * detj;
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;
302 if(val != val) have_nan =
true;
308 printf(
"\n vol: %g \n\n", vol);
309 buff.
disp(
"Mass Matrix:");
void assign(const dmat< S > &m)
copy a mtrix.
void set_size(const short irows, const short icols)
set the matrix dimensions
void disp(const char *name)
Comfort class. Provides getter functions to access the mesh member variables more comfortably.
Point fiber() const
Get element fiber direction.
void integration_points(const short order, Point *ip, double *w, int &nint) const
elem_t type() const
Getter function for the element type.
const T & global_node(short nidx) const
Access the connectivity information.
size_t global_element_index() const
Get currently selected element index.
size_t element_index() const
Get currently selected element index.
bool has_sheet() const
Check if a sheet direction is present.
T num_nodes() const
Getter function for the number of nodes.
Point sheet() const
Get element sheet direction.
size_t size() const
The current size of the vector.
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.
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 ...
double inner_prod(const Point &a, const Point &b)
void outer_prod(const Point &a, const Point &b, const double s, double *buff, const bool add=false)
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.
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.
@ NBR_ELEM_REF
The element numbering of the reference mesh (the one stored on HD).
@ NBR_REF
The nodal numbering of the reference mesh (the one stored on HD).
void invert_jacobian_matrix(const elem_t type, double *J, double &detJ)
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)
vec3< V > cross(const vec3< V > &a, const vec3< V > &b)
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)
#define PHYSREG_INTRA_ELEC
SF::vector< RegionSpecs > regions
array with region params
SF::vector< double > el_scale
optionally provided per-element params scale
SF::vector< int > regionIDs
elemental region IDs (multiple tags are grouped in one region)
region based variations of arbitrary material parameters
physMaterial * material
material parameter description
double ExVal[3]
extracellular conductivity eigenvalues
cond_t g
rule to build conductivity tensor
double InVal[3]
intracellular conductivity eigenvalues
double BathVal[3]
bath conductivity eigenvalues
physMat_t material_type
ID of physics material.