void HCurlFETransformation<OutputShape>::map_curl( const unsigned int dim, const Elem* const, const std::vector<Point>&, const FEGenericBase<OutputShape>& fe, std::vector<std::vector<OutputShape> >& curl_phi ) const { switch(dim) { // These element transformations only make sense in 2D and 3D case 0: case 1: { libmesh_error(); } case 2: { const std::vector<std::vector<OutputShape> >& dphi_dxi = fe.get_dphidxi(); const std::vector<std::vector<OutputShape> >& dphi_deta = fe.get_dphideta(); const std::vector<Real>& J = fe.get_fe_map().get_jacobian(); // FIXME: I don't think this is valid for 2D elements in 3D space /* In 2D: curl(phi) = J^{-1} * curl(\hat{phi}) */ for (unsigned int i=0; i<curl_phi.size(); i++) for (unsigned int p=0; p<curl_phi[i].size(); p++) { curl_phi[i][p].slice(0) = curl_phi[i][p].slice(1) = 0.0; curl_phi[i][p].slice(2) = ( dphi_dxi[i][p].slice(1) - dphi_deta[i][p].slice(0) )/J[p]; } break; } case 3: { const std::vector<std::vector<OutputShape> >& dphi_dxi = fe.get_dphidxi(); const std::vector<std::vector<OutputShape> >& dphi_deta = fe.get_dphideta(); const std::vector<std::vector<OutputShape> >& dphi_dzeta = fe.get_dphidzeta(); const std::vector<RealGradient>& dxyz_dxi = fe.get_fe_map().get_dxyzdxi(); const std::vector<RealGradient>& dxyz_deta = fe.get_fe_map().get_dxyzdeta(); const std::vector<RealGradient>& dxyz_dzeta = fe.get_fe_map().get_dxyzdzeta(); const std::vector<Real>& J = fe.get_fe_map().get_jacobian(); for (unsigned int i=0; i<curl_phi.size(); i++) for (unsigned int p=0; p<curl_phi[i].size(); p++) { Real dx_dxi = dxyz_dxi[p](0); Real dx_deta = dxyz_deta[p](0); Real dx_dzeta = dxyz_dzeta[p](0); Real dy_dxi = dxyz_dxi[p](1); Real dy_deta = dxyz_deta[p](1); Real dy_dzeta = dxyz_dzeta[p](1); Real dz_dxi = dxyz_dxi[p](2); Real dz_deta = dxyz_deta[p](2); Real dz_dzeta = dxyz_dzeta[p](2); const Real inv_jac = 1.0/J[p]; /* In 3D: curl(phi) = J^{-1} dx/dxi * curl(\hat{phi}) dx/dxi = [ dx/dxi dx/deta dx/dzeta dy/dxi dy/deta dy/dzeta dz/dxi dz/deta dz/dzeta ] curl(u) = [ du_z/deta - du_y/dzeta du_x/dzeta - du_z/dxi du_y/dxi - du_x/deta ] */ curl_phi[i][p].slice(0) = inv_jac*( dx_dxi*( dphi_deta[i][p].slice(2) - dphi_dzeta[i][p].slice(1) ) + dx_deta*( dphi_dzeta[i][p].slice(0) - dphi_dxi[i][p].slice(2) ) + dx_dzeta*( dphi_dxi[i][p].slice(1) - dphi_deta[i][p].slice(0) ) ); curl_phi[i][p].slice(1) = inv_jac*( dy_dxi*( dphi_deta[i][p].slice(2) - dphi_dzeta[i][p].slice(1) ) + dy_deta*( dphi_dzeta[i][p].slice(0)- dphi_dxi[i][p].slice(2) ) + dy_dzeta*( dphi_dxi[i][p].slice(1) - dphi_deta[i][p].slice(0) ) ); curl_phi[i][p].slice(2) = inv_jac*( dz_dxi*( dphi_deta[i][p].slice(2) - dphi_dzeta[i][p].slice(1) ) + dz_deta*( dphi_dzeta[i][p].slice(0) - dphi_dxi[i][p].slice(2) ) + dz_dzeta*( dphi_dxi[i][p].slice(1) - dphi_deta[i][p].slice(0) ) ); } break; } default: libmesh_error(); } // switch(dim) return; }
void H1FETransformation<OutputShape>::map_dphi( const unsigned int dim, const Elem* const, const std::vector<Point>&, const FEGenericBase<OutputShape>& fe, std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> >& dphi, std::vector<std::vector<OutputShape> >& dphidx, std::vector<std::vector<OutputShape> >& dphidy, std::vector<std::vector<OutputShape> >& dphidz ) const { switch(dim) { case 0: // No derivatives in 0D { for (unsigned int i=0; i<dphi.size(); i++) for (unsigned int p=0; p<dphi[i].size(); p++) { dphi[i][p] = 0.; } break; } case 1: { const std::vector<std::vector<OutputShape> >& dphidxi = fe.get_dphidxi(); const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx(); #if LIBMESH_DIM>1 const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy(); #endif #if LIBMESH_DIM>2 const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz(); #endif for (unsigned int i=0; i<dphi.size(); i++) for (unsigned int p=0; p<dphi[i].size(); p++) { // dphi/dx = (dphi/dxi)*(dxi/dx) dphi[i][p].slice(0) = dphidx[i][p] = dphidxi[i][p]*dxidx_map[p]; #if LIBMESH_DIM>1 dphi[i][p].slice(1) = dphidy[i][p] = dphidxi[i][p]*dxidy_map[p]; #endif #if LIBMESH_DIM>2 dphi[i][p].slice(2) = dphidz[i][p] = dphidxi[i][p]*dxidz_map[p]; #endif } break; } case 2: { const std::vector<std::vector<OutputShape> >& dphidxi = fe.get_dphidxi(); const std::vector<std::vector<OutputShape> >& dphideta = fe.get_dphideta(); const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx(); const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy(); #if LIBMESH_DIM > 2 const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz(); #endif const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx(); const std::vector<Real>& detady_map = fe.get_fe_map().get_detady(); #if LIBMESH_DIM > 2 const std::vector<Real>& detadz_map = fe.get_fe_map().get_detadz(); #endif for (unsigned int i=0; i<dphi.size(); i++) for (unsigned int p=0; p<dphi[i].size(); p++) { // dphi/dx = (dphi/dxi)*(dxi/dx) + (dphi/deta)*(deta/dx) dphi[i][p].slice(0) = dphidx[i][p] = (dphidxi[i][p]*dxidx_map[p] + dphideta[i][p]*detadx_map[p]); // dphi/dy = (dphi/dxi)*(dxi/dy) + (dphi/deta)*(deta/dy) dphi[i][p].slice(1) = dphidy[i][p] = (dphidxi[i][p]*dxidy_map[p] + dphideta[i][p]*detady_map[p]); #if LIBMESH_DIM > 2 // dphi/dz = (dphi/dxi)*(dxi/dz) + (dphi/deta)*(deta/dz) dphi[i][p].slice(2) = dphidz[i][p] = (dphidxi[i][p]*dxidz_map[p] + dphideta[i][p]*detadz_map[p]); #endif } break; } case 3: { const std::vector<std::vector<OutputShape> >& dphidxi = fe.get_dphidxi(); const std::vector<std::vector<OutputShape> >& dphideta = fe.get_dphideta(); const std::vector<std::vector<OutputShape> >& dphidzeta = fe.get_dphidzeta(); const std::vector<Real>& dxidx_map = fe.get_fe_map().get_dxidx(); const std::vector<Real>& dxidy_map = fe.get_fe_map().get_dxidy(); const std::vector<Real>& dxidz_map = fe.get_fe_map().get_dxidz(); const std::vector<Real>& detadx_map = fe.get_fe_map().get_detadx(); const std::vector<Real>& detady_map = fe.get_fe_map().get_detady(); const std::vector<Real>& detadz_map = fe.get_fe_map().get_detadz(); const std::vector<Real>& dzetadx_map = fe.get_fe_map().get_dzetadx(); const std::vector<Real>& dzetady_map = fe.get_fe_map().get_dzetady(); const std::vector<Real>& dzetadz_map = fe.get_fe_map().get_dzetadz(); for (unsigned int i=0; i<dphi.size(); i++) for (unsigned int p=0; p<dphi[i].size(); p++) { // dphi/dx = (dphi/dxi)*(dxi/dx) + (dphi/deta)*(deta/dx) + (dphi/dzeta)*(dzeta/dx); dphi[i][p].slice(0) = dphidx[i][p] = (dphidxi[i][p]*dxidx_map[p] + dphideta[i][p]*detadx_map[p] + dphidzeta[i][p]*dzetadx_map[p]); // dphi/dy = (dphi/dxi)*(dxi/dy) + (dphi/deta)*(deta/dy) + (dphi/dzeta)*(dzeta/dy); dphi[i][p].slice(1) = dphidy[i][p] = (dphidxi[i][p]*dxidy_map[p] + dphideta[i][p]*detady_map[p] + dphidzeta[i][p]*dzetady_map[p]); // dphi/dz = (dphi/dxi)*(dxi/dz) + (dphi/deta)*(deta/dz) + (dphi/dzeta)*(dzeta/dz); dphi[i][p].slice(2) = dphidz[i][p] = (dphidxi[i][p]*dxidz_map[p] + dphideta[i][p]*detadz_map[p] + dphidzeta[i][p]*dzetadz_map[p]); } break; } default: libmesh_error(); } // switch(dim) return; }