$extrastylesheet
fe_xyz_map.C
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 #include "libmesh/fe_xyz_map.h"
00019 
00020 namespace libMesh {
00021 
00022 void FEXYZMap::compute_face_map(int dim, const std::vector<Real>& qw, const Elem* side)
00023 {
00024   libmesh_assert(side);
00025 
00026   START_LOG("compute_face_map()", "FEXYZMap");
00027 
00028   // The number of quadrature points.
00029   const unsigned int n_qp = cast_int<unsigned int>(qw.size());
00030 
00031   switch(dim)
00032     {
00033     case 2:
00034       {
00035 
00036         // Resize the vectors to hold data at the quadrature points
00037         {
00038           this->xyz.resize(n_qp);
00039           this->dxyzdxi_map.resize(n_qp);
00040           this->d2xyzdxi2_map.resize(n_qp);
00041           this->tangents.resize(n_qp);
00042           this->normals.resize(n_qp);
00043           this->curvatures.resize(n_qp);
00044 
00045           this->JxW.resize(n_qp);
00046         }
00047 
00048         // Clear the entities that will be summed
00049         // Compute the tangent & normal at the quadrature point
00050         for (unsigned int p=0; p<n_qp; p++)
00051           {
00052             this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D
00053             this->xyz[p].zero();
00054             this->dxyzdxi_map[p].zero();
00055             this->d2xyzdxi2_map[p].zero();
00056           }
00057 
00058         // compute x, dxdxi at the quadrature points
00059         for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes
00060           {
00061             const Point& side_point = side->point(i);
00062 
00063             for (unsigned int p=0; p<n_qp; p++) // for each quadrature point...
00064               {
00065                 this->xyz[p].add_scaled          (side_point, this->psi_map[i][p]);
00066                 this->dxyzdxi_map[p].add_scaled  (side_point, this->dpsidxi_map[i][p]);
00067                 this->d2xyzdxi2_map[p].add_scaled(side_point, this->d2psidxi2_map[i][p]);
00068               }
00069           }
00070 
00071         // Compute the tangent & normal at the quadrature point
00072         for (unsigned int p=0; p<n_qp; p++)
00073           {
00074             const Point n(this->dxyzdxi_map[p](1), -this->dxyzdxi_map[p](0), 0.);
00075 
00076             this->normals[p]     = n.unit();
00077             this->tangents[p][0] = this->dxyzdxi_map[p].unit();
00078 #if LIBMESH_DIM == 3  // Only good in 3D space
00079             this->tangents[p][1] = this->dxyzdxi_map[p].cross(n).unit();
00080 #endif
00081             // The curvature is computed via the familiar Frenet formula:
00082             // curvature = [d^2(x) / d (xi)^2] dot [normal]
00083             // For a reference, see:
00084             // F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill, p. 310
00085             //
00086             // Note: The sign convention here is different from the
00087             // 3D case.  Concave-upward curves (smiles) have a positive
00088             // curvature.  Concave-downward curves (frowns) have a
00089             // negative curvature.  Be sure to take that into account!
00090             const Real numerator   = this->d2xyzdxi2_map[p] * this->normals[p];
00091             const Real denominator = this->dxyzdxi_map[p].size_sq();
00092             libmesh_assert_not_equal_to (denominator, 0);
00093             this->curvatures[p] = numerator / denominator;
00094           }
00095 
00096         // compute the jacobian at the quadrature points
00097         for (unsigned int p=0; p<n_qp; p++)
00098           {
00099             const Real the_jac = std::sqrt(this->dxdxi_map(p)*this->dxdxi_map(p) +
00100                                            this->dydxi_map(p)*this->dydxi_map(p));
00101 
00102             libmesh_assert_greater (the_jac, 0.);
00103 
00104             this->JxW[p] = the_jac*qw[p];
00105           }
00106 
00107         break;
00108       }
00109 
00110     case 3:
00111       {
00112         // Resize the vectors to hold data at the quadrature points
00113         {
00114           this->xyz.resize(n_qp);
00115           this->dxyzdxi_map.resize(n_qp);
00116           this->dxyzdeta_map.resize(n_qp);
00117           this->d2xyzdxi2_map.resize(n_qp);
00118           this->d2xyzdxideta_map.resize(n_qp);
00119           this->d2xyzdeta2_map.resize(n_qp);
00120           this->tangents.resize(n_qp);
00121           this->normals.resize(n_qp);
00122           this->curvatures.resize(n_qp);
00123 
00124           this->JxW.resize(n_qp);
00125         }
00126 
00127         // Clear the entities that will be summed
00128         for (unsigned int p=0; p<n_qp; p++)
00129           {
00130             this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D
00131             this->xyz[p].zero();
00132             this->dxyzdxi_map[p].zero();
00133             this->dxyzdeta_map[p].zero();
00134             this->d2xyzdxi2_map[p].zero();
00135             this->d2xyzdxideta_map[p].zero();
00136             this->d2xyzdeta2_map[p].zero();
00137           }
00138 
00139         // compute x, dxdxi at the quadrature points
00140         for (unsigned int i=0; i<this->psi_map.size(); i++) // sum over the nodes
00141           {
00142             const Point& side_point = side->point(i);
00143 
00144             for (unsigned int p=0; p<n_qp; p++) // for each quadrature point...
00145               {
00146                 this->xyz[p].add_scaled             (side_point, this->psi_map[i][p]);
00147                 this->dxyzdxi_map[p].add_scaled     (side_point, this->dpsidxi_map[i][p]);
00148                 this->dxyzdeta_map[p].add_scaled    (side_point, this->dpsideta_map[i][p]);
00149                 this->d2xyzdxi2_map[p].add_scaled   (side_point, this->d2psidxi2_map[i][p]);
00150                 this->d2xyzdxideta_map[p].add_scaled(side_point, this->d2psidxideta_map[i][p]);
00151                 this->d2xyzdeta2_map[p].add_scaled  (side_point, this->d2psideta2_map[i][p]);
00152               }
00153           }
00154 
00155         // Compute the tangents, normal, and curvature at the quadrature point
00156         for (unsigned int p=0; p<n_qp; p++)
00157           {
00158             const Point n  = this->dxyzdxi_map[p].cross(this->dxyzdeta_map[p]);
00159             this->normals[p]     = n.unit();
00160             this->tangents[p][0] = this->dxyzdxi_map[p].unit();
00161             this->tangents[p][1] = n.cross(this->dxyzdxi_map[p]).unit();
00162 
00163             // Compute curvature using the typical nomenclature
00164             // of the first and second fundamental forms.
00165             // For reference, see:
00166             // 1) http://mathworld.wolfram.com/MeanCurvature.html
00167             //    (note -- they are using inward normal)
00168             // 2) F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill
00169             const Real L  = -this->d2xyzdxi2_map[p]    * this->normals[p];
00170             const Real M  = -this->d2xyzdxideta_map[p] * this->normals[p];
00171             const Real N  = -this->d2xyzdeta2_map[p]   * this->normals[p];
00172             const Real E  =  this->dxyzdxi_map[p].size_sq();
00173             const Real F  =  this->dxyzdxi_map[p]      * this->dxyzdeta_map[p];
00174             const Real G  =  this->dxyzdeta_map[p].size_sq();
00175 
00176             const Real numerator   = E*N -2.*F*M + G*L;
00177             const Real denominator = E*G - F*F;
00178             libmesh_assert_not_equal_to (denominator, 0.);
00179             this->curvatures[p] = 0.5*numerator/denominator;
00180           }
00181 
00182         // compute the jacobian at the quadrature points, see
00183         // http://sp81.msi.umn.edu:999/fluent/fidap/help/theory/thtoc.htm
00184         for (unsigned int p=0; p<n_qp; p++)
00185           {
00186             const Real g11 = (this->dxdxi_map(p)*this->dxdxi_map(p) +
00187                               this->dydxi_map(p)*this->dydxi_map(p) +
00188                               this->dzdxi_map(p)*this->dzdxi_map(p));
00189 
00190             const Real g12 = (this->dxdxi_map(p)*this->dxdeta_map(p) +
00191                               this->dydxi_map(p)*this->dydeta_map(p) +
00192                               this->dzdxi_map(p)*this->dzdeta_map(p));
00193 
00194             const Real g21 = g12;
00195 
00196             const Real g22 = (this->dxdeta_map(p)*this->dxdeta_map(p) +
00197                               this->dydeta_map(p)*this->dydeta_map(p) +
00198                               this->dzdeta_map(p)*this->dzdeta_map(p));
00199 
00200 
00201             const Real the_jac = std::sqrt(g11*g22 - g12*g21);
00202 
00203             libmesh_assert_greater (the_jac, 0.);
00204 
00205             this->JxW[p] = the_jac*qw[p];
00206           }
00207 
00208         break;
00209       }
00210     default:
00211       libmesh_error_msg("Invalid dim = " << dim);
00212 
00213     } // switch(dim)
00214 
00215   STOP_LOG("compute_face_map()", "FEXYZMap");
00216 
00217   return;
00218 }
00219 
00220 } // namespace libMesh