26 #include "libmesh/fe.h"
27 #include "libmesh/elem.h"
28 #include "libmesh/libmesh_logging.h"
29 #include "libmesh/fe_interface.h"
30 #include "libmesh/fe_macro.h"
31 #include "libmesh/fe_map.h"
32 #include "libmesh/fe_xyz_map.h"
33 #include "libmesh/inf_fe_map.h"
34 #include "libmesh/mesh_subdivision_support.h"
35 #include "libmesh/dense_matrix.h"
36 #include "libmesh/dense_vector.h"
37 #include "libmesh/tensor_value.h"
38 #include "libmesh/auto_ptr.h"
39 #include "libmesh/enum_elem_type.h"
40 #include "libmesh/int_range.h"
56 libmesh_error_msg(
"Unknown mapping type " << elem.
mapping_type());
65 calculations_started(false),
67 calculate_dxyz(false),
68 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
69 calculate_d2xyz(false),
71 jacobian_tolerance(jtol)
81 return libmesh_make_unique<FEXYZMap>();
84 return libmesh_make_unique<FEMap>();
90 template<
unsigned int Dim>
95 LOG_SCOPE(
"init_reference_to_physical_map()",
"FEMap");
100 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
105 libmesh_not_implemented();
110 const std::size_t n_qp = qp.size();
122 const unsigned int n_mapping_shape_functions =
127 this->
phi_map.resize (n_mapping_shape_functions);
131 this->
dphidxi_map.resize (n_mapping_shape_functions);
132 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
135 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
142 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
148 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
155 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
162 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
166 for (
unsigned int i=0; i<n_mapping_shape_functions; i++)
169 this->
phi_map[i].resize (n_qp);
174 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
190 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
209 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
212 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
221 for (
unsigned int i=0; i<n_mapping_shape_functions; i++)
222 for (std::size_t p=0; p<n_qp; p++)
224 shape_ptr(elem, mapping_order, i, qp[p],
false);
237 for (
unsigned int i=0; i<n_mapping_shape_functions; i++)
241 shape_ptr(elem, mapping_order, i, qp[0],
false);
245 shape_deriv_ptr(elem, mapping_order, i, 0, qp[0],
false);
247 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
250 shape_second_deriv_ptr(elem, mapping_order, i, 0, qp[0],
false);
251 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
252 for (std::size_t p=1; p<n_qp; p++)
256 shape_ptr(elem, mapping_order, i, qp[p],
false);
259 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
262 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
267 for (
unsigned int i=0; i<n_mapping_shape_functions; i++)
268 for (std::size_t p=0; p<n_qp; p++)
272 shape_ptr (elem, mapping_order, i, qp[p],
false);
274 this->
dphidxi_map[i][p] = shape_deriv_ptr (elem, mapping_order, i, 0, qp[p],
false);
275 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
277 this->
d2phidxi2_map[i][p] = shape_second_deriv_ptr (elem, mapping_order, i, 0, qp[p],
false);
278 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
291 for (
unsigned int i=0; i<n_mapping_shape_functions; i++)
295 shape_ptr (elem, mapping_order, i, qp[0],
false);
298 this->
dphidxi_map[i][0] = shape_deriv_ptr (elem, mapping_order, i, 0, qp[0],
false);
299 this->
dphideta_map[i][0] = shape_deriv_ptr (elem, mapping_order, i, 1, qp[0],
false);
301 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
304 this->
d2phidxi2_map[i][0] = shape_second_deriv_ptr (elem, mapping_order, i, 0, qp[0],
false);
305 this->
d2phidxideta_map[i][0] = shape_second_deriv_ptr (elem, mapping_order, i, 1, qp[0],
false);
306 this->
d2phideta2_map[i][0] = shape_second_deriv_ptr (elem, mapping_order, i, 2, qp[0],
false);
308 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
309 for (std::size_t p=1; p<n_qp; p++)
313 shape_ptr (elem, mapping_order, i, qp[p],
false);
319 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
326 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
331 for (
unsigned int i=0; i<n_mapping_shape_functions; i++)
332 for (std::size_t p=0; p<n_qp; p++)
336 shape_ptr(elem, mapping_order, i, qp[p],
false);
339 this->
dphidxi_map[i][p] = shape_deriv_ptr (elem, mapping_order, i, 0, qp[p],
false);
340 this->
dphideta_map[i][p] = shape_deriv_ptr (elem, mapping_order, i, 1, qp[p],
false);
342 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
345 this->
d2phidxi2_map[i][p] = shape_second_deriv_ptr (elem, mapping_order, i, 0, qp[p],
false);
346 this->
d2phidxideta_map[i][p] = shape_second_deriv_ptr (elem, mapping_order, i, 1, qp[p],
false);
347 this->
d2phideta2_map[i][p] = shape_second_deriv_ptr (elem, mapping_order, i, 2, qp[p],
false);
349 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
363 for (
unsigned int i=0; i<n_mapping_shape_functions; i++)
367 shape_ptr (elem, mapping_order, i, qp[0],
false);
370 this->
dphidxi_map[i][0] = shape_deriv_ptr (elem, mapping_order, i, 0, qp[0],
false);
371 this->
dphideta_map[i][0] = shape_deriv_ptr (elem, mapping_order, i, 1, qp[0],
false);
372 this->
dphidzeta_map[i][0] = shape_deriv_ptr (elem, mapping_order, i, 2, qp[0],
false);
374 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
377 this->
d2phidxi2_map[i][0] = shape_second_deriv_ptr (elem, mapping_order, i, 0, qp[0],
false);
378 this->
d2phidxideta_map[i][0] = shape_second_deriv_ptr (elem, mapping_order, i, 1, qp[0],
false);
379 this->
d2phideta2_map[i][0] = shape_second_deriv_ptr (elem, mapping_order, i, 2, qp[0],
false);
380 this->
d2phidxidzeta_map[i][0] = shape_second_deriv_ptr (elem, mapping_order, i, 3, qp[0],
false);
381 this->
d2phidetadzeta_map[i][0] = shape_second_deriv_ptr (elem, mapping_order, i, 4, qp[0],
false);
382 this->
d2phidzeta2_map[i][0] = shape_second_deriv_ptr (elem, mapping_order, i, 5, qp[0],
false);
384 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
385 for (std::size_t p=1; p<n_qp; p++)
389 shape_ptr (elem, mapping_order, i, qp[p],
false);
396 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
406 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
411 for (
unsigned int i=0; i<n_mapping_shape_functions; i++)
412 for (std::size_t p=0; p<n_qp; p++)
416 shape_ptr(elem, mapping_order, i, qp[p],
false);
419 this->
dphidxi_map[i][p] = shape_deriv_ptr (elem, mapping_order, i, 0, qp[p],
false);
420 this->
dphideta_map[i][p] = shape_deriv_ptr (elem, mapping_order, i, 1, qp[p],
false);
421 this->
dphidzeta_map[i][p] = shape_deriv_ptr (elem, mapping_order, i, 2, qp[p],
false);
423 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
426 this->
d2phidxi2_map[i][p] = shape_second_deriv_ptr (elem, mapping_order, i, 0, qp[p],
false);
427 this->
d2phidxideta_map[i][p] = shape_second_deriv_ptr (elem, mapping_order, i, 1, qp[p],
false);
428 this->
d2phideta2_map[i][p] = shape_second_deriv_ptr (elem, mapping_order, i, 2, qp[p],
false);
429 this->
d2phidxidzeta_map[i][p] = shape_second_deriv_ptr (elem, mapping_order, i, 3, qp[p],
false);
430 this->
d2phidetadzeta_map[i][p] = shape_second_deriv_ptr (elem, mapping_order, i, 4, qp[p],
false);
431 this->
d2phidzeta2_map[i][p] = shape_second_deriv_ptr (elem, mapping_order, i, 5, qp[p],
false);
433 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
440 libmesh_error_msg(
"Invalid Dim = " << Dim);
447 const std::vector<Real> & qw,
450 const std::vector<const Node *> & elem_nodes,
451 bool compute_second_derivatives)
455 #ifndef LIBMESH_ENABLE_SECOND_DERIVATIVES
460 libmesh_assert_equal_to(
phi_map.size(), elem_nodes.size());
470 xyz[p] = *elem_nodes[0];
488 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
503 const Point & elem_point = *elem_nodes[i];
506 xyz[p].add_scaled (elem_point,
phi_map[i][p] );
509 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
542 static bool failing =
false;
550 libmesh_error_msg(
"ERROR: negative Jacobian " \
561 libmesh_error_msg(
"ERROR: negative Jacobian " \
563 <<
" at point index " \
593 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
600 #elif LIBMESH_DIM == 2
619 static bool failing =
false;
625 libmesh_error_msg(
"Encountered invalid 1D element!");
645 #elif LIBMESH_DIM == 3
666 static bool failing =
false;
672 libmesh_error_msg(
"Encountered invalid 1D element!");
700 #endif //LIBMESH_DIM == 3
703 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
726 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
745 const Point & elem_point = *elem_nodes[i];
748 xyz[p].add_scaled (elem_point,
phi_map[i][p] );
755 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
783 jac[p] = (dx_dxi*dy_deta - dx_deta*dy_dxi);
790 static bool failing =
false;
798 libmesh_error_msg(
"ERROR: negative Jacobian " \
809 libmesh_error_msg(
"ERROR: negative Jacobian " \
811 <<
" at point index " \
830 const Real inv_jac = 1./
jac[p];
839 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
840 if (compute_second_derivatives)
843 #else // LIBMESH_DIM == 3
870 const Real g11 = (dx_dxi*dx_dxi +
874 const Real g12 = (dx_dxi*dx_deta +
878 const Real g21 = g12;
880 const Real g22 = (dx_deta*dx_deta +
884 const Real det = (g11*g22 - g12*g21);
891 static bool failing =
false;
899 libmesh_error_msg(
"ERROR: negative Jacobian " \
910 libmesh_error_msg(
"ERROR: negative Jacobian " \
912 <<
" at point index " \
927 const Real inv_det = 1./det;
932 const Real g11inv = g22*inv_det;
933 const Real g12inv = -g12*inv_det;
934 const Real g21inv = -g21*inv_det;
935 const Real g22inv = g11*inv_det;
937 dxidx_map[p] = g11inv*dx_dxi + g12inv*dx_deta;
938 dxidy_map[p] = g11inv*dy_dxi + g12inv*dy_deta;
939 dxidz_map[p] = g11inv*dz_dxi + g12inv*dz_deta;
941 detadx_map[p] = g21inv*dx_dxi + g22inv*dx_deta;
942 detady_map[p] = g21inv*dy_dxi + g22inv*dy_deta;
943 detadz_map[p] = g21inv*dz_dxi + g22inv*dz_deta;
945 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
964 JT(0,0) = dx_dxi; JT(0,1) = dy_dxi; JT(0,2) = dz_dxi;
965 JT(1,0) = dx_deta; JT(1,1) = dy_deta; JT(1,2) = dz_deta;
969 JTJinv(0,0) = g11inv; JTJinv(0,1) = g12inv;
970 JTJinv(1,0) = g21inv; JTJinv(1,1) = g22inv;
985 for (
unsigned s=0; s<3; ++s)
986 for (
unsigned t=s; t<3; ++t)
989 tmp1(0) = dxi(s)*dxi(t);
990 tmp1(1) = deta(s)*deta(t);
993 A.vector_mult(tmp2, tmp1);
996 Real alpha = dxi(s)*deta(t) + deta(s)*dxi(t);
999 for (
unsigned i=0; i<3; ++i)
1003 JT.vector_mult(tmp3, tmp2);
1017 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
1019 #endif // LIBMESH_DIM == 3
1044 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1070 const Point & elem_point = *elem_nodes[i];
1073 xyz[p].add_scaled (elem_point,
phi_map[i][p] );
1080 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1117 jac[p] = (dx_dxi*(dy_deta*dz_dzeta - dz_deta*dy_dzeta) +
1118 dy_dxi*(dz_deta*dx_dzeta - dx_deta*dz_dzeta) +
1119 dz_dxi*(dx_deta*dy_dzeta - dy_deta*dx_dzeta));
1126 static bool failing =
false;
1134 libmesh_error_msg(
"ERROR: negative Jacobian " \
1145 libmesh_error_msg(
"ERROR: negative Jacobian " \
1147 <<
" at point index " \
1166 const Real inv_jac = 1./
jac[p];
1168 dxidx_map[p] = (dy_deta*dz_dzeta - dz_deta*dy_dzeta)*inv_jac;
1169 dxidy_map[p] = (dz_deta*dx_dzeta - dx_deta*dz_dzeta)*inv_jac;
1170 dxidz_map[p] = (dx_deta*dy_dzeta - dy_deta*dx_dzeta)*inv_jac;
1172 detadx_map[p] = (dz_dxi*dy_dzeta - dy_dxi*dz_dzeta )*inv_jac;
1173 detady_map[p] = (dx_dxi*dz_dzeta - dz_dxi*dx_dzeta )*inv_jac;
1174 detadz_map[p] = (dy_dxi*dx_dzeta - dx_dxi*dy_dzeta )*inv_jac;
1176 dzetadx_map[p] = (dy_dxi*dz_deta - dz_dxi*dy_deta )*inv_jac;
1177 dzetady_map[p] = (dz_dxi*dx_deta - dx_dxi*dz_deta )*inv_jac;
1178 dzetadz_map[p] = (dx_dxi*dy_deta - dy_dxi*dx_deta )*inv_jac;
1181 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1182 if (compute_second_derivatives)
1190 libmesh_error_msg(
"Invalid dim = " <<
dim);
1211 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1231 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1252 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1278 const std::vector<Real> & qw,
1282 LOG_SCOPE(
"compute_affine_map()",
"FEMap");
1286 const unsigned int n_qp = cast_int<unsigned int>(qw.size());
1294 for (
unsigned int i=0; i<
n_nodes; i++)
1302 for (
unsigned int p=1; p<n_qp; p++)
1311 for (
unsigned int p=1; p<n_qp; p++)
1317 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1328 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1341 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1352 JxW[p] =
JxW[0] / qw[0] * qw[p];
1359 const std::vector<Real> & qw)
1362 LOG_SCOPE(
"compute_null_map()",
"FEMap");
1364 const unsigned int n_qp = cast_int<unsigned int>(qw.size());
1370 for (
unsigned int p=1; p<n_qp; p++)
1382 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1397 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1413 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1434 const std::vector<Real> & qw,
1436 bool calculate_d2phi)
1449 #ifndef LIBMESH_ENABLE_SECOND_DERIVATIVES
1454 LOG_SCOPE(
"compute_map()",
"FEMap");
1458 const unsigned int n_qp = cast_int<unsigned int>(qw.size());
1467 libmesh_assert_equal_to (
dim, 2);
1468 const Tri3Subdivision * sd_elem = static_cast<const Tri3Subdivision *>(elem);
1480 for (
unsigned int p=0; p!=n_qp; p++)
1489 os <<
" [" << i <<
"]: " <<
JxW[i] << std::endl;
1497 os <<
" [" << i <<
"]: " <<
xyz[i];
1504 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
1507 std::set<unsigned> valid_indices;
1531 valid_indices.insert(0);
1533 #elif LIBMESH_DIM==2
1553 const unsigned tmp[3] = {0,1,3};
1554 valid_indices.insert(tmp, tmp+3);
1556 #elif LIBMESH_DIM==3
1576 const unsigned tmp[6] = {0,1,2,3,4,5};
1577 valid_indices.insert(tmp, tmp+6);
1584 for (
unsigned s=0; s<3; ++s)
1585 for (
unsigned t=s; t<3; ++t)
1587 if (valid_indices.count(ctr))
1594 v2(dxi(s)*deta(t) + deta(s)*dxi(t),
1595 dxi(s)*dzeta(t) + dzeta(s)*dxi(t),
1596 deta(s)*dzeta(t) + dzeta(s)*deta(t));
1604 if (LIBMESH_DIM > 1)
1607 if (LIBMESH_DIM > 2)
1617 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
1624 const Point & physical_point,
1625 const Real tolerance,
1629 libmesh_assert_greater_equal (tolerance, 0.);
1631 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
1638 LOG_SCOPE(
"inverse_map()",
"FEMap");
1642 Real inverse_map_error = 0.;
1655 unsigned int cnt = 0;
1659 const unsigned int max_cnt = 10;
1674 const Point delta = physical_point - physical_guess;
1717 const Real G = dxi*dxi;
1720 libmesh_assert_greater (G, 0.);
1722 const Real Ginv = 1./G;
1724 const Real dxidelta = dxi*delta;
1726 dp(0) = Ginv*dxidelta;
1766 G11 = dxi*dxi, G12 = dxi*deta,
1767 G21 = dxi*deta, G22 = deta*deta;
1770 const Real det = (G11*G22 - G12*G21);
1773 libmesh_assert_not_equal_to (det, 0.);
1775 const Real inv_det = 1./det;
1778 Ginv11 = G22*inv_det,
1779 Ginv12 = -G12*inv_det,
1781 Ginv21 = -G21*inv_det,
1782 Ginv22 = G11*inv_det;
1785 const Real dxidelta = dxi*delta;
1786 const Real detadelta = deta*delta;
1788 dp(0) = (Ginv11*dxidelta + Ginv12*detadelta);
1789 dp(1) = (Ginv21*dxidelta + Ginv22*detadelta);
1830 dxi(1), deta(1), dzeta(1),
1831 dxi(2), deta(2), dzeta(2)).solve(delta, dp);
1846 libMesh::err <<
"ERROR: Newton scheme encountered a singular Jacobian in element: "
1852 libmesh_error_msg(
"Exiting...");
1856 for (
unsigned int i=0; i !=
dim; ++i)
1873 libmesh_error_msg(
"Invalid dim = " <<
dim);
1879 inverse_map_error = dp.
norm();
1908 libMesh::err <<
"WARNING: Newton scheme has not converged in "
1909 << cnt <<
" iterations:" << std::endl
1910 <<
" physical_point="
1912 <<
" physical_guess="
1918 <<
" error=" << inverse_map_error
1919 <<
" in element " << elem->
id()
1926 libmesh_do_once(
libMesh::err <<
"WARNING: At least one element took more than "
1928 <<
" iterations to converge in inverse_map()...\n"
1929 <<
"Rerun in devel/dbg mode for more details."
1934 if (cnt > 2*max_cnt)
1936 libMesh::err <<
"ERROR: Newton scheme FAILED to converge in "
1938 <<
" iterations in element "
1940 <<
" for physical point = "
1946 libmesh_error_msg(
"Exiting...");
1954 for (
unsigned int i=0; i !=
dim; ++i)
1961 while (inverse_map_error > tolerance);
1974 const Point diff = physical_point - check;
1976 if (diff.
norm() > tolerance)
1996 libMesh::err <<
"WARNING: inverse_map of physical point "
1998 <<
" is not on element." <<
'\n';
2012 const std::vector<Point> & physical_points,
2013 std::vector<Point> & reference_points,
2014 const Real tolerance,
2017 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2028 const std::size_t n_points = physical_points.size();
2032 reference_points.resize(n_points);
2036 for (std::size_t p=0; p<n_points; p++)
2037 reference_points[p] =
2045 const Point & reference_point)
2049 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2059 const FEType fe_type (order, mapping_family);
2067 for (
unsigned int i=0; i<n_sf; i++)
2069 shape_ptr(elem, order, i, reference_point,
false));
2078 const unsigned int j,
2079 const Point & reference_point)
2083 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
2085 libmesh_not_implemented();
2093 const FEType fe_type (order, mapping_family);
2100 for (
unsigned int i=0; i<n_sf; i++)
2102 shape_deriv_ptr(elem, order, i, j, reference_point,
2111 template void FEMap::init_reference_to_physical_map<0>(
const std::vector<Point> &,
const Elem *);
2112 template void FEMap::init_reference_to_physical_map<1>(
const std::vector<Point> &,
const Elem *);
2113 template void FEMap::init_reference_to_physical_map<2>(
const std::vector<Point> &,
const Elem *);
2114 template void FEMap::init_reference_to_physical_map<3>(
const std::vector<Point> &,
const Elem *);