20 static PetscErrorCode
ierr;
31 det_jac *= jac[3 *
i +
i];
35 if ((
j - (
j / 2) * 2) != 0)
47 SETERRQ1(PETSC_COMM_SELF, 1,
"info = %d", info);
50 SETERRQ1(PETSC_COMM_SELF, 1,
"info = %d", info);
78 ierr = PetscMalloc(point_num *
sizeof(
double), &
w);
80 ierr = PetscMalloc(dim_num * point_num *
sizeof(
double), &x);
85 for (point = 0; point < point_num; point++) {
86 G_TRI_X[point] = x[0 + point * dim_num];
87 G_TRI_W[point] =
w[point];
122 ierr = PetscMalloc(point_num *
sizeof(
double), &
w);
124 ierr = PetscMalloc(dim_num * point_num *
sizeof(
double), &x);
129 for (point = 0; point < point_num; point++) {
130 G_TRI_X[point] = x[0 + point * dim_num];
131 G_TRI_Y[point] = x[1 + point * dim_num];
132 G_TRI_W[point] =
w[point];
162 ierr = PetscMalloc(point_num *
sizeof(
double), &
w);
164 ierr = PetscMalloc(dim_num * point_num *
sizeof(
double), &x);
168 for (point = 0; point < point_num; point++) {
169 G_TET_X[point] = x[0 + point * dim_num];
170 G_TET_Y[point] = x[1 + point * dim_num];
171 G_TET_Z[point] = x[2 + point * dim_num];
172 G_TET_W[point] =
w[point];
182 PetscErrorCode
ShapeMBTRI(
double *
N,
const double *X,
const double *Y,
186 for (; ii < G_DIM; ii++) {
187 double x = X[ii], y = Y[ii];
205 double *normal,
double *s1,
double *s2) {
211 for (; ii < 3; ii++) {
212 diffX_ksi[ii] = cblas_ddot(3, &coords[ii], 3, &diffN[0], 2);
213 diffX_eta[ii] = cblas_ddot(3, &coords[ii], 3, &diffN[1], 2);
216 cblas_dcopy(3, diffX_ksi, 1, s1, 1);
219 cblas_dcopy(3, diffX_eta, 1, s2, 1);
221 double Spin_diffX_ksi[9];
222 ierr =
Spin(Spin_diffX_ksi, diffX_ksi);
224 cblas_dgemv(CblasRowMajor, CblasNoTrans, 3, 3, 1., Spin_diffX_ksi, 3,
225 diffX_eta, 1, 0., normal, 1);
238 double *diff_normal) {
245 for (; ii < 3; ii++) {
246 diffX_ksi[ii] = cblas_ddot(3, &coords[ii], 3, &diffN[0], 2);
247 diffX_eta[ii] = cblas_ddot(3, &coords[ii], 3, &diffN[1], 2);
249 double Spin_diffX_ksi[9];
250 ierr =
Spin(Spin_diffX_ksi, diffX_ksi);
252 double Spin_diffX_eta[9];
253 ierr =
Spin(Spin_diffX_eta, diffX_eta);
256 bzero(B_ksi, 3 * 9 *
sizeof(
double));
258 bzero(B_eta, 3 * 9 *
sizeof(
double));
276 for (; ii < 3; ii++) {
277 cblas_dcopy(3, &diffN[0], 2, &B_ksi[ii * 9 + ii], 3);
278 cblas_dcopy(3, &diffN[1], 2, &B_eta[ii * 9 + ii], 3);
280 cblas_dgemm(CblasRowMajor, CblasNoTrans, CblasNoTrans, 3, 9, 3, +1.,
281 Spin_diffX_ksi, 3, B_eta, 9, 0., diff_normal, 9);
282 cblas_dgemm(CblasRowMajor, CblasNoTrans, CblasNoTrans, 3, 9, 3, -1.,
283 Spin_diffX_eta, 3, B_ksi, 9, 1., diff_normal, 9);
288 PetscErrorCode
ShapeJacMBTET(
double *diffN,
const double *coords,
double *jac) {
291 bzero(jac,
sizeof(
double) * 9);
292 for (ii = 0; ii < 4; ii++)
293 for (jj = 0; jj < 3; jj++)
294 for (kk = 0; kk < 3; kk++)
295 jac[jj * 3 + kk] += diffN[ii * 3 + kk] * coords[ii * 3 + jj];
306 PetscErrorCode
ShapeMBTET(
double *
N,
const double *G_X,
const double *G_Y,
307 const double *G_Z,
int DIM) {
310 for (; ii < DIM; ii++) {
311 double x = G_X[ii], y = G_Y[ii], z = G_Z[ii];
336 const double *elem_coords,
337 const double *glob_coords,
338 double *loc_coords) {
345 cblas_ddot(4, &diffN[0 * 3 + 0], 3, &elem_coords[0 * 3 + 0], 3);
347 cblas_ddot(4, &diffN[0 * 3 + 1], 3, &elem_coords[0 * 3 + 0], 3);
349 cblas_ddot(4, &diffN[0 * 3 + 2], 3, &elem_coords[0 * 3 + 0], 3);
351 glob_coords[0] - cblas_ddot(4, &
N[0], 1, &elem_coords[0 * 3 + 0], 3);
355 cblas_ddot(4, &diffN[0 * 3 + 0], 3, &elem_coords[0 * 3 + 1], 3);
357 cblas_ddot(4, &diffN[0 * 3 + 1], 3, &elem_coords[0 * 3 + 1], 3);
359 cblas_ddot(4, &diffN[0 * 3 + 2], 3, &elem_coords[0 * 3 + 1], 3);
361 glob_coords[1] - cblas_ddot(4, &
N[0], 1, &elem_coords[0 * 3 + 1], 3);
365 cblas_ddot(4, &diffN[0 * 3 + 0], 3, &elem_coords[0 * 3 + 2], 3);
367 cblas_ddot(4, &diffN[0 * 3 + 1], 3, &elem_coords[0 * 3 + 2], 3);
369 cblas_ddot(4, &diffN[0 * 3 + 2], 3, &elem_coords[0 * 3 + 2], 3);
371 glob_coords[2] - cblas_ddot(4, &
N[0], 1, &elem_coords[0 * 3 + 2], 3);
376 SETERRQ1(PETSC_COMM_SELF, 1,
"info == %d", info);
381 const double *elem_coords,
382 const double *glob_coords,
383 double *loc_coords) {
388 A[0] = cblas_ddot(3, &diffN[0], 2, &elem_coords[0], 2);
389 A[1] = cblas_ddot(3, &diffN[1], 2, &elem_coords[0], 2);
390 loc_coords[0] = glob_coords[0] - cblas_ddot(3, &
N[0], 1, &elem_coords[0], 2);
393 A[2] = cblas_ddot(3, &diffN[0], 2, &elem_coords[1], 2);
394 A[3] = cblas_ddot(3, &diffN[1], 2, &elem_coords[1], 2);
395 loc_coords[1] = glob_coords[1] - cblas_ddot(3, &
N[0], 1, &elem_coords[1], 2);
398 double invA[2 * 2], detA;
399 detA =
A[0] *
A[3] -
A[1] *
A[2];
401 invA[0] =
A[3] * detA;
402 invA[1] = -1.0 *
A[1] * detA;
403 invA[2] = -1.0 *
A[2] * detA;
404 invA[3] =
A[0] * detA;
406 double loc_coords_new[2];
407 loc_coords_new[0] = invA[0] * loc_coords[0] + invA[1] * loc_coords[1];
408 loc_coords_new[1] = invA[2] * loc_coords[0] + invA[3] * loc_coords[1];
410 loc_coords[0] = loc_coords_new[0];
411 loc_coords[1] = loc_coords_new[1];
416 double *diffNinvJac) {
419 for (; ii < 4; ii++) {
420 cblas_dgemv(CblasRowMajor, CblasTrans, 3, 3, 1.,
invJac, 3, &diffN[ii * 3],
421 1, 0., &diffNinvJac[ii * 3], 1);
428 for (; row < 3; row++)
429 for (col = 0; col < 3; col++) {
430 F[3 * row + col] = cblas_ddot(4, &diffN[col], 3, &dofs[row], 3);
439 enum CBLAS_TRANSPOSE Trans) {
442 for (; ii < 4; ii++) {
444 for (jj = 0; jj < 3; jj++) {
445 tmp3[jj].
r = diffN[ii * 3 + jj];
448 cblas_zgemv(CblasRowMajor, Trans, 3, 3, &tmp1,
invJac, 3, tmp3, 1, &tmp2,
449 &diffNinvJac[ii * 3], 1);
456 double complex diffX_x, diffX_y, diffX_z;
457 double complex diffY_x, diffY_y, diffY_z;
458 diffX_x = diffX_y = diffX_z = 0.;
459 diffY_x = diffY_y = diffY_z = 0.;
461 for (ii = 0; ii < 3; ii++) {
463 (xcoords[3 * ii + 0].
r +
I * xcoords[3 * ii + 0].
i) * diffN[2 * ii + 0];
465 (xcoords[3 * ii + 1].
r +
I * xcoords[3 * ii + 1].
i) * diffN[2 * ii + 0];
467 (xcoords[3 * ii + 2].
r +
I * xcoords[3 * ii + 2].
i) * diffN[2 * ii + 0];
469 (xcoords[3 * ii + 0].
r +
I * xcoords[3 * ii + 0].
i) * diffN[2 * ii + 1];
471 (xcoords[3 * ii + 1].
r +
I * xcoords[3 * ii + 1].
i) * diffN[2 * ii + 1];
473 (xcoords[3 * ii + 2].
r +
I * xcoords[3 * ii + 2].
i) * diffN[2 * ii + 1];
476 tmp = diffX_y * diffY_z - diffX_z * diffY_y;
477 xnormal[0].
r = creal(tmp);
478 xnormal[0].
i = cimag(tmp);
479 tmp = diffX_z * diffY_x - diffX_x * diffY_z;
480 xnormal[1].
r = creal(tmp);
481 xnormal[1].
i = cimag(tmp);
482 tmp = diffX_x * diffY_y - diffX_y * diffY_x;
483 xnormal[2].
r = creal(tmp);
484 xnormal[2].
i = cimag(tmp);
491 for (; ii < 3; ii++) {
492 for (jj = 0; jj < 3; jj++) {
493 xA[3 * ii + jj].
r = reA[3 * ii + jj];
494 xA[3 * ii + jj].
i = imA[3 * ii + jj];
507 SETERRQ(PETSC_COMM_SELF, 1,
"info == 0");
510 SETERRQ(PETSC_COMM_SELF, 1,
"info == 0");
518 SETERRQ(PETSC_COMM_SELF, 1,
"info == 0");
522 SETERRQ(PETSC_COMM_SELF, 1,
"info == 0");
531 SETERRQ(PETSC_COMM_SELF, 1,
"lapack_zgetrf(3,3,xF,3,IPIV) != 0");
533 double complex det = 1;
536 det *= xF[3 *
i +
i].
r +
I * xF[3 *
i +
i].
i;
537 if (IPIV[
i] !=
i + 1)
540 if ((
j - (
j / 2) * 2) != 0)
542 (*det_xF).r = creal(det);
543 (*det_xF).i = cimag(det);
546 PetscErrorCode
Spin(
double *spinOmega,
double *vecOmega) {
548 bzero(spinOmega, 9 *
sizeof(
double));
549 spinOmega[0 * 3 + 1] = -vecOmega[2];
550 spinOmega[0 * 3 + 2] = +vecOmega[1];
551 spinOmega[1 * 3 + 0] = +vecOmega[2];
552 spinOmega[1 * 3 + 2] = -vecOmega[0];
553 spinOmega[2 * 3 + 0] = -vecOmega[1];
554 spinOmega[2 * 3 + 1] = +vecOmega[0];
561 for (; ii < 3; ii++) {
562 for (jj = 0; jj < 3; jj++) {
563 xA[3 * ii + jj].
r = reA[3 * ii + jj];
564 xA[3 * ii + jj].
i = imA[3 * ii + jj];
570 int order_approx,
int *order_edge_approx,
int order,
int *order_edge,
571 double *diffN,
double *diffN_face,
double *diffN_edge[],
double *dofs,
572 double *dofs_edge[],
double *dofs_face,
double *idofs,
double *idofs_edge[],
578 double complex diffX_x_node, diffX_y_node, diffX_z_node;
579 double complex diffY_x_node, diffY_y_node, diffY_z_node;
586 if (dofs != NULL || idofs != NULL) {
588 for (; nn < 3; nn++) {
590 diffX_x_node += dofs[3 * nn + 0] * diffN[2 * nn + 0];
591 diffX_y_node += dofs[3 * nn + 1] * diffN[2 * nn + 0];
592 diffX_z_node += dofs[3 * nn + 2] * diffN[2 * nn + 0];
593 diffY_x_node += dofs[3 * nn + 0] * diffN[2 * nn + 1];
594 diffY_y_node += dofs[3 * nn + 1] * diffN[2 * nn + 1];
595 diffY_z_node += dofs[3 * nn + 2] * diffN[2 * nn + 1];
598 diffX_x_node +=
I * idofs[3 * nn + 0] * diffN[2 * nn + 0];
599 diffX_y_node +=
I * idofs[3 * nn + 1] * diffN[2 * nn + 0];
600 diffX_z_node +=
I * idofs[3 * nn + 2] * diffN[2 * nn + 0];
601 diffY_x_node +=
I * idofs[3 * nn + 0] * diffN[2 * nn + 1];
602 diffY_y_node +=
I * idofs[3 * nn + 1] * diffN[2 * nn + 1];
603 diffY_z_node +=
I * idofs[3 * nn + 2] * diffN[2 * nn + 1];
607 double complex diffX_x, diffX_y, diffX_z;
608 double complex diffY_x, diffY_y, diffY_z;
609 diffX_x = diffX_x_node;
610 diffX_y = diffX_y_node;
611 diffX_z = diffX_z_node;
612 diffY_x = diffY_x_node;
613 diffY_y = diffY_y_node;
614 diffY_z = diffY_z_node;
615 if (dofs_face != NULL || idofs_face != NULL) {
618 if (nb_dofs_face > 0) {
619 if (dofs_face != NULL) {
620 diffX_x += cblas_ddot(nb_dofs_face, &dofs_face[0], 3,
621 &diffN_face[gg * 2 * nb_dofs_approx_face + 0], 2);
622 diffX_y += cblas_ddot(nb_dofs_face, &dofs_face[1], 3,
623 &diffN_face[gg * 2 * nb_dofs_approx_face + 0], 2);
624 diffX_z += cblas_ddot(nb_dofs_face, &dofs_face[2], 3,
625 &diffN_face[gg * 2 * nb_dofs_approx_face + 0], 2);
626 diffY_x += cblas_ddot(nb_dofs_face, &dofs_face[0], 3,
627 &diffN_face[gg * 2 * nb_dofs_approx_face + 1], 2);
628 diffY_y += cblas_ddot(nb_dofs_face, &dofs_face[1], 3,
629 &diffN_face[gg * 2 * nb_dofs_approx_face + 1], 2);
630 diffY_z += cblas_ddot(nb_dofs_face, &dofs_face[2], 3,
631 &diffN_face[gg * 2 * nb_dofs_approx_face + 1], 2);
633 if (idofs_face != NULL) {
635 I * cblas_ddot(nb_dofs_face, &idofs_face[0], 3,
636 &diffN_face[gg * 2 * nb_dofs_approx_face + 0], 2);
638 I * cblas_ddot(nb_dofs_face, &idofs_face[1], 3,
639 &diffN_face[gg * 2 * nb_dofs_approx_face + 0], 2);
641 I * cblas_ddot(nb_dofs_face, &idofs_face[2], 3,
642 &diffN_face[gg * 2 * nb_dofs_approx_face + 0], 2);
644 I * cblas_ddot(nb_dofs_face, &idofs_face[0], 3,
645 &diffN_face[gg * 2 * nb_dofs_approx_face + 1], 2);
647 I * cblas_ddot(nb_dofs_face, &idofs_face[1], 3,
648 &diffN_face[gg * 2 * nb_dofs_approx_face + 1], 2);
650 I * cblas_ddot(nb_dofs_face, &idofs_face[2], 3,
651 &diffN_face[gg * 2 * nb_dofs_approx_face + 1], 2);
656 if (dofs_edge != NULL || idofs_edge != NULL) {
657 for (; ee < 3; ee++) {
658 int nb_dofs_edge =
NBEDGE_H1(order_edge[ee]);
659 int nb_dofs_approx_edge =
NBEDGE_H1(order_edge_approx[ee]);
660 if (nb_dofs_edge > 0) {
661 if (dofs_edge != NULL) {
662 if (dofs_edge[ee] != NULL) {
663 diffX_x += cblas_ddot(
664 nb_dofs_edge, &(dofs_edge[ee])[0], 3,
665 &(diffN_edge[ee])[gg * 2 * nb_dofs_approx_edge + 0], 2);
666 diffX_y += cblas_ddot(
667 nb_dofs_edge, &(dofs_edge[ee])[1], 3,
668 &(diffN_edge[ee])[gg * 2 * nb_dofs_approx_edge + 0], 2);
669 diffX_z += cblas_ddot(
670 nb_dofs_edge, &(dofs_edge[ee])[2], 3,
671 &(diffN_edge[ee])[gg * 2 * nb_dofs_approx_edge + 0], 2);
672 diffY_x += cblas_ddot(
673 nb_dofs_edge, &(dofs_edge[ee])[0], 3,
674 &(diffN_edge[ee])[gg * 2 * nb_dofs_approx_edge + 1], 2);
675 diffY_y += cblas_ddot(
676 nb_dofs_edge, &(dofs_edge[ee])[1], 3,
677 &(diffN_edge[ee])[gg * 2 * nb_dofs_approx_edge + 1], 2);
678 diffY_z += cblas_ddot(
679 nb_dofs_edge, &(dofs_edge[ee])[2], 3,
680 &(diffN_edge[ee])[gg * 2 * nb_dofs_approx_edge + 1], 2);
683 if (idofs_edge != NULL) {
684 if (idofs_edge[ee] == NULL)
688 nb_dofs_edge, &(idofs_edge[ee])[0], 3,
689 &(diffN_edge[ee])[gg * 2 * nb_dofs_approx_edge + 0], 2);
692 nb_dofs_edge, &(idofs_edge[ee])[1], 3,
693 &(diffN_edge[ee])[gg * 2 * nb_dofs_approx_edge + 0], 2);
696 nb_dofs_edge, &(idofs_edge[ee])[2], 3,
697 &(diffN_edge[ee])[gg * 2 * nb_dofs_approx_edge + 0], 2);
700 nb_dofs_edge, &(idofs_edge[ee])[0], 3,
701 &(diffN_edge[ee])[gg * 2 * nb_dofs_approx_edge + 1], 2);
704 nb_dofs_edge, &(idofs_edge[ee])[1], 3,
705 &(diffN_edge[ee])[gg * 2 * nb_dofs_approx_edge + 1], 2);
708 nb_dofs_edge, &(idofs_edge[ee])[2], 3,
709 &(diffN_edge[ee])[gg * 2 * nb_dofs_approx_edge + 1], 2);
714 double complex normal[3];
715 normal[0] = diffX_y * diffY_z - diffX_z * diffY_y;
716 normal[1] = diffX_z * diffY_x - diffX_x * diffY_z;
717 normal[2] = diffX_x * diffY_y - diffX_y * diffY_x;
719 for (;
dd < 3;
dd++) {
720 xnormal[
dd].
r = creal(normal[
dd]);
721 xnormal[
dd].
i = cimag(normal[
dd]);
724 xs1[0].
r = creal(diffX_x);
725 xs1[0].
i = cimag(diffX_x);
726 xs1[1].
r = creal(diffX_y);
727 xs1[1].
i = cimag(diffX_y);
728 xs1[2].
r = creal(diffX_z);
729 xs1[2].
i = cimag(diffX_z);
732 xs2[0].
r = creal(diffY_x);
733 xs2[0].
i = cimag(diffY_x);
734 xs2[1].
r = creal(diffY_y);
735 xs2[1].
i = cimag(diffY_y);
736 xs2[2].
r = creal(diffY_z);
737 xs2[2].
i = cimag(diffY_z);
745 complex
double xnrm2_normal = csqrt(cpow(xnormal[0].
r +
I * xnormal[0].
i, 2) +
746 cpow(xnormal[1].
r +
I * xnormal[1].
i, 2) +
747 cpow(xnormal[2].
r +
I * xnormal[2].
i, 2));
749 for (;
dd < 3;
dd++) {
750 complex
double s1 = (xs1[
dd].
r +
I * xs1[
dd].
i) * xnrm2_normal;
751 complex
double s2 = (xs2[
dd].
r +
I * xs2[
dd].
i) * xnrm2_normal;
752 xs1[
dd].
r = creal(s1);
753 xs1[
dd].
i = cimag(s1);
754 xs2[
dd].
r = creal(s2);
755 xs2[
dd].
i = cimag(s2);
764 for (; ii < DIM; ii++) {
780 PetscErrorCode
ShapeMBTRIQ(
double *
N,
const double *X,
const double *Y,
784 for (; ii < G_DIM; ii++) {
785 double x = X[ii], y = Y[ii];
799 for (; ii < G_DIM; ii++) {
800 double x = X[ii], y = Y[ii];
818 #define N_MBTETQ0(x, y, z) ((2. * (1. - x - y - z) - 1.) * (1. - x - y - z))
819 #define N_MBTETQ1(x, y, z) ((2. * x - 1.) * x)
820 #define N_MBTETQ2(x, y, z) ((2. * y - 1.) * y)
821 #define N_MBTETQ3(x, y, z) ((2. * z - 1.) * z)
822 #define N_MBTETQ4(x, y, z) (4. * (1. - x - y - z) * x)
823 #define N_MBTETQ5(x, y, z) (4. * x * y)
824 #define N_MBTETQ6(x, y, z) (4. * (1. - x - y - z) * y)
825 #define N_MBTETQ7(x, y, z) (4. * (1. - x - y - z) * z)
826 #define N_MBTETQ8(x, y, z) (4. * x * z)
827 #define N_MBTETQ9(x, y, z) (4. * y * z)
828 #define diffN_MBTETQ0x(x, y, z) (-3. + 4. * x + 4. * y + 4. * z)
829 #define diffN_MBTETQ0y(x, y, z) (-3. + 4. * x + 4. * y + 4. * z)
830 #define diffN_MBTETQ0z(x, y, z) (-3. + 4. * x + 4. * y + 4. * z)
831 #define diffN_MBTETQ1x(x, y, z) (4. * x - 1.)
832 #define diffN_MBTETQ1y(x, y, z) (0.)
833 #define diffN_MBTETQ1z(x, y, z) (0.)
834 #define diffN_MBTETQ2x(x, y, z) (0.)
835 #define diffN_MBTETQ2y(x, y, z) (4. * y - 1.)
836 #define diffN_MBTETQ2z(x, y, z) (0.)
837 #define diffN_MBTETQ3x(x, y, z) (0.)
838 #define diffN_MBTETQ3y(x, y, z) (0.)
839 #define diffN_MBTETQ3z(x, y, z) (4. * z - 1.)
840 #define diffN_MBTETQ4x(x, y, z) (-8. * x + 4. - 4. * y - 4. * z)
841 #define diffN_MBTETQ4y(x, y, z) (-4. * x)
842 #define diffN_MBTETQ4z(x, y, z) (-4. * x)
843 #define diffN_MBTETQ5x(x, y, z) (4. * y)
844 #define diffN_MBTETQ5y(x, y, z) (4. * x)
845 #define diffN_MBTETQ5z(x, y, z) (0.)
846 #define diffN_MBTETQ6x(x, y, z) (-4. * y)
847 #define diffN_MBTETQ6y(x, y, z) (-8. * y + 4. - 4. * x - 4. * z)
848 #define diffN_MBTETQ6z(x, y, z) (-4. * y)
849 #define diffN_MBTETQ7x(x, y, z) (-4. * z)
850 #define diffN_MBTETQ7y(x, y, z) (-4. * z)
851 #define diffN_MBTETQ7z(x, y, z) (-8. * z + 4. - 4. * x - 4. * y)
852 #define diffN_MBTETQ8x(x, y, z) (4. * z)
853 #define diffN_MBTETQ8y(x, y, z) (0.)
854 #define diffN_MBTETQ8z(x, y, z) (4. * x)
855 #define diffN_MBTETQ9x(x, y, z) (0.)
856 #define diffN_MBTETQ9y(x, y, z) (4. * z)
857 #define diffN_MBTETQ9z(x, y, z) (4. * y)
909 const double *Z,
const int G_DIM) {
912 for (; ii < G_DIM; ii++) {
913 double x = X[ii], y = Y[ii], z = Z[ii];
928 const double *Y,
const double *Z,
932 for (; ii < G_DIM; ii++) {
933 double x = X[ii], y = Y[ii], z = Z[ii];
971 bzero(Jac,
sizeof(
double) * 9);
972 for (ii = 0; ii < 10; ii++)
973 for (jj = 0; jj < 3; jj++)
974 for (kk = 0; kk < 3; kk++)
975 Jac[jj * 3 + kk] += diffN[ii * 3 + kk] * coords[ii * 3 + jj];
980 const double *diffN,
const double *coords,
986 for (; ii < G_DIM; ii++) {
998 double detJac_at_Gauss_Points[G_DIM];
1002 for (; ii < G_DIM; ii++) {
1003 vol += G_TET_W[ii] * (detJac_at_Gauss_Points[ii]) / 6;
1008 const double *elem_coords,
1009 const double *glob_coords,
1010 double *loc_coords,
const double eps) {
1015 float NORM_dR = 1000.;
1019 R[0] = glob_coords[0] - cblas_ddot(10, &
N[0], 1, &elem_coords[0], 3);
1020 R[1] = glob_coords[1] - cblas_ddot(10, &
N[0], 1, &elem_coords[1], 3);
1021 R[2] = glob_coords[2] - cblas_ddot(10, &
N[0], 1, &elem_coords[2], 3);
1022 NORM_R0 = cblas_dnrm2(3, &
R[0], 1);
1023 while ((NORM_dR / NORM_R0) >
eps) {
1026 A[0 + 3 * 0] = cblas_ddot(10, &diffN[0 * 3 + 0], 3, &elem_coords[0], 3);
1027 A[0 + 3 * 1] = cblas_ddot(10, &diffN[0 * 3 + 1], 3, &elem_coords[0], 3);
1028 A[0 + 3 * 2] = cblas_ddot(10, &diffN[0 * 3 + 2], 3, &elem_coords[0], 3);
1029 R[0] = glob_coords[0] - cblas_ddot(10, &
N[0], 1, &elem_coords[0], 3);
1031 A[1 + 3 * 0] = cblas_ddot(10, &diffN[0 * 3 + 0], 3, &elem_coords[1], 3);
1032 A[1 + 3 * 1] = cblas_ddot(10, &diffN[0 * 3 + 1], 3, &elem_coords[1], 3);
1033 A[1 + 3 * 2] = cblas_ddot(10, &diffN[0 * 3 + 2], 3, &elem_coords[1], 3);
1034 R[1] = glob_coords[1] - cblas_ddot(10, &
N[0], 1, &elem_coords[1], 3);
1037 cblas_ddot(10, &diffN[0 * 3 + 0], 3, &elem_coords[0 * 3 + 2], 3);
1039 cblas_ddot(10, &diffN[0 * 3 + 1], 3, &elem_coords[0 * 3 + 2], 3);
1041 cblas_ddot(10, &diffN[0 * 3 + 2], 3, &elem_coords[0 * 3 + 2], 3);
1042 R[2] = glob_coords[2] - cblas_ddot(10, &
N[0], 1, &elem_coords[2], 3);
1046 cblas_daxpy(3, 1.,
R, 1, loc_coords, 1);
1047 NORM_dR = cblas_dnrm2(3, &
R[0], 1);
1048 ShapeMBTETQ(
N, loc_coords[0], loc_coords[1], loc_coords[2]);