libMesh
fe_xyz_map.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 #include "libmesh/fe_xyz_map.h"
19 #include "libmesh/elem.h"
20 #include "libmesh/libmesh_logging.h"
21 
22 namespace libMesh
23 {
24 
25 void FEXYZMap::compute_face_map(int dim, const std::vector<Real> & qw, const Elem * side)
26 {
27  libmesh_assert(side);
28 
29  LOG_SCOPE("compute_face_map()", "FEXYZMap");
30 
31  // The number of quadrature points.
32  const unsigned int n_qp = cast_int<unsigned int>(qw.size());
33 
34  switch(dim)
35  {
36  case 2:
37  {
38 
39  // Resize the vectors to hold data at the quadrature points
40  {
41  this->xyz.resize(n_qp);
42  this->dxyzdxi_map.resize(n_qp);
43 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
44  this->d2xyzdxi2_map.resize(n_qp);
45  this->curvatures.resize(n_qp);
46 #endif
47  this->tangents.resize(n_qp);
48  this->normals.resize(n_qp);
49 
50  this->JxW.resize(n_qp);
51  }
52 
53  // Clear the entities that will be summed
54  // Compute the tangent & normal at the quadrature point
55  for (unsigned int p=0; p<n_qp; p++)
56  {
57  this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D
58  this->xyz[p].zero();
59  this->dxyzdxi_map[p].zero();
60 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
61  this->d2xyzdxi2_map[p].zero();
62 #endif
63  }
64 
65  // compute x, dxdxi at the quadrature points
66  for (unsigned int i=0,
67  n_psi_map = cast_int<unsigned int>(this->psi_map.size());
68  i != n_psi_map; i++) // sum over the nodes
69  {
70  const Point & side_point = side->point(i);
71 
72  for (unsigned int p=0; p<n_qp; p++) // for each quadrature point...
73  {
74  this->xyz[p].add_scaled (side_point, this->psi_map[i][p]);
75  this->dxyzdxi_map[p].add_scaled (side_point, this->dpsidxi_map[i][p]);
76 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
77  this->d2xyzdxi2_map[p].add_scaled(side_point, this->d2psidxi2_map[i][p]);
78 #endif
79  }
80  }
81 
82  // Compute the tangent & normal at the quadrature point
83  for (unsigned int p=0; p<n_qp; p++)
84  {
85  const Point n(this->dxyzdxi_map[p](1), -this->dxyzdxi_map[p](0), 0.);
86 
87  this->normals[p] = n.unit();
88  this->tangents[p][0] = this->dxyzdxi_map[p].unit();
89 #if LIBMESH_DIM == 3 // Only good in 3D space
90  this->tangents[p][1] = this->dxyzdxi_map[p].cross(n).unit();
91 
92 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
93  // The curvature is computed via the familiar Frenet formula:
94  // curvature = [d^2(x) / d (xi)^2] dot [normal]
95  // For a reference, see:
96  // F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill, p. 310
97  //
98  // Note: The sign convention here is different from the
99  // 3D case. Concave-upward curves (smiles) have a positive
100  // curvature. Concave-downward curves (frowns) have a
101  // negative curvature. Be sure to take that into account!
102  const Real numerator = this->d2xyzdxi2_map[p] * this->normals[p];
103  const Real denominator = this->dxyzdxi_map[p].norm_sq();
104  libmesh_assert_not_equal_to (denominator, 0);
105  this->curvatures[p] = numerator / denominator;
106 #endif
107 #endif
108  }
109 
110  // compute the jacobian at the quadrature points
111  for (unsigned int p=0; p<n_qp; p++)
112  {
113  const Real the_jac = std::sqrt(this->dxdxi_map(p)*this->dxdxi_map(p) +
114  this->dydxi_map(p)*this->dydxi_map(p));
115 
116  libmesh_assert_greater (the_jac, 0.);
117 
118  this->JxW[p] = the_jac*qw[p];
119  }
120 
121  break;
122  }
123 
124  case 3:
125  {
126  // Resize the vectors to hold data at the quadrature points
127  {
128  this->xyz.resize(n_qp);
129  this->dxyzdxi_map.resize(n_qp);
130  this->dxyzdeta_map.resize(n_qp);
131 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
132  this->d2xyzdxi2_map.resize(n_qp);
133  this->d2xyzdxideta_map.resize(n_qp);
134  this->d2xyzdeta2_map.resize(n_qp);
135  this->curvatures.resize(n_qp);
136 #endif
137  this->tangents.resize(n_qp);
138  this->normals.resize(n_qp);
139 
140  this->JxW.resize(n_qp);
141  }
142 
143  // Clear the entities that will be summed
144  for (unsigned int p=0; p<n_qp; p++)
145  {
146  this->tangents[p].resize(LIBMESH_DIM-1); // 1 Tangent in 2D, 2 in 3D
147  this->xyz[p].zero();
148  this->dxyzdxi_map[p].zero();
149  this->dxyzdeta_map[p].zero();
150 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
151  this->d2xyzdxi2_map[p].zero();
152  this->d2xyzdxideta_map[p].zero();
153  this->d2xyzdeta2_map[p].zero();
154 #endif
155  }
156 
157  // compute x, dxdxi at the quadrature points
158  for (unsigned int i=0,
159  n_psi_map = cast_int<unsigned int>(this->psi_map.size());
160  i != n_psi_map; i++) // sum over the nodes
161  {
162  const Point & side_point = side->point(i);
163 
164  for (unsigned int p=0; p<n_qp; p++) // for each quadrature point...
165  {
166  this->xyz[p].add_scaled (side_point, this->psi_map[i][p]);
167  this->dxyzdxi_map[p].add_scaled (side_point, this->dpsidxi_map[i][p]);
168  this->dxyzdeta_map[p].add_scaled (side_point, this->dpsideta_map[i][p]);
169 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
170  this->d2xyzdxi2_map[p].add_scaled (side_point, this->d2psidxi2_map[i][p]);
171  this->d2xyzdxideta_map[p].add_scaled(side_point, this->d2psidxideta_map[i][p]);
172  this->d2xyzdeta2_map[p].add_scaled (side_point, this->d2psideta2_map[i][p]);
173 #endif
174  }
175  }
176 
177  // Compute the tangents, normal, and curvature at the quadrature point
178  for (unsigned int p=0; p<n_qp; p++)
179  {
180  const Point n = this->dxyzdxi_map[p].cross(this->dxyzdeta_map[p]);
181  this->normals[p] = n.unit();
182  this->tangents[p][0] = this->dxyzdxi_map[p].unit();
183  this->tangents[p][1] = n.cross(this->dxyzdxi_map[p]).unit();
184 
185 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
186  // Compute curvature using the typical nomenclature
187  // of the first and second fundamental forms.
188  // For reference, see:
189  // 1) http://mathworld.wolfram.com/MeanCurvature.html
190  // (note -- they are using inward normal)
191  // 2) F.S. Merritt, Mathematics Manual, 1962, McGraw-Hill
192  const Real L = -this->d2xyzdxi2_map[p] * this->normals[p];
193  const Real M = -this->d2xyzdxideta_map[p] * this->normals[p];
194  const Real N = -this->d2xyzdeta2_map[p] * this->normals[p];
195  const Real E = this->dxyzdxi_map[p].norm_sq();
196  const Real F = this->dxyzdxi_map[p] * this->dxyzdeta_map[p];
197  const Real G = this->dxyzdeta_map[p].norm_sq();
198 
199  const Real numerator = E*N -2.*F*M + G*L;
200  const Real denominator = E*G - F*F;
201  libmesh_assert_not_equal_to (denominator, 0.);
202  this->curvatures[p] = 0.5*numerator/denominator;
203 #endif
204  }
205 
206  // compute the jacobian at the quadrature points, see
207  // http://sp81.msi.umn.edu:999/fluent/fidap/help/theory/thtoc.htm
208  for (unsigned int p=0; p<n_qp; p++)
209  {
210  const Real g11 = (this->dxdxi_map(p)*this->dxdxi_map(p) +
211  this->dydxi_map(p)*this->dydxi_map(p) +
212  this->dzdxi_map(p)*this->dzdxi_map(p));
213 
214  const Real g12 = (this->dxdxi_map(p)*this->dxdeta_map(p) +
215  this->dydxi_map(p)*this->dydeta_map(p) +
216  this->dzdxi_map(p)*this->dzdeta_map(p));
217 
218  const Real g21 = g12;
219 
220  const Real g22 = (this->dxdeta_map(p)*this->dxdeta_map(p) +
221  this->dydeta_map(p)*this->dydeta_map(p) +
222  this->dzdeta_map(p)*this->dzdeta_map(p));
223 
224 
225  const Real the_jac = std::sqrt(g11*g22 - g12*g21);
226 
227  libmesh_assert_greater (the_jac, 0.);
228 
229  this->JxW[p] = the_jac*qw[p];
230  }
231 
232  break;
233  }
234  default:
235  libmesh_error_msg("Invalid dim = " << dim);
236 
237  } // switch(dim)
238 }
239 
240 } // namespace libMesh
libMesh::FEMap::d2psidxi2_map
std::vector< std::vector< Real > > d2psidxi2_map
Map for the second derivatives (in xi) of the side shape functions.
Definition: fe_map.h:926
libMesh::FEMap::dydeta_map
Real dydeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
Definition: fe_map.h:675
libMesh::FEMap::dpsidxi_map
std::vector< std::vector< Real > > dpsidxi_map
Map for the derivative of the side functions, d(psi)/d(xi).
Definition: fe_map.h:911
libMesh::FEMap::d2xyzdxi2_map
std::vector< RealGradient > d2xyzdxi2_map
Vector of second partial derivatives in xi: d^2(x)/d(xi)^2, d^2(y)/d(xi)^2, d^2(z)/d(xi)^2.
Definition: fe_map.h:738
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::FEXYZMap::compute_face_map
virtual void compute_face_map(int dim, const std::vector< Real > &qw, const Elem *side) override
Special implementation for XYZ finite elements.
Definition: fe_xyz_map.C:25
libMesh::FEMap::JxW
std::vector< Real > JxW
Jacobian*Weight values at quadrature points.
Definition: fe_map.h:972
libMesh::FEMap::dxdxi_map
Real dxdxi_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
Definition: fe_map.h:643
std::sqrt
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::FEMap::dxyzdeta_map
std::vector< RealGradient > dxyzdeta_map
Vector of partial derivatives: d(x)/d(eta), d(y)/d(eta), d(z)/d(eta)
Definition: fe_map.h:724
libMesh::FEMap::dydxi_map
Real dydxi_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
Definition: fe_map.h:651
libMesh::FEMap::normals
std::vector< Point > normals
Normal vectors on boundary at quadrature points.
Definition: fe_map.h:952
libMesh::FEMap::xyz
std::vector< Point > xyz
The spatial locations of the quadrature points.
Definition: fe_map.h:712
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::Elem::point
const Point & point(const unsigned int i) const
Definition: elem.h:1955
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::FEMap::dpsideta_map
std::vector< std::vector< Real > > dpsideta_map
Map for the derivative of the side function, d(psi)/d(eta).
Definition: fe_map.h:917
libMesh::FEMap::tangents
std::vector< std::vector< Point > > tangents
Tangent vectors on boundary at quadrature points.
Definition: fe_map.h:947
libMesh::FEMap::psi_map
std::vector< std::vector< Real > > psi_map
Map for the side shape functions, psi.
Definition: fe_map.h:905
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::FEMap::d2psideta2_map
std::vector< std::vector< Real > > d2psideta2_map
Map for the second derivatives (in eta) of the side shape functions.
Definition: fe_map.h:940
libMesh::FEMap::dzdeta_map
Real dzdeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
Definition: fe_map.h:683
libMesh::TypeVector::unit
TypeVector< T > unit() const
Definition: type_vector.h:1158
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::FEMap::dzdxi_map
Real dzdxi_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
Definition: fe_map.h:659
libMesh::TypeVector::cross
TypeVector< typename CompareTypes< T, T2 >::supertype > cross(const TypeVector< T2 > &v) const
Definition: type_vector.h:920
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FEMap::dxdeta_map
Real dxdeta_map(const unsigned int p) const
Used in FEMap::compute_map(), which should be be usable in derived classes, and therefore protected.
Definition: fe_map.h:667
libMesh::FEMap::dxyzdxi_map
std::vector< RealGradient > dxyzdxi_map
Vector of partial derivatives: d(x)/d(xi), d(y)/d(xi), d(z)/d(xi)
Definition: fe_map.h:718
libMesh::FEMap::d2xyzdxideta_map
std::vector< RealGradient > d2xyzdxideta_map
Vector of mixed second partial derivatives in xi-eta: d^2(x)/d(xi)d(eta) d^2(y)/d(xi)d(eta) d^2(z)/d(...
Definition: fe_map.h:744
libMesh::FEMap::d2xyzdeta2_map
std::vector< RealGradient > d2xyzdeta2_map
Vector of second partial derivatives in eta: d^2(x)/d(eta)^2.
Definition: fe_map.h:750
libMesh::FEMap::d2psidxideta_map
std::vector< std::vector< Real > > d2psidxideta_map
Map for the second (cross) derivatives in xi, eta of the side shape functions.
Definition: fe_map.h:933
libMesh::FEMap::curvatures
std::vector< Real > curvatures
The mean curvature (= one half the sum of the principal curvatures) on the boundary at the quadrature...
Definition: fe_map.h:960