libMesh
fe_map.h
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 
19 
20 #ifndef LIBMESH_FE_MAP_H
21 #define LIBMESH_FE_MAP_H
22 
23 // libMesh includes
24 #include "libmesh/reference_counted_object.h"
25 #include "libmesh/point.h"
26 #include "libmesh/vector_value.h"
27 #include "libmesh/fe_type.h"
28 
29 // C++ includes
30 #include <memory>
31 
32 namespace libMesh
33 {
34 
35 // forward declarations
36 class Elem;
37 class Node;
38 
47 class FEMap
48 {
49 public:
50 
51  FEMap(Real jtol = 0);
52  virtual ~FEMap(){}
53 
54  static std::unique_ptr<FEMap> build(FEType fe_type);
55 
56  static FEFamily map_fe_type(const Elem & elem);
57 
58  template<unsigned int Dim>
59  void init_reference_to_physical_map(const std::vector<Point> & qp,
60  const Elem * elem);
61 
69  void compute_single_point_map(const unsigned int dim,
70  const std::vector<Real> & qw,
71  const Elem * elem,
72  unsigned int p,
73  const std::vector<const Node *> & elem_nodes,
74  bool compute_second_derivatives);
75 
82  virtual void compute_affine_map(const unsigned int dim,
83  const std::vector<Real> & qw,
84  const Elem * elem);
85 
91  virtual void compute_null_map(const unsigned int dim,
92  const std::vector<Real> & qw);
93 
94 
103  virtual void compute_map(const unsigned int dim,
104  const std::vector<Real> & qw,
105  const Elem * elem,
106  bool calculate_d2phi);
107 
111  virtual void compute_face_map(int dim,
112  const std::vector<Real> & qw,
113  const Elem * side);
114 
118  void compute_edge_map(int dim,
119  const std::vector<Real> & qw,
120  const Elem * side);
121 
126  template<unsigned int Dim>
127  void init_face_shape_functions(const std::vector<Point> & qp,
128  const Elem * side);
129 
134  template<unsigned int Dim>
135  void init_edge_shape_functions(const std::vector<Point> & qp,
136  const Elem * edge);
137 
142  static Point map (const unsigned int dim,
143  const Elem * elem,
144  const Point & reference_point);
145 
150  static Point map_deriv (const unsigned int dim,
151  const Elem * elem,
152  const unsigned int j,
153  const Point & reference_point);
154 
166  static Point inverse_map (const unsigned int dim,
167  const Elem * elem,
168  const Point & p,
169  const Real tolerance = TOLERANCE,
170  const bool secure = true);
171 
184  static void inverse_map (unsigned int dim,
185  const Elem * elem,
186  const std::vector<Point> & physical_points,
187  std::vector<Point> & reference_points,
188  const Real tolerance = TOLERANCE,
189  const bool secure = true);
190 
195  const std::vector<Point> & get_xyz() const
197  calculate_xyz = true; return xyz; }
198 
202  const std::vector<Real> & get_jacobian() const
204  calculate_dxyz = true; return jac; }
205 
210  const std::vector<Real> & get_JxW() const
212  calculate_dxyz = true; return JxW; }
213 
218  const std::vector<RealGradient> & get_dxyzdxi() const
220  calculate_dxyz = true; return dxyzdxi_map; }
221 
226  const std::vector<RealGradient> & get_dxyzdeta() const
228  calculate_dxyz = true; return dxyzdeta_map; }
229 
234  const std::vector<RealGradient> & get_dxyzdzeta() const
236  calculate_dxyz = true; return dxyzdzeta_map; }
237 
238 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
239 
243  const std::vector<RealGradient> & get_d2xyzdxi2() const
245  calculate_d2xyz = true; return d2xyzdxi2_map; }
246 
250  const std::vector<RealGradient> & get_d2xyzdeta2() const
252  calculate_d2xyz = true; return d2xyzdeta2_map; }
253 
257  const std::vector<RealGradient> & get_d2xyzdzeta2() const
259  calculate_d2xyz = true; return d2xyzdzeta2_map; }
260 
264  const std::vector<RealGradient> & get_d2xyzdxideta() const
266  calculate_d2xyz = true; return d2xyzdxideta_map; }
267 
271  const std::vector<RealGradient> & get_d2xyzdxidzeta() const
273  calculate_d2xyz = true; return d2xyzdxidzeta_map; }
274 
278  const std::vector<RealGradient> & get_d2xyzdetadzeta() const
280  calculate_d2xyz = true; return d2xyzdetadzeta_map; }
281 
282 #endif
283 
288  const std::vector<Real> & get_dxidx() const
290  calculate_dxyz = true; return dxidx_map; }
291 
296  const std::vector<Real> & get_dxidy() const
298  calculate_dxyz = true; return dxidy_map; }
299 
304  const std::vector<Real> & get_dxidz() const
306  calculate_dxyz = true; return dxidz_map; }
307 
312  const std::vector<Real> & get_detadx() const
314  calculate_dxyz = true; return detadx_map; }
315 
320  const std::vector<Real> & get_detady() const
322  calculate_dxyz = true; return detady_map; }
323 
328  const std::vector<Real> & get_detadz() const
330  calculate_dxyz = true; return detadz_map; }
331 
336  const std::vector<Real> & get_dzetadx() const
338  calculate_dxyz = true; return dzetadx_map; }
339 
344  const std::vector<Real> & get_dzetady() const
346  calculate_dxyz = true; return dzetady_map; }
347 
352  const std::vector<Real> & get_dzetadz() const
354  calculate_dxyz = true; return dzetadz_map; }
355 
356 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
357 
360  const std::vector<std::vector<Real>> & get_d2xidxyz2() const
362  calculate_d2xyz = true; return d2xidxyz2_map; }
363 
367  const std::vector<std::vector<Real>> & get_d2etadxyz2() const
369  calculate_d2xyz = true; return d2etadxyz2_map; }
370 
374  const std::vector<std::vector<Real>> & get_d2zetadxyz2() const
376  calculate_d2xyz = true; return d2zetadxyz2_map; }
377 #endif
378 
382  const std::vector<std::vector<Real>> & get_psi() const
383  { return psi_map; }
384 
388  const std::vector<std::vector<Real>> & get_phi_map() const
390  calculate_xyz = true; return phi_map; }
391 
395  const std::vector<std::vector<Real>> & get_dphidxi_map() const
397  calculate_dxyz = true; return dphidxi_map; }
398 
402  const std::vector<std::vector<Real>> & get_dphideta_map() const
404  calculate_dxyz = true; return dphideta_map; }
405 
409  const std::vector<std::vector<Real>> & get_dphidzeta_map() const
411  calculate_dxyz = true; return dphidzeta_map; }
412 
416  const std::vector<std::vector<Point>> & get_tangents() const
418  calculate_dxyz = true; return tangents; }
419 
423  const std::vector<Point> & get_normals() const
425  calculate_dxyz = true; return normals; }
426 
427 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
428 
432  const std::vector<Real> & get_curvatures() const
434  calculate_d2xyz = true; return curvatures;}
435 
436 #endif
437 
441  void print_JxW(std::ostream & os) const;
442 
447  void print_xyz(std::ostream & os) const;
448 
449  /* FIXME: PB: The public functions below break encapsulation! I'm not happy
450  about it, but I don't have time to redo the infinite element routines.
451  These are needed in inf_fe_boundary.C and inf_fe.C. A proper implementation
452  would subclass this class and hide all the radial function business. */
456  std::vector<std::vector<Real>> & get_psi()
457  { return psi_map; }
458 
462  std::vector<std::vector<Real>> & get_dpsidxi()
464  calculate_dxyz = true; return dpsidxi_map; }
465 
466  const std::vector<std::vector<Real>> & get_dpsidxi() const
468  calculate_dxyz = true; return dpsidxi_map; }
469 
473  std::vector<std::vector<Real>> & get_dpsideta()
475  calculate_dxyz = true; return dpsideta_map; }
476 
477  const std::vector<std::vector<Real>> & get_dpsideta() const
479  calculate_dxyz = true; return dpsideta_map; }
480 
481 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
482 
486  std::vector<std::vector<Real>> & get_d2psidxi2()
488  calculate_d2xyz = true; return d2psidxi2_map; }
489 
493  const std::vector<std::vector<Real>> & get_d2psidxi2() const
495  calculate_d2xyz = true; return d2psidxi2_map; }
496 
500  std::vector<std::vector<Real>> & get_d2psidxideta()
502  calculate_d2xyz = true; return d2psidxideta_map; }
503 
507  const std::vector<std::vector<Real>> & get_d2psidxideta() const
509  calculate_d2xyz = true; return d2psidxideta_map; }
510 
514  std::vector<std::vector<Real>> & get_d2psideta2()
516  calculate_d2xyz = true; return d2psideta2_map; }
517 
518 
522  const std::vector<std::vector<Real>> & get_d2psideta2() const
524  calculate_d2xyz = true; return d2psideta2_map; }
525 
526 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
527 
531  std::vector<std::vector<Real>> & get_phi_map()
533  calculate_xyz = true; return phi_map; }
534 
538  std::vector<std::vector<Real>> & get_dphidxi_map()
540  calculate_dxyz = true; return dphidxi_map; }
541 
545  std::vector<std::vector<Real>> & get_dphideta_map()
547  calculate_dxyz = true; return dphideta_map; }
548 
552  std::vector<std::vector<Real>> & get_dphidzeta_map()
554  calculate_dxyz = true; return dphidzeta_map; }
555 
556 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
557 
560  std::vector<std::vector<Real>> & get_d2phidxi2_map()
562  calculate_d2xyz = true; return d2phidxi2_map; }
563 
567  std::vector<std::vector<Real>> & get_d2phidxideta_map()
569  calculate_d2xyz = true; return d2phidxideta_map; }
570 
574  std::vector<std::vector<Real>> & get_d2phidxidzeta_map()
576  calculate_d2xyz = true; return d2phidxidzeta_map; }
577 
581  std::vector<std::vector<Real>> & get_d2phideta2_map()
583  calculate_d2xyz = true; return d2phideta2_map; }
584 
588  std::vector<std::vector<Real>> & get_d2phidetadzeta_map()
590  calculate_d2xyz = true; return d2phidetadzeta_map; }
591 
595  std::vector<std::vector<Real>> & get_d2phidzeta2_map()
597  calculate_d2xyz = true; return d2phidzeta2_map; }
598 #endif
599 
600  /* FIXME: PB: This function breaks encapsulation! Needed in FE<>::reinit and
601  InfFE<>::reinit. Not sure yet if the algorithm can be redone to avoid this. */
606  std::vector<Real> & get_JxW()
608  calculate_dxyz = true; return JxW; }
609 
615 
616 protected:
617 
622  calculations_started = true;
623 
624 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
625  // Second derivative calculations currently have first derivative
626  // calculations as a prerequisite
627  if (calculate_d2xyz)
628  calculate_dxyz = true;
629 #endif
630  }
631 
635  void resize_quadrature_map_vectors(const unsigned int dim, unsigned int n_qp);
636 
643  Real dxdxi_map(const unsigned int p) const { return dxyzdxi_map[p](0); }
644 
651  Real dydxi_map(const unsigned int p) const { return dxyzdxi_map[p](1); }
652 
659  Real dzdxi_map(const unsigned int p) const { return dxyzdxi_map[p](2); }
660 
667  Real dxdeta_map(const unsigned int p) const { return dxyzdeta_map[p](0); }
668 
675  Real dydeta_map(const unsigned int p) const { return dxyzdeta_map[p](1); }
676 
683  Real dzdeta_map(const unsigned int p) const { return dxyzdeta_map[p](2); }
684 
691  Real dxdzeta_map(const unsigned int p) const { return dxyzdzeta_map[p](0); }
692 
699  Real dydzeta_map(const unsigned int p) const { return dxyzdzeta_map[p](1); }
700 
707  Real dzdzeta_map(const unsigned int p) const { return dxyzdzeta_map[p](2); }
708 
712  std::vector<Point> xyz;
713 
718  std::vector<RealGradient> dxyzdxi_map;
719 
724  std::vector<RealGradient> dxyzdeta_map;
725 
730  std::vector<RealGradient> dxyzdzeta_map;
731 
732 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
733 
738  std::vector<RealGradient> d2xyzdxi2_map;
739 
744  std::vector<RealGradient> d2xyzdxideta_map;
745 
750  std::vector<RealGradient> d2xyzdeta2_map;
751 
756  std::vector<RealGradient> d2xyzdxidzeta_map;
757 
762  std::vector<RealGradient> d2xyzdetadzeta_map;
763 
768  std::vector<RealGradient> d2xyzdzeta2_map;
769 
770 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
771 
776  std::vector<Real> dxidx_map;
777 
782  std::vector<Real> dxidy_map;
783 
788  std::vector<Real> dxidz_map;
789 
790 
795  std::vector<Real> detadx_map;
796 
801  std::vector<Real> detady_map;
802 
807  std::vector<Real> detadz_map;
808 
809 
814  std::vector<Real> dzetadx_map;
815 
820  std::vector<Real> dzetady_map;
821 
826  std::vector<Real> dzetadz_map;
827 
828 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
829 
833  std::vector<std::vector<Real>> d2xidxyz2_map;
834 
839  std::vector<std::vector<Real>> d2etadxyz2_map;
840 
845  std::vector<std::vector<Real>> d2zetadxyz2_map;
846 #endif
847 
851  std::vector<std::vector<Real>> phi_map;
852 
856  std::vector<std::vector<Real>> dphidxi_map;
857 
861  std::vector<std::vector<Real>> dphideta_map;
862 
866  std::vector<std::vector<Real>> dphidzeta_map;
867 
868 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
869 
873  std::vector<std::vector<Real>> d2phidxi2_map;
874 
878  std::vector<std::vector<Real>> d2phidxideta_map;
879 
883  std::vector<std::vector<Real>> d2phidxidzeta_map;
884 
888  std::vector<std::vector<Real>> d2phideta2_map;
889 
893  std::vector<std::vector<Real>> d2phidetadzeta_map;
894 
898  std::vector<std::vector<Real>> d2phidzeta2_map;
899 
900 #endif
901 
905  std::vector<std::vector<Real>> psi_map;
906 
911  std::vector<std::vector<Real>> dpsidxi_map;
912 
917  std::vector<std::vector<Real>> dpsideta_map;
918 
919 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
920 
926  std::vector<std::vector<Real>> d2psidxi2_map;
927 
933  std::vector<std::vector<Real>> d2psidxideta_map;
934 
940  std::vector<std::vector<Real>> d2psideta2_map;
941 
942 #endif
943 
947  std::vector<std::vector<Point>> tangents;
948 
952  std::vector<Point> normals;
953 
954 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
955 
960  std::vector<Real> curvatures;
961 
962 #endif
963 
967  std::vector<Real> jac;
968 
972  std::vector<Real> JxW;
973 
978  mutable bool calculations_started;
979 
983  mutable bool calculate_xyz;
984 
988  mutable bool calculate_dxyz;
989 
990 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
991 
995  mutable bool calculate_d2xyz;
996 
997 #endif
998 
1003  template <unsigned int Dim, FEFamily T>
1004  friend class FE;
1005 
1012 
1013 private:
1018  void compute_inverse_map_second_derivs(unsigned p);
1019 
1023  std::vector<const Node *> _elem_nodes;
1024 };
1025 
1026 }
1027 
1028 #endif // LIBMESH_FE_MAP_H
libMesh::FEMap::get_jacobian
const std::vector< Real > & get_jacobian() const
Definition: fe_map.h:202
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::get_psi
const std::vector< std::vector< Real > > & get_psi() const
Definition: fe_map.h:382
libMesh::FEMap::dzetadx_map
std::vector< Real > dzetadx_map
Map for partial derivatives: d(zeta)/d(x).
Definition: fe_map.h:814
libMesh::FEMap::calculate_d2xyz
bool calculate_d2xyz
Should we calculate mapping hessians?
Definition: fe_map.h:995
libMesh::FEMap::d2xyzdetadzeta_map
std::vector< RealGradient > d2xyzdetadzeta_map
Vector of mixed second partial derivatives in eta-zeta: d^2(x)/d(eta)d(zeta) d^2(y)/d(eta)d(zeta) d^2...
Definition: fe_map.h:762
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::get_d2psidxideta
std::vector< std::vector< Real > > & get_d2psidxideta()
Definition: fe_map.h:500
libMesh::FEMap::get_d2zetadxyz2
const std::vector< std::vector< Real > > & get_d2zetadxyz2() const
Second derivatives of "zeta" reference coordinate wrt physical coordinates.
Definition: fe_map.h:374
libMesh::FEMap::get_d2xyzdzeta2
const std::vector< RealGradient > & get_d2xyzdzeta2() const
Definition: fe_map.h:257
libMesh::FEMap::get_d2phidxidzeta_map
std::vector< std::vector< Real > > & get_d2phidxidzeta_map()
Definition: fe_map.h:574
libMesh::FEMap::get_d2xyzdxi2
const std::vector< RealGradient > & get_d2xyzdxi2() const
Definition: fe_map.h:243
libMesh::FEMap::set_jacobian_tolerance
void set_jacobian_tolerance(Real tol)
Set the Jacobian tolerance used for determining when the mapping fails.
Definition: fe_map.h:614
libMesh::FEMap::get_d2xyzdeta2
const std::vector< RealGradient > & get_d2xyzdeta2() const
Definition: fe_map.h:250
libMesh::FEMap::get_dzetadx
const std::vector< Real > & get_dzetadx() const
Definition: fe_map.h:336
libMesh::FEMap::get_dxidx
const std::vector< Real > & get_dxidx() const
Definition: fe_map.h:288
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::get_d2phidzeta2_map
std::vector< std::vector< Real > > & get_d2phidzeta2_map()
Definition: fe_map.h:595
libMesh::FEMap::map_deriv
static Point map_deriv(const unsigned int dim, const Elem *elem, const unsigned int j, const Point &reference_point)
Definition: fe_map.C:2076
libMesh::FEMap::dzetadz_map
std::vector< Real > dzetadz_map
Map for partial derivatives: d(zeta)/d(z).
Definition: fe_map.h:826
libMesh::FEMap::get_dxyzdxi
const std::vector< RealGradient > & get_dxyzdxi() const
Definition: fe_map.h:218
libMesh::FEMap::dxidy_map
std::vector< Real > dxidy_map
Map for partial derivatives: d(xi)/d(y).
Definition: fe_map.h:782
libMesh::FEMap::get_d2psidxideta
const std::vector< std::vector< Real > > & get_d2psidxideta() const
Definition: fe_map.h:507
libMesh::FEMap::get_dphidxi_map
const std::vector< std::vector< Real > > & get_dphidxi_map() const
Definition: fe_map.h:395
libMesh::FEMap::FEMap
FEMap(Real jtol=0)
Definition: fe_map.C:64
libMesh::FEMap::calculate_xyz
bool calculate_xyz
Should we calculate physical point locations?
Definition: fe_map.h:983
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::FEMap::get_dzetadz
const std::vector< Real > & get_dzetadz() const
Definition: fe_map.h:352
libMesh::FEMap::jac
std::vector< Real > jac
Jacobian values at quadrature points.
Definition: fe_map.h:967
libMesh::FEMap::get_dzetady
const std::vector< Real > & get_dzetady() const
Definition: fe_map.h:344
libMesh::FEMap::dxdzeta_map
Real dxdzeta_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:691
libMesh::FEMap::d2phidxideta_map
std::vector< std::vector< Real > > d2phidxideta_map
Map for the second derivative, d^2(phi)/d(xi)d(eta).
Definition: fe_map.h:878
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
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::FEMap::get_psi
std::vector< std::vector< Real > > & get_psi()
Definition: fe_map.h:456
libMesh::FEMap::print_JxW
void print_JxW(std::ostream &os) const
Prints the Jacobian times the weight for each quadrature point.
Definition: fe_map.C:1486
libMesh::FEMap::map
static Point map(const unsigned int dim, const Elem *elem, const Point &reference_point)
Definition: fe_map.C:2043
libMesh::FEMap::get_JxW
const std::vector< Real > & get_JxW() const
Definition: fe_map.h:210
libMesh::FEMap::_elem_nodes
std::vector< const Node * > _elem_nodes
Work vector for compute_affine_map()
Definition: fe_map.h:1023
libMesh::FEMap::dydzeta_map
Real dydzeta_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:699
libMesh::FEMap::get_xyz
const std::vector< Point > & get_xyz() const
Definition: fe_map.h:195
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::get_d2psidxi2
std::vector< std::vector< Real > > & get_d2psidxi2()
Definition: fe_map.h:486
libMesh::FEMap::get_d2xidxyz2
const std::vector< std::vector< Real > > & get_d2xidxyz2() const
Second derivatives of "xi" reference coordinate wrt physical coordinates.
Definition: fe_map.h:360
libMesh::FEMap::normals
std::vector< Point > normals
Normal vectors on boundary at quadrature points.
Definition: fe_map.h:952
libMesh::FEMap::get_dxidy
const std::vector< Real > & get_dxidy() const
Definition: fe_map.h:296
libMesh::FEMap::xyz
std::vector< Point > xyz
The spatial locations of the quadrature points.
Definition: fe_map.h:712
libMesh::FEMap::init_face_shape_functions
void init_face_shape_functions(const std::vector< Point > &qp, const Elem *side)
Initializes the reference to physical element map for a side.
Definition: fe_boundary.C:381
libMesh::FEMap::~FEMap
virtual ~FEMap()
Definition: fe_map.h:52
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::FEMap::compute_affine_map
virtual void compute_affine_map(const unsigned int dim, const std::vector< Real > &qw, const Elem *elem)
Compute the jacobian and some other additional data fields.
Definition: fe_map.C:1277
libMesh::FEMap::get_dphideta_map
std::vector< std::vector< Real > > & get_dphideta_map()
Definition: fe_map.h:545
libMesh::FEMap::d2phidxi2_map
std::vector< std::vector< Real > > d2phidxi2_map
Map for the second derivative, d^2(phi)/d(xi)^2.
Definition: fe_map.h:873
libMesh::FEMap::compute_null_map
virtual void compute_null_map(const unsigned int dim, const std::vector< Real > &qw)
Assign a fake jacobian and some other additional data fields.
Definition: fe_map.C:1358
libMesh::FEMap::detady_map
std::vector< Real > detady_map
Map for partial derivatives: d(eta)/d(y).
Definition: fe_map.h:801
libMesh::FEMap::compute_single_point_map
void compute_single_point_map(const unsigned int dim, const std::vector< Real > &qw, const Elem *elem, unsigned int p, const std::vector< const Node * > &elem_nodes, bool compute_second_derivatives)
Compute the jacobian and some other additional data fields at the single point with index p.
Definition: fe_map.C:446
libMesh::FEMap::print_xyz
void print_xyz(std::ostream &os) const
Prints the spatial location of each quadrature point (on the physical element).
Definition: fe_map.C:1494
libMesh::FEMap::get_dxyzdeta
const std::vector< RealGradient > & get_dxyzdeta() const
Definition: fe_map.h:226
libMesh::FEMap::compute_face_map
virtual void compute_face_map(int dim, const std::vector< Real > &qw, const Elem *side)
Same as compute_map, but for a side.
Definition: fe_boundary.C:587
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::FEMap::get_d2psidxi2
const std::vector< std::vector< Real > > & get_d2psidxi2() const
Definition: fe_map.h:493
libMesh::FEMap::get_detady
const std::vector< Real > & get_detady() const
Definition: fe_map.h:320
libMesh::FEMap::get_d2xyzdxidzeta
const std::vector< RealGradient > & get_d2xyzdxidzeta() const
Definition: fe_map.h:271
libMesh::FEMap::determine_calculations
void determine_calculations()
Determine which values are to be calculated.
Definition: fe_map.h:621
libMesh::FEMap::get_phi_map
std::vector< std::vector< Real > > & get_phi_map()
Definition: fe_map.h:531
libMesh::FEMap::d2xyzdxidzeta_map
std::vector< RealGradient > d2xyzdxidzeta_map
Vector of second partial derivatives in xi-zeta: d^2(x)/d(xi)d(zeta), d^2(y)/d(xi)d(zeta),...
Definition: fe_map.h:756
libMesh::FEMap::inverse_map
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: fe_map.C:1622
libMesh::FEMap::resize_quadrature_map_vectors
void resize_quadrature_map_vectors(const unsigned int dim, unsigned int n_qp)
A utility function for use by compute_*_map.
Definition: fe_map.C:1196
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::get_dpsideta
std::vector< std::vector< Real > > & get_dpsideta()
Definition: fe_map.h:473
libMesh::FEMap::detadz_map
std::vector< Real > detadz_map
Map for partial derivatives: d(eta)/d(z).
Definition: fe_map.h:807
libMesh::FEMap::tangents
std::vector< std::vector< Point > > tangents
Tangent vectors on boundary at quadrature points.
Definition: fe_map.h:947
libMesh::FEMap::get_d2phidxideta_map
std::vector< std::vector< Real > > & get_d2phidxideta_map()
Definition: fe_map.h:567
libMesh::FE
A specific instantiation of the FEBase class.
Definition: fe.h:95
libMesh::FEMap::psi_map
std::vector< std::vector< Real > > psi_map
Map for the side shape functions, psi.
Definition: fe_map.h:905
libMesh::FEMap::d2phidzeta2_map
std::vector< std::vector< Real > > d2phidzeta2_map
Map for the second derivative, d^2(phi)/d(zeta)^2.
Definition: fe_map.h:898
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::FEMap::calculate_dxyz
bool calculate_dxyz
Should we calculate mapping gradients?
Definition: fe_map.h:988
libMesh::FEMap::get_d2xyzdxideta
const std::vector< RealGradient > & get_d2xyzdxideta() const
Definition: fe_map.h:264
libMesh::FEMap::d2xyzdzeta2_map
std::vector< RealGradient > d2xyzdzeta2_map
Vector of second partial derivatives in zeta: d^2(x)/d(zeta)^2.
Definition: fe_map.h:768
libMesh::FEMap::dphideta_map
std::vector< std::vector< Real > > dphideta_map
Map for the derivative, d(phi)/d(eta).
Definition: fe_map.h:861
libMesh::FEMap::build
static std::unique_ptr< FEMap > build(FEType fe_type)
Definition: fe_map.C:76
libMesh::FEMap::get_dphidzeta_map
const std::vector< std::vector< Real > > & get_dphidzeta_map() const
Definition: fe_map.h:409
libMesh::FEMap::get_detadx
const std::vector< Real > & get_detadx() const
Definition: fe_map.h:312
libMesh::FEMap::d2xidxyz2_map
std::vector< std::vector< Real > > d2xidxyz2_map
Second derivatives of "xi" reference coordinate wrt physical coordinates.
Definition: fe_map.h:833
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::compute_edge_map
void compute_edge_map(int dim, const std::vector< Real > &qw, const Elem *side)
Same as before, but for an edge.
Definition: fe_boundary.C:950
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::FEMap::map_fe_type
static FEFamily map_fe_type(const Elem &elem)
Definition: fe_map.C:47
libMesh::FEMap::d2etadxyz2_map
std::vector< std::vector< Real > > d2etadxyz2_map
Second derivatives of "eta" reference coordinate wrt physical coordinates.
Definition: fe_map.h:839
libMesh::FEMap::detadx_map
std::vector< Real > detadx_map
Map for partial derivatives: d(eta)/d(x).
Definition: fe_map.h:795
libMesh::FEMap::get_d2psideta2
const std::vector< std::vector< Real > > & get_d2psideta2() const
Definition: fe_map.h:522
libMesh::FEMap::init_reference_to_physical_map
void init_reference_to_physical_map(const std::vector< Point > &qp, const Elem *elem)
Definition: fe_map.C:91
libMesh::FEMap::jacobian_tolerance
Real jacobian_tolerance
The Jacobian tolerance used for determining when the mapping fails.
Definition: fe_map.h:1011
libMesh::FEMap::init_edge_shape_functions
void init_edge_shape_functions(const std::vector< Point > &qp, const Elem *edge)
Same as before, but for an edge.
Definition: fe_boundary.C:510
libMesh::FEMap::dphidxi_map
std::vector< std::vector< Real > > dphidxi_map
Map for the derivative, d(phi)/d(xi).
Definition: fe_map.h:856
libMesh::FEMap::get_dphidzeta_map
std::vector< std::vector< Real > > & get_dphidzeta_map()
Definition: fe_map.h:552
libMesh::FEMap::compute_inverse_map_second_derivs
void compute_inverse_map_second_derivs(unsigned p)
A helper function used by FEMap::compute_single_point_map() to compute second derivatives of the inve...
Definition: fe_map.C:1502
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::FEMap::get_JxW
std::vector< Real > & get_JxW()
Definition: fe_map.h:606
libMesh::FEMap::d2phidetadzeta_map
std::vector< std::vector< Real > > d2phidetadzeta_map
Map for the second derivative, d^2(phi)/d(eta)d(zeta).
Definition: fe_map.h:893
libMesh::FEFamily
FEFamily
Definition: enum_fe_family.h:34
libMesh::FEMap::phi_map
std::vector< std::vector< Real > > phi_map
Map for the shape function phi.
Definition: fe_map.h:851
libMesh::FEMap::compute_map
virtual void compute_map(const unsigned int dim, const std::vector< Real > &qw, const Elem *elem, bool calculate_d2phi)
Compute the jacobian and some other additional data fields.
Definition: fe_map.C:1433
libMesh::FEMap::get_dpsidxi
std::vector< std::vector< Real > > & get_dpsidxi()
Definition: fe_map.h:462
libMesh::FEMap::get_d2phideta2_map
std::vector< std::vector< Real > > & get_d2phideta2_map()
Definition: fe_map.h:581
libMesh::FEMap::dxyzdzeta_map
std::vector< RealGradient > dxyzdzeta_map
Vector of partial derivatives: d(x)/d(zeta), d(y)/d(zeta), d(z)/d(zeta)
Definition: fe_map.h:730
libMesh::FEMap::dzetady_map
std::vector< Real > dzetady_map
Map for partial derivatives: d(zeta)/d(y).
Definition: fe_map.h:820
libMesh::FEMap::dphidzeta_map
std::vector< std::vector< Real > > dphidzeta_map
Map for the derivative, d(phi)/d(zeta).
Definition: fe_map.h:866
libMesh::FEMap::get_d2phidxi2_map
std::vector< std::vector< Real > > & get_d2phidxi2_map()
Definition: fe_map.h:560
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::FEMap::get_curvatures
const std::vector< Real > & get_curvatures() const
Definition: fe_map.h:432
libMesh::FEMap::get_normals
const std::vector< Point > & get_normals() const
Definition: fe_map.h:423
libMesh::FEMap::get_detadz
const std::vector< Real > & get_detadz() const
Definition: fe_map.h:328
libMesh::FEMap::get_phi_map
const std::vector< std::vector< Real > > & get_phi_map() const
Definition: fe_map.h:388
libMesh::FEMap::get_dpsidxi
const std::vector< std::vector< Real > > & get_dpsidxi() const
Definition: fe_map.h:466
libMesh::FEMap::d2phideta2_map
std::vector< std::vector< Real > > d2phideta2_map
Map for the second derivative, d^2(phi)/d(eta)^2.
Definition: fe_map.h:888
libMesh::FEMap
Class contained in FE that encapsulates mapping (i.e.
Definition: fe_map.h:47
libMesh::FEMap::d2zetadxyz2_map
std::vector< std::vector< Real > > d2zetadxyz2_map
Second derivatives of "zeta" reference coordinate wrt physical coordinates.
Definition: fe_map.h:845
libMesh::FEMap::get_d2psideta2
std::vector< std::vector< Real > > & get_d2psideta2()
Definition: fe_map.h:514
libMesh::FEMap::d2phidxidzeta_map
std::vector< std::vector< Real > > d2phidxidzeta_map
Map for the second derivative, d^2(phi)/d(xi)d(zeta).
Definition: fe_map.h:883
libMesh::FEMap::get_dpsideta
const std::vector< std::vector< Real > > & get_dpsideta() const
Definition: fe_map.h:477
libMesh::FEMap::get_dxyzdzeta
const std::vector< RealGradient > & get_dxyzdzeta() const
Definition: fe_map.h:234
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::FEMap::dzdzeta_map
Real dzdzeta_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:707
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::get_tangents
const std::vector< std::vector< Point > > & get_tangents() const
Definition: fe_map.h:416
libMesh::FEMap::get_dphidxi_map
std::vector< std::vector< Real > > & get_dphidxi_map()
Definition: fe_map.h:538
libMesh::FEMap::dxidx_map
std::vector< Real > dxidx_map
Map for partial derivatives: d(xi)/d(x).
Definition: fe_map.h:776
libMesh::FEMap::get_dxidz
const std::vector< Real > & get_dxidz() const
Definition: fe_map.h:304
libMesh::FEMap::get_d2phidetadzeta_map
std::vector< std::vector< Real > > & get_d2phidetadzeta_map()
Definition: fe_map.h:588
libMesh::FEMap::dxidz_map
std::vector< Real > dxidz_map
Map for partial derivatives: d(xi)/d(z).
Definition: fe_map.h:788
libMesh::FEMap::calculations_started
bool calculations_started
Have calculations with this object already been started? Then all get_* functions should already have...
Definition: fe_map.h:978
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::get_dphideta_map
const std::vector< std::vector< Real > > & get_dphideta_map() const
Definition: fe_map.h:402
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::get_d2etadxyz2
const std::vector< std::vector< Real > > & get_d2etadxyz2() const
Second derivatives of "eta" reference coordinate wrt physical coordinates.
Definition: fe_map.h:367
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::get_d2xyzdetadzeta
const std::vector< RealGradient > & get_d2xyzdetadzeta() const
Definition: fe_map.h:278
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