void FEMap::compute_face_map(int dim, const std::vector<Real>& qw, const Elem* side) { libmesh_assert (side != NULL); START_LOG("compute_face_map()", "FEMap"); // The number of quadrature points. const unsigned int n_qp = qw.size(); switch (dim) { case 1: { // A 1D finite element, currently assumed to be in 1D space // This means the boundary is a "0D finite element", a // NODEELEM. // Resize the vectors to hold data at the quadrature points { this->xyz.resize(n_qp); normals.resize(n_qp); this->JxW.resize(n_qp); } // If we have no quadrature points, there's nothing else to do if (!n_qp) break; // We need to look back at the full edge to figure out the normal // vector const Elem *elem = side->parent(); libmesh_assert (elem); if (side->node(0) == elem->node(0)) normals[0] = Point(-1.); else { libmesh_assert (side->node(0) == elem->node(1)); normals[0] = Point(1.); } // Calculate x at the point libmesh_assert (this->psi_map.size() == 1); // In the unlikely event we have multiple quadrature // points, they'll be in the same place for (unsigned int p=0; p<n_qp; p++) { this->xyz[p].zero(); this->xyz[p].add_scaled (side->point(0), this->psi_map[0][p]); normals[p] = normals[0]; this->JxW[p] = 1.0*qw[p]; } // done computing the map break; } case 2: { // A 2D finite element living in either 2D or 3D space. // This means the boundary is a 1D finite element, i.e. // and EDGE2 or EDGE3. // Resize the vectors to hold data at the quadrature points { this->xyz.resize(n_qp); this->dxyzdxi_map.resize(n_qp); this->d2xyzdxi2_map.resize(n_qp); this->tangents.resize(n_qp); this->normals.resize(n_qp); this->curvatures.resize(n_qp); this->JxW.resize(n_qp); } // Clear the entities that will be summed // Compute the tangent & normal at the quadrature point for (unsigned int p=0; p<n_qp; p++) { this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D this->xyz[p].zero(); this->dxyzdxi_map[p].zero(); this->d2xyzdxi2_map[p].zero(); } // compute x, dxdxi at the quadrature points for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes { const Point& side_point = side->point(i); for (unsigned int p=0; p<n_qp; p++) // for each quadrature point... { this->xyz[p].add_scaled (side_point, this->psi_map[i][p]); this->dxyzdxi_map[p].add_scaled (side_point, this->dpsidxi_map[i][p]); this->d2xyzdxi2_map[p].add_scaled(side_point, this->d2psidxi2_map[i][p]); } } // Compute the tangent & normal at the quadrature point for (unsigned int p=0; p<n_qp; p++) { // The first tangent comes from just the edge's Jacobian this->tangents[p][0] = this->dxyzdxi_map[p].unit(); #if LIBMESH_DIM == 2 // For a 2D element living in 2D, the normal is given directly // from the entries in the edge Jacobian. this->normals[p] = (Point(this->dxyzdxi_map[p](1), -this->dxyzdxi_map[p](0), 0.)).unit(); #elif LIBMESH_DIM == 3 // For a 2D element living in 3D, there is a second tangent. // For the second tangent, we need to refer to the full // element's (not just the edge's) Jacobian. const Elem *elem = side->parent(); libmesh_assert (elem != NULL); // Inverse map xyz[p] to a reference point on the parent... Point reference_point = FE<2,LAGRANGE>::inverse_map(elem, this->xyz[p]); // Get dxyz/dxi and dxyz/deta from the parent map. Point dx_dxi = FE<2,LAGRANGE>::map_xi (elem, reference_point); Point dx_deta = FE<2,LAGRANGE>::map_eta(elem, reference_point); // The second tangent vector is formed by crossing these vectors. tangents[p][1] = dx_dxi.cross(dx_deta).unit(); // Finally, the normal in this case is given by crossing these // two tangents. normals[p] = tangents[p][0].cross(tangents[p][1]).unit(); #endif // The curvature is computed via the familiar Frenet formula: // curvature = [d^2(x) / d (xi)^2] dot [normal] // For a reference, see: // F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill, p. 310 // // Note: The sign convention here is different from the // 3D case. Concave-upward curves (smiles) have a positive // curvature. Concave-downward curves (frowns) have a // negative curvature. Be sure to take that into account! const Real numerator = this->d2xyzdxi2_map[p] * this->normals[p]; const Real denominator = this->dxyzdxi_map[p].size_sq(); libmesh_assert (denominator != 0); curvatures[p] = numerator / denominator; } // compute the jacobian at the quadrature points for (unsigned int p=0; p<n_qp; p++) { const Real jac = this->dxyzdxi_map[p].size(); libmesh_assert (jac > 0.); this->JxW[p] = jac*qw[p]; } // done computing the map break; } case 3: { // A 3D finite element living in 3D space. // Resize the vectors to hold data at the quadrature points { this->xyz.resize(n_qp); this->dxyzdxi_map.resize(n_qp); this->dxyzdeta_map.resize(n_qp); this->d2xyzdxi2_map.resize(n_qp); this->d2xyzdxideta_map.resize(n_qp); this->d2xyzdeta2_map.resize(n_qp); this->tangents.resize(n_qp); this->normals.resize(n_qp); this->curvatures.resize(n_qp); this->JxW.resize(n_qp); } // Clear the entities that will be summed for (unsigned int p=0; p<n_qp; p++) { this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D this->xyz[p].zero(); this->dxyzdxi_map[p].zero(); this->dxyzdeta_map[p].zero(); this->d2xyzdxi2_map[p].zero(); this->d2xyzdxideta_map[p].zero(); this->d2xyzdeta2_map[p].zero(); } // compute x, dxdxi at the quadrature points for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes { const Point& side_point = side->point(i); for (unsigned int p=0; p<n_qp; p++) // for each quadrature point... { this->xyz[p].add_scaled (side_point, this->psi_map[i][p]); this->dxyzdxi_map[p].add_scaled (side_point, this->dpsidxi_map[i][p]); this->dxyzdeta_map[p].add_scaled(side_point, this->dpsideta_map[i][p]); this->d2xyzdxi2_map[p].add_scaled (side_point, this->d2psidxi2_map[i][p]); this->d2xyzdxideta_map[p].add_scaled(side_point, this->d2psidxideta_map[i][p]); this->d2xyzdeta2_map[p].add_scaled (side_point, this->d2psideta2_map[i][p]); } } // Compute the tangents, normal, and curvature at the quadrature point for (unsigned int p=0; p<n_qp; p++) { const Point n = this->dxyzdxi_map[p].cross(this->dxyzdeta_map[p]); this->normals[p] = n.unit(); this->tangents[p][0] = this->dxyzdxi_map[p].unit(); this->tangents[p][1] = n.cross(this->dxyzdxi_map[p]).unit(); // Compute curvature using the typical nomenclature // of the first and second fundamental forms. // For reference, see: // 1) http://mathworld.wolfram.com/MeanCurvature.html // (note -- they are using inward normal) // 2) F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill const Real L = -this->d2xyzdxi2_map[p] * this->normals[p]; const Real M = -this->d2xyzdxideta_map[p] * this->normals[p]; const Real N = -this->d2xyzdeta2_map[p] * this->normals[p]; const Real E = this->dxyzdxi_map[p].size_sq(); const Real F = this->dxyzdxi_map[p] * this->dxyzdeta_map[p]; const Real G = this->dxyzdeta_map[p].size_sq(); const Real numerator = E*N -2.*F*M + G*L; const Real denominator = E*G - F*F; libmesh_assert (denominator != 0.); curvatures[p] = 0.5*numerator/denominator; } // compute the jacobian at the quadrature points, see // http://sp81.msi.umn.edu:999/fluent/fidap/help/theory/thtoc.htm for (unsigned int p=0; p<n_qp; p++) { const Real g11 = (dxdxi_map(p)*dxdxi_map(p) + dydxi_map(p)*dydxi_map(p) + dzdxi_map(p)*dzdxi_map(p)); const Real g12 = (dxdxi_map(p)*dxdeta_map(p) + dydxi_map(p)*dydeta_map(p) + dzdxi_map(p)*dzdeta_map(p)); const Real g21 = g12; const Real g22 = (dxdeta_map(p)*dxdeta_map(p) + dydeta_map(p)*dydeta_map(p) + dzdeta_map(p)*dzdeta_map(p)); const Real jac = std::sqrt(g11*g22 - g12*g21); libmesh_assert (jac > 0.); this->JxW[p] = jac*qw[p]; } // done computing the map break; } default: libmesh_error(); } STOP_LOG("compute_face_map()", "FEMap"); }
void FEMap::compute_single_point_map(const unsigned int dim, const std::vector<Real>& qw, const Elem* elem, unsigned int p) { libmesh_assert(elem); switch (dim) { //-------------------------------------------------------------------- // 0D case 0: { xyz[p] = elem->point(0); jac[p] = 1.0; JxW[p] = qw[p]; break; } //-------------------------------------------------------------------- // 1D case 1: { // Clear the entities that will be summed xyz[p].zero(); dxyzdxi_map[p].zero(); #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES d2xyzdxi2_map[p].zero(); #endif // compute x, dx, d2x at the quadrature point for (unsigned int i=0; i<phi_map.size(); i++) // sum over the nodes { // Reference to the point, helps eliminate // exessive temporaries in the inner loop const Point& elem_point = elem->point(i); xyz[p].add_scaled (elem_point, phi_map[i][p] ); dxyzdxi_map[p].add_scaled (elem_point, dphidxi_map[i][p]); #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES d2xyzdxi2_map[p].add_scaled(elem_point, d2phidxi2_map[i][p]); #endif } // Compute the jacobian // // 1D elements can live in 2D or 3D space. // The transformation matrix from local->global // coordinates is // // T = | dx/dxi | // | dy/dxi | // | dz/dxi | // // The generalized determinant of T (from the // so-called "normal" eqns.) is // jac = "det(T)" = sqrt(det(T'T)) // // where T'= transpose of T, so // // jac = sqrt( (dx/dxi)^2 + (dy/dxi)^2 + (dz/dxi)^2 ) jac[p] = dxyzdxi_map[p].size(); if (jac[p] <= 0.) { libMesh::err << "ERROR: negative Jacobian: " << jac[p] << " in element " << elem->id() << std::endl; libmesh_error(); } // The inverse Jacobian entries also come from the // generalized inverse of T (see also the 2D element // living in 3D code). const Real jacm2 = 1./jac[p]/jac[p]; dxidx_map[p] = jacm2*dxdxi_map(p); #if LIBMESH_DIM > 1 dxidy_map[p] = jacm2*dydxi_map(p); #endif #if LIBMESH_DIM > 2 dxidz_map[p] = jacm2*dzdxi_map(p); #endif JxW[p] = jac[p]*qw[p]; // done computing the map break; } //-------------------------------------------------------------------- // 2D case 2: { //------------------------------------------------------------------ // Compute the (x,y) values at the quadrature points, // the Jacobian at the quadrature points xyz[p].zero(); dxyzdxi_map[p].zero(); dxyzdeta_map[p].zero(); #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES d2xyzdxi2_map[p].zero(); d2xyzdxideta_map[p].zero(); d2xyzdeta2_map[p].zero(); #endif // compute (x,y) at the quadrature points, derivatives once for (unsigned int i=0; i<phi_map.size(); i++) // sum over the nodes { // Reference to the point, helps eliminate // exessive temporaries in the inner loop const Point& elem_point = elem->point(i); xyz[p].add_scaled (elem_point, phi_map[i][p] ); dxyzdxi_map[p].add_scaled (elem_point, dphidxi_map[i][p] ); dxyzdeta_map[p].add_scaled (elem_point, dphideta_map[i][p]); #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES d2xyzdxi2_map[p].add_scaled (elem_point, d2phidxi2_map[i][p]); d2xyzdxideta_map[p].add_scaled (elem_point, d2phidxideta_map[i][p]); d2xyzdeta2_map[p].add_scaled (elem_point, d2phideta2_map[i][p]); #endif } // compute the jacobian once const Real dx_dxi = dxdxi_map(p), dx_deta = dxdeta_map(p), dy_dxi = dydxi_map(p), dy_deta = dydeta_map(p); #if LIBMESH_DIM == 2 // Compute the Jacobian. This assumes the 2D face // lives in 2D space // // Symbolically, the matrix determinant is // // | dx/dxi dx/deta | // jac = | dy/dxi dy/deta | // // jac = dx/dxi*dy/deta - dx/deta*dy/dxi jac[p] = (dx_dxi*dy_deta - dx_deta*dy_dxi); if (jac[p] <= 0.) { libMesh::err << "ERROR: negative Jacobian: " << jac[p] << " in element " << elem->id() << std::endl; libmesh_error(); } JxW[p] = jac[p]*qw[p]; // Compute the shape function derivatives wrt x,y at the // quadrature points const Real inv_jac = 1./jac[p]; dxidx_map[p] = dy_deta*inv_jac; //dxi/dx = (1/J)*dy/deta dxidy_map[p] = -dx_deta*inv_jac; //dxi/dy = -(1/J)*dx/deta detadx_map[p] = -dy_dxi* inv_jac; //deta/dx = -(1/J)*dy/dxi detady_map[p] = dx_dxi* inv_jac; //deta/dy = (1/J)*dx/dxi dxidz_map[p] = detadz_map[p] = 0.; #else const Real dz_dxi = dzdxi_map(p), dz_deta = dzdeta_map(p); // Compute the Jacobian. This assumes a 2D face in // 3D space. // // The transformation matrix T from local to global // coordinates is // // | dx/dxi dx/deta | // T = | dy/dxi dy/deta | // | dz/dxi dz/deta | // note det(T' T) = det(T')det(T) = det(T)det(T) // so det(T) = std::sqrt(det(T' T)) // //---------------------------------------------- // Notes: // // dX = R dXi -> R'dX = R'R dXi // (R^-1)dX = dXi [(R'R)^-1 R']dX = dXi // // so R^-1 = (R'R)^-1 R' // // and R^-1 R = (R'R)^-1 R'R = I. // const Real g11 = (dx_dxi*dx_dxi + dy_dxi*dy_dxi + dz_dxi*dz_dxi); const Real g12 = (dx_dxi*dx_deta + dy_dxi*dy_deta + dz_dxi*dz_deta); const Real g21 = g12; const Real g22 = (dx_deta*dx_deta + dy_deta*dy_deta + dz_deta*dz_deta); const Real det = (g11*g22 - g12*g21); if (det <= 0.) { libMesh::err << "ERROR: negative Jacobian! " << " in element " << elem->id() << std::endl; libmesh_error(); } const Real inv_det = 1./det; jac[p] = std::sqrt(det); JxW[p] = jac[p]*qw[p]; const Real g11inv = g22*inv_det; const Real g12inv = -g12*inv_det; const Real g21inv = -g21*inv_det; const Real g22inv = g11*inv_det; dxidx_map[p] = g11inv*dx_dxi + g12inv*dx_deta; dxidy_map[p] = g11inv*dy_dxi + g12inv*dy_deta; dxidz_map[p] = g11inv*dz_dxi + g12inv*dz_deta; detadx_map[p] = g21inv*dx_dxi + g22inv*dx_deta; detady_map[p] = g21inv*dy_dxi + g22inv*dy_deta; detadz_map[p] = g21inv*dz_dxi + g22inv*dz_deta; #endif // done computing the map break; } //-------------------------------------------------------------------- // 3D case 3: { //------------------------------------------------------------------ // Compute the (x,y,z) values at the quadrature points, // the Jacobian at the quadrature point // Clear the entities that will be summed xyz[p].zero (); dxyzdxi_map[p].zero (); dxyzdeta_map[p].zero (); dxyzdzeta_map[p].zero (); #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES d2xyzdxi2_map[p].zero(); d2xyzdxideta_map[p].zero(); d2xyzdxidzeta_map[p].zero(); d2xyzdeta2_map[p].zero(); d2xyzdetadzeta_map[p].zero(); d2xyzdzeta2_map[p].zero(); #endif // compute (x,y,z) at the quadrature points, // dxdxi, dydxi, dzdxi, // dxdeta, dydeta, dzdeta, // dxdzeta, dydzeta, dzdzeta all once for (unsigned int i=0; i<phi_map.size(); i++) // sum over the nodes { // Reference to the point, helps eliminate // exessive temporaries in the inner loop const Point& elem_point = elem->point(i); xyz[p].add_scaled (elem_point, phi_map[i][p] ); dxyzdxi_map[p].add_scaled (elem_point, dphidxi_map[i][p] ); dxyzdeta_map[p].add_scaled (elem_point, dphideta_map[i][p] ); dxyzdzeta_map[p].add_scaled (elem_point, dphidzeta_map[i][p]); #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES d2xyzdxi2_map[p].add_scaled (elem_point, d2phidxi2_map[i][p]); d2xyzdxideta_map[p].add_scaled (elem_point, d2phidxideta_map[i][p]); d2xyzdxidzeta_map[p].add_scaled (elem_point, d2phidxidzeta_map[i][p]); d2xyzdeta2_map[p].add_scaled (elem_point, d2phideta2_map[i][p]); d2xyzdetadzeta_map[p].add_scaled (elem_point, d2phidetadzeta_map[i][p]); d2xyzdzeta2_map[p].add_scaled (elem_point, d2phidzeta2_map[i][p]); #endif } // compute the jacobian const Real dx_dxi = dxdxi_map(p), dy_dxi = dydxi_map(p), dz_dxi = dzdxi_map(p), dx_deta = dxdeta_map(p), dy_deta = dydeta_map(p), dz_deta = dzdeta_map(p), dx_dzeta = dxdzeta_map(p), dy_dzeta = dydzeta_map(p), dz_dzeta = dzdzeta_map(p); // Symbolically, the matrix determinant is // // | dx/dxi dy/dxi dz/dxi | // jac = | dx/deta dy/deta dz/deta | // | dx/dzeta dy/dzeta dz/dzeta | // // jac = dx/dxi*(dy/deta*dz/dzeta - dz/deta*dy/dzeta) + // dy/dxi*(dz/deta*dx/dzeta - dx/deta*dz/dzeta) + // dz/dxi*(dx/deta*dy/dzeta - dy/deta*dx/dzeta) jac[p] = (dx_dxi*(dy_deta*dz_dzeta - dz_deta*dy_dzeta) + dy_dxi*(dz_deta*dx_dzeta - dx_deta*dz_dzeta) + dz_dxi*(dx_deta*dy_dzeta - dy_deta*dx_dzeta)); if (jac[p] <= 0.) { libMesh::err << "ERROR: negative Jacobian: " << jac[p] << " in element " << elem->id() << std::endl; libmesh_error(); } JxW[p] = jac[p]*qw[p]; // Compute the shape function derivatives wrt x,y at the // quadrature points const Real inv_jac = 1./jac[p]; dxidx_map[p] = (dy_deta*dz_dzeta - dz_deta*dy_dzeta)*inv_jac; dxidy_map[p] = (dz_deta*dx_dzeta - dx_deta*dz_dzeta)*inv_jac; dxidz_map[p] = (dx_deta*dy_dzeta - dy_deta*dx_dzeta)*inv_jac; detadx_map[p] = (dz_dxi*dy_dzeta - dy_dxi*dz_dzeta )*inv_jac; detady_map[p] = (dx_dxi*dz_dzeta - dz_dxi*dx_dzeta )*inv_jac; detadz_map[p] = (dy_dxi*dx_dzeta - dx_dxi*dy_dzeta )*inv_jac; dzetadx_map[p] = (dy_dxi*dz_deta - dz_dxi*dy_deta )*inv_jac; dzetady_map[p] = (dz_dxi*dx_deta - dx_dxi*dz_deta )*inv_jac; dzetadz_map[p] = (dx_dxi*dy_deta - dy_dxi*dx_deta )*inv_jac; // done computing the map break; } default: libmesh_error(); } }