fe_map.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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
25 #include "libmesh/point.h"
26 #include "libmesh/vector_value.h"
27 #include "libmesh/fe_type.h"
28 #include "libmesh/auto_ptr.h" // deprecated
29 
30 // C++ includes
31 #include <memory>
32 
33 namespace libMesh
34 {
35 
36 // forward declarations
37 class Elem;
38 class Node;
39 
48 class FEMap
49 {
50 public:
51 
52  FEMap();
53  virtual ~FEMap(){}
54 
55  static std::unique_ptr<FEMap> build(FEType fe_type);
56 
57  template<unsigned int Dim>
58  void init_reference_to_physical_map(const std::vector<Point> & qp,
59  const Elem * elem);
60 
68  void compute_single_point_map(const unsigned int dim,
69  const std::vector<Real> & qw,
70  const Elem * elem,
71  unsigned int p,
72  const std::vector<const Node *> & elem_nodes,
73  bool compute_second_derivatives);
74 
81  virtual void compute_affine_map(const unsigned int dim,
82  const std::vector<Real> & qw,
83  const Elem * elem);
84 
90  virtual void compute_null_map(const unsigned int dim,
91  const std::vector<Real> & qw);
92 
93 
102  virtual void compute_map(const unsigned int dim,
103  const std::vector<Real> & qw,
104  const Elem * elem,
105  bool calculate_d2phi);
106 
110  virtual void compute_face_map(int dim,
111  const std::vector<Real> & qw,
112  const Elem * side);
113 
117  void compute_edge_map(int dim,
118  const std::vector<Real> & qw,
119  const Elem * side);
120 
125  template<unsigned int Dim>
126  void init_face_shape_functions(const std::vector<Point> & qp,
127  const Elem * side);
128 
133  template<unsigned int Dim>
134  void init_edge_shape_functions(const std::vector<Point> & qp,
135  const Elem * edge);
136 
141  const std::vector<Point> & get_xyz() const
142  { libmesh_assert(!calculations_started || calculate_xyz);
143  calculate_xyz = true; return xyz; }
144 
148  const std::vector<Real> & get_jacobian() const
149  { libmesh_assert(!calculations_started || calculate_dxyz);
150  calculate_dxyz = true; return jac; }
151 
156  const std::vector<Real> & get_JxW() const
157  { libmesh_assert(!calculations_started || calculate_dxyz);
158  calculate_dxyz = true; return JxW; }
159 
164  const std::vector<RealGradient> & get_dxyzdxi() const
165  { libmesh_assert(!calculations_started || calculate_dxyz);
166  calculate_dxyz = true; return dxyzdxi_map; }
167 
172  const std::vector<RealGradient> & get_dxyzdeta() const
173  { libmesh_assert(!calculations_started || calculate_dxyz);
174  calculate_dxyz = true; return dxyzdeta_map; }
175 
180  const std::vector<RealGradient> & get_dxyzdzeta() const
181  { libmesh_assert(!calculations_started || calculate_dxyz);
182  calculate_dxyz = true; return dxyzdzeta_map; }
183 
187  const std::vector<RealGradient> & get_d2xyzdxi2() const
188  { libmesh_assert(!calculations_started || calculate_d2xyz);
189  calculate_d2xyz = true; return d2xyzdxi2_map; }
190 
194  const std::vector<RealGradient> & get_d2xyzdeta2() const
195  { libmesh_assert(!calculations_started || calculate_d2xyz);
196  calculate_d2xyz = true; return d2xyzdeta2_map; }
197 
198 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
199 
203  const std::vector<RealGradient> & get_d2xyzdzeta2() const
204  { libmesh_assert(!calculations_started || calculate_d2xyz);
205  calculate_d2xyz = true; return d2xyzdzeta2_map; }
206 
207 #endif
208 
212  const std::vector<RealGradient> & get_d2xyzdxideta() const
213  { libmesh_assert(!calculations_started || calculate_d2xyz);
214  calculate_d2xyz = true; return d2xyzdxideta_map; }
215 
216 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
217 
221  const std::vector<RealGradient> & get_d2xyzdxidzeta() const
222  { libmesh_assert(!calculations_started || calculate_d2xyz);
223  calculate_d2xyz = true; return d2xyzdxidzeta_map; }
224 
228  const std::vector<RealGradient> & get_d2xyzdetadzeta() const
229  { libmesh_assert(!calculations_started || calculate_d2xyz);
230  calculate_d2xyz = true; return d2xyzdetadzeta_map; }
231 
232 #endif
233 
238  const std::vector<Real> & get_dxidx() const
239  { libmesh_assert(!calculations_started || calculate_dxyz);
240  calculate_dxyz = true; return dxidx_map; }
241 
246  const std::vector<Real> & get_dxidy() const
247  { libmesh_assert(!calculations_started || calculate_dxyz);
248  calculate_dxyz = true; return dxidy_map; }
249 
254  const std::vector<Real> & get_dxidz() const
255  { libmesh_assert(!calculations_started || calculate_dxyz);
256  calculate_dxyz = true; return dxidz_map; }
257 
262  const std::vector<Real> & get_detadx() const
263  { libmesh_assert(!calculations_started || calculate_dxyz);
264  calculate_dxyz = true; return detadx_map; }
265 
270  const std::vector<Real> & get_detady() const
271  { libmesh_assert(!calculations_started || calculate_dxyz);
272  calculate_dxyz = true; return detady_map; }
273 
278  const std::vector<Real> & get_detadz() const
279  { libmesh_assert(!calculations_started || calculate_dxyz);
280  calculate_dxyz = true; return detadz_map; }
281 
286  const std::vector<Real> & get_dzetadx() const
287  { libmesh_assert(!calculations_started || calculate_dxyz);
288  calculate_dxyz = true; return dzetadx_map; }
289 
294  const std::vector<Real> & get_dzetady() const
295  { libmesh_assert(!calculations_started || calculate_dxyz);
296  calculate_dxyz = true; return dzetady_map; }
297 
302  const std::vector<Real> & get_dzetadz() const
303  { libmesh_assert(!calculations_started || calculate_dxyz);
304  calculate_dxyz = true; return dzetadz_map; }
305 
306 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
307 
310  const std::vector<std::vector<Real>> & get_d2xidxyz2() const
311  { libmesh_assert(!calculations_started || calculate_d2xyz);
312  calculate_d2xyz = true; return d2xidxyz2_map; }
313 
317  const std::vector<std::vector<Real>> & get_d2etadxyz2() const
318  { libmesh_assert(!calculations_started || calculate_d2xyz);
319  calculate_d2xyz = true; return d2etadxyz2_map; }
320 
324  const std::vector<std::vector<Real>> & get_d2zetadxyz2() const
325  { libmesh_assert(!calculations_started || calculate_d2xyz);
326  calculate_d2xyz = true; return d2zetadxyz2_map; }
327 #endif
328 
332  const std::vector<std::vector<Real>> & get_psi() const
333  { return psi_map; }
334 
338  const std::vector<std::vector<Real>> & get_phi_map() const
339  { libmesh_assert(!calculations_started || calculate_xyz);
340  calculate_xyz = true; return phi_map; }
341 
345  const std::vector<std::vector<Real>> & get_dphidxi_map() const
346  { libmesh_assert(!calculations_started || calculate_dxyz);
347  calculate_dxyz = true; return dphidxi_map; }
348 
352  const std::vector<std::vector<Real>> & get_dphideta_map() const
353  { libmesh_assert(!calculations_started || calculate_dxyz);
354  calculate_dxyz = true; return dphideta_map; }
355 
359  const std::vector<std::vector<Real>> & get_dphidzeta_map() const
360  { libmesh_assert(!calculations_started || calculate_dxyz);
361  calculate_dxyz = true; return dphidzeta_map; }
362 
366  const std::vector<std::vector<Point>> & get_tangents() const
367  { libmesh_assert(!calculations_started || calculate_dxyz);
368  calculate_dxyz = true; return tangents; }
369 
373  const std::vector<Point> & get_normals() const
374  { libmesh_assert(!calculations_started || calculate_dxyz);
375  calculate_dxyz = true; return normals; }
376 
380  const std::vector<Real> & get_curvatures() const
381  { libmesh_assert(!calculations_started || calculate_d2xyz);
382  calculate_d2xyz = true; return curvatures;}
383 
387  void print_JxW(std::ostream & os) const;
388 
393  void print_xyz(std::ostream & os) const;
394 
395  /* FIXME: PB: The public functions below break encapsulation! I'm not happy
396  about it, but I don't have time to redo the infinite element routines.
397  These are needed in inf_fe_boundary.C and inf_fe.C. A proper implementation
398  would subclass this class and hide all the radial function business. */
402  std::vector<std::vector<Real>> & get_psi()
403  { return psi_map; }
404 
408  std::vector<std::vector<Real>> & get_dpsidxi()
409  { libmesh_assert(!calculations_started || calculate_dxyz);
410  calculate_dxyz = true; return dpsidxi_map; }
411 
415  std::vector<std::vector<Real>> & get_dpsideta()
416  { libmesh_assert(!calculations_started || calculate_dxyz);
417  calculate_dxyz = true; return dpsideta_map; }
418 
422  std::vector<std::vector<Real>> & get_d2psidxi2()
423  { libmesh_assert(!calculations_started || calculate_d2xyz);
424  calculate_d2xyz = true; return d2psidxi2_map; }
425 
429  std::vector<std::vector<Real>> & get_d2psidxideta()
430  { libmesh_assert(!calculations_started || calculate_d2xyz);
431  calculate_d2xyz = true; return d2psidxideta_map; }
432 
436  std::vector<std::vector<Real>> & get_d2psideta2()
437  { libmesh_assert(!calculations_started || calculate_d2xyz);
438  calculate_d2xyz = true; return d2psideta2_map; }
439 
443  std::vector<std::vector<Real>> & get_phi_map()
444  { libmesh_assert(!calculations_started || calculate_xyz);
445  calculate_xyz = true; return phi_map; }
446 
450  std::vector<std::vector<Real>> & get_dphidxi_map()
451  { libmesh_assert(!calculations_started || calculate_dxyz);
452  calculate_dxyz = true; return dphidxi_map; }
453 
457  std::vector<std::vector<Real>> & get_dphideta_map()
458  { libmesh_assert(!calculations_started || calculate_dxyz);
459  calculate_dxyz = true; return dphideta_map; }
460 
464  std::vector<std::vector<Real>> & get_dphidzeta_map()
465  { libmesh_assert(!calculations_started || calculate_dxyz);
466  calculate_dxyz = true; return dphidzeta_map; }
467 
468 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
469 
472  std::vector<std::vector<Real>> & get_d2phidxi2_map()
473  { libmesh_assert(!calculations_started || calculate_d2xyz);
474  calculate_d2xyz = true; return d2phidxi2_map; }
475 
479  std::vector<std::vector<Real>> & get_d2phidxideta_map()
480  { libmesh_assert(!calculations_started || calculate_d2xyz);
481  calculate_d2xyz = true; return d2phidxideta_map; }
482 
486  std::vector<std::vector<Real>> & get_d2phidxidzeta_map()
487  { libmesh_assert(!calculations_started || calculate_d2xyz);
488  calculate_d2xyz = true; return d2phidxidzeta_map; }
489 
493  std::vector<std::vector<Real>> & get_d2phideta2_map()
494  { libmesh_assert(!calculations_started || calculate_d2xyz);
495  calculate_d2xyz = true; return d2phideta2_map; }
496 
500  std::vector<std::vector<Real>> & get_d2phidetadzeta_map()
501  { libmesh_assert(!calculations_started || calculate_d2xyz);
502  calculate_d2xyz = true; return d2phidetadzeta_map; }
503 
507  std::vector<std::vector<Real>> & get_d2phidzeta2_map()
508  { libmesh_assert(!calculations_started || calculate_d2xyz);
509  calculate_d2xyz = true; return d2phidzeta2_map; }
510 #endif
511 
512  /* FIXME: PB: This function breaks encapsulation! Needed in FE<>::reinit and
513  InfFE<>::reinit. Not sure yet if the algorithm can be redone to avoid this. */
518  std::vector<Real> & get_JxW()
519  { libmesh_assert(!calculations_started || calculate_dxyz);
520  calculate_dxyz = true; return JxW; }
521 
522 protected:
523 
528  calculations_started = true;
529 
530 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
531  // Second derivative calculations currently have first derivative
532  // calculations as a prerequisite
533  if (calculate_d2xyz)
534  calculate_dxyz = true;
535 #endif
536  }
537 
541  void resize_quadrature_map_vectors(const unsigned int dim, unsigned int n_qp);
542 
549  Real dxdxi_map(const unsigned int p) const { return dxyzdxi_map[p](0); }
550 
557  Real dydxi_map(const unsigned int p) const { return dxyzdxi_map[p](1); }
558 
565  Real dzdxi_map(const unsigned int p) const { return dxyzdxi_map[p](2); }
566 
573  Real dxdeta_map(const unsigned int p) const { return dxyzdeta_map[p](0); }
574 
581  Real dydeta_map(const unsigned int p) const { return dxyzdeta_map[p](1); }
582 
589  Real dzdeta_map(const unsigned int p) const { return dxyzdeta_map[p](2); }
590 
597  Real dxdzeta_map(const unsigned int p) const { return dxyzdzeta_map[p](0); }
598 
605  Real dydzeta_map(const unsigned int p) const { return dxyzdzeta_map[p](1); }
606 
613  Real dzdzeta_map(const unsigned int p) const { return dxyzdzeta_map[p](2); }
614 
618  std::vector<Point> xyz;
619 
624  std::vector<RealGradient> dxyzdxi_map;
625 
630  std::vector<RealGradient> dxyzdeta_map;
631 
636  std::vector<RealGradient> dxyzdzeta_map;
637 
642  std::vector<RealGradient> d2xyzdxi2_map;
643 
648  std::vector<RealGradient> d2xyzdxideta_map;
649 
654  std::vector<RealGradient> d2xyzdeta2_map;
655 
656 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
657 
662  std::vector<RealGradient> d2xyzdxidzeta_map;
663 
668  std::vector<RealGradient> d2xyzdetadzeta_map;
669 
674  std::vector<RealGradient> d2xyzdzeta2_map;
675 
676 #endif
677 
682  std::vector<Real> dxidx_map;
683 
688  std::vector<Real> dxidy_map;
689 
694  std::vector<Real> dxidz_map;
695 
696 
701  std::vector<Real> detadx_map;
702 
707  std::vector<Real> detady_map;
708 
713  std::vector<Real> detadz_map;
714 
715 
720  std::vector<Real> dzetadx_map;
721 
726  std::vector<Real> dzetady_map;
727 
732  std::vector<Real> dzetadz_map;
733 
734 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
735 
739  std::vector<std::vector<Real>> d2xidxyz2_map;
740 
745  std::vector<std::vector<Real>> d2etadxyz2_map;
746 
751  std::vector<std::vector<Real>> d2zetadxyz2_map;
752 #endif
753 
757  std::vector<std::vector<Real>> phi_map;
758 
762  std::vector<std::vector<Real>> dphidxi_map;
763 
767  std::vector<std::vector<Real>> dphideta_map;
768 
772  std::vector<std::vector<Real>> dphidzeta_map;
773 
774 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
775 
779  std::vector<std::vector<Real>> d2phidxi2_map;
780 
784  std::vector<std::vector<Real>> d2phidxideta_map;
785 
789  std::vector<std::vector<Real>> d2phidxidzeta_map;
790 
794  std::vector<std::vector<Real>> d2phideta2_map;
795 
799  std::vector<std::vector<Real>> d2phidetadzeta_map;
800 
804  std::vector<std::vector<Real>> d2phidzeta2_map;
805 
806 #endif
807 
811  std::vector<std::vector<Real>> psi_map;
812 
817  std::vector<std::vector<Real>> dpsidxi_map;
818 
823  std::vector<std::vector<Real>> dpsideta_map;
824 
830  std::vector<std::vector<Real>> d2psidxi2_map;
831 
837  std::vector<std::vector<Real>> d2psidxideta_map;
838 
844  std::vector<std::vector<Real>> d2psideta2_map;
845 
849  std::vector<std::vector<Point>> tangents;
850 
854  std::vector<Point> normals;
855 
861  std::vector<Real> curvatures;
862 
866  std::vector<Real> jac;
867 
871  std::vector<Real> JxW;
872 
877  mutable bool calculations_started;
878 
882  mutable bool calculate_xyz;
883 
887  mutable bool calculate_dxyz;
888 
892  mutable bool calculate_d2xyz;
893 
898  template <unsigned int Dim, FEFamily T>
899  friend class FE;
900 
901 private:
906  void compute_inverse_map_second_derivs(unsigned p);
907 
911  std::vector<const Node *> _elem_nodes;
912 };
913 
914 }
915 
916 #endif // LIBMESH_FE_MAP_H
Manages the family, order, etc. parameters for a given FE.
Definition: fe_type.h:179
virtual void compute_affine_map(const unsigned int dim, const std::vector< Real > &qw, const Elem *elem)
Definition: fe_map.C:1232
std::vector< std::vector< Real > > & get_dphidzeta_map()
Definition: fe_map.h:464
void init_edge_shape_functions(const std::vector< Point > &qp, const Elem *edge)
Definition: fe_boundary.C:514
Real dzdeta_map(const unsigned int p) const
Definition: fe_map.h:589
std::vector< std::vector< Real > > dpsideta_map
Definition: fe_map.h:823
std::vector< std::vector< Real > > psi_map
Definition: fe_map.h:811
std::vector< std::vector< Real > > d2psideta2_map
Definition: fe_map.h:844
std::vector< std::vector< Real > > & get_d2phidxideta_map()
Definition: fe_map.h:479
std::vector< std::vector< Real > > & get_d2phideta2_map()
Definition: fe_map.h:493
std::vector< std::vector< Real > > dphidzeta_map
Definition: fe_map.h:772
void compute_inverse_map_second_derivs(unsigned p)
Definition: fe_map.C:1454
virtual void compute_map(const unsigned int dim, const std::vector< Real > &qw, const Elem *elem, bool calculate_d2phi)
Definition: fe_map.C:1388
std::vector< std::vector< Real > > dphidxi_map
Definition: fe_map.h:762
std::vector< std::vector< Real > > d2etadxyz2_map
Definition: fe_map.h:745
bool calculate_dxyz
Definition: fe_map.h:887
void compute_edge_map(int dim, const std::vector< Real > &qw, const Elem *side)
Definition: fe_boundary.C:922
std::vector< std::vector< Real > > & get_d2phidxi2_map()
Definition: fe_map.h:472
static std::unique_ptr< FEMap > build(FEType fe_type)
Definition: fe_map.C:53
std::vector< RealGradient > d2xyzdzeta2_map
Definition: fe_map.h:674
std::vector< Real > dzetady_map
Definition: fe_map.h:726
const std::vector< std::vector< Real > > & get_dphidzeta_map() const
Definition: fe_map.h:359
std::vector< std::vector< Real > > d2xidxyz2_map
Definition: fe_map.h:739
void init_reference_to_physical_map(const std::vector< Point > &qp, const Elem *elem)
Definition: fe_map.C:68
Real dxdeta_map(const unsigned int p) const
Definition: fe_map.h:573
unsigned short int side
Definition: xdr_io.C:50
const std::vector< Real > & get_detadz() const
Definition: fe_map.h:278
Real dzdxi_map(const unsigned int p) const
Definition: fe_map.h:565
const std::vector< Real > & get_curvatures() const
Definition: fe_map.h:380
The base class for all geometric element types.
Definition: elem.h:100
std::vector< std::vector< Real > > d2phideta2_map
Definition: fe_map.h:794
const std::vector< Point > & get_normals() const
Definition: fe_map.h:373
std::vector< std::vector< Real > > d2phidxidzeta_map
Definition: fe_map.h:789
std::vector< Real > dxidz_map
Definition: fe_map.h:694
const std::vector< std::vector< Real > > & get_phi_map() const
Definition: fe_map.h:338
std::vector< Real > & get_JxW()
Definition: fe_map.h:518
std::vector< std::vector< Real > > phi_map
Definition: fe_map.h:757
const std::vector< RealGradient > & get_dxyzdzeta() const
Definition: fe_map.h:180
const std::vector< std::vector< Real > > & get_d2etadxyz2() const
Definition: fe_map.h:317
std::vector< std::vector< Real > > d2phidetadzeta_map
Definition: fe_map.h:799
const std::vector< Real > & get_dxidz() const
Definition: fe_map.h:254
std::vector< std::vector< Real > > d2psidxideta_map
Definition: fe_map.h:837
std::vector< RealGradient > d2xyzdxideta_map
Definition: fe_map.h:648
const std::vector< std::vector< Real > > & get_dphideta_map() const
Definition: fe_map.h:352
const std::vector< std::vector< Point > > & get_tangents() const
Definition: fe_map.h:366
std::vector< RealGradient > dxyzdzeta_map
Definition: fe_map.h:636
std::vector< std::vector< Real > > & get_dpsidxi()
Definition: fe_map.h:408
std::vector< std::vector< Real > > & get_d2phidxidzeta_map()
Definition: fe_map.h:486
std::vector< Real > dzetadx_map
Definition: fe_map.h:720
std::vector< RealGradient > dxyzdxi_map
Definition: fe_map.h:624
Template class which generates the different FE families and orders.
Definition: fe.h:89
std::vector< std::vector< Real > > & get_d2phidetadzeta_map()
Definition: fe_map.h:500
const std::vector< RealGradient > & get_d2xyzdeta2() const
Definition: fe_map.h:194
std::vector< std::vector< Real > > & get_d2psideta2()
Definition: fe_map.h:436
std::vector< RealGradient > d2xyzdeta2_map
Definition: fe_map.h:654
const std::vector< std::vector< Real > > & get_d2zetadxyz2() const
Definition: fe_map.h:324
const std::vector< RealGradient > & get_dxyzdxi() const
Definition: fe_map.h:164
std::vector< std::vector< Real > > & get_d2psidxideta()
Definition: fe_map.h:429
std::vector< RealGradient > d2xyzdxi2_map
Definition: fe_map.h:642
std::vector< std::vector< Real > > & get_dphidxi_map()
Definition: fe_map.h:450
std::vector< Real > dzetadz_map
Definition: fe_map.h:732
const std::vector< std::vector< Real > > & get_psi() const
Definition: fe_map.h:332
Real dzdzeta_map(const unsigned int p) const
Definition: fe_map.h:613
bool calculate_d2xyz
Definition: fe_map.h:892
std::vector< std::vector< Real > > d2zetadxyz2_map
Definition: fe_map.h:751
const std::vector< Real > & get_dzetadx() const
Definition: fe_map.h:286
bool calculations_started
Definition: fe_map.h:877
std::vector< RealGradient > d2xyzdetadzeta_map
Definition: fe_map.h:668
std::vector< Real > dxidx_map
Definition: fe_map.h:682
bool calculate_xyz
Definition: fe_map.h:882
Real dxdzeta_map(const unsigned int p) const
Definition: fe_map.h:597
std::vector< Real > dxidy_map
Definition: fe_map.h:688
const std::vector< Real > & get_jacobian() const
Definition: fe_map.h:148
Real dydzeta_map(const unsigned int p) const
Definition: fe_map.h:605
const std::vector< Real > & get_dxidx() const
Definition: fe_map.h:238
std::vector< Real > curvatures
Definition: fe_map.h:861
const std::vector< RealGradient > & get_d2xyzdxi2() const
Definition: fe_map.h:187
const std::vector< Real > & get_dzetady() const
Definition: fe_map.h:294
const std::vector< Real > & get_JxW() const
Definition: fe_map.h:156
std::vector< std::vector< Real > > d2psidxi2_map
Definition: fe_map.h:830
const std::vector< Real > & get_dxidy() const
Definition: fe_map.h:246
std::vector< const Node * > _elem_nodes
Definition: fe_map.h:911
std::vector< std::vector< Real > > & get_d2psidxi2()
Definition: fe_map.h:422
std::vector< Real > JxW
Definition: fe_map.h:871
std::vector< std::vector< Real > > d2phidxideta_map
Definition: fe_map.h:784
const std::vector< Point > & get_xyz() const
Definition: fe_map.h:141
const std::vector< RealGradient > & get_d2xyzdetadzeta() const
Definition: fe_map.h:228
std::vector< Point > normals
Definition: fe_map.h:854
std::vector< Point > xyz
Definition: fe_map.h:618
const std::vector< Real > & get_dzetadz() const
Definition: fe_map.h:302
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real dydeta_map(const unsigned int p) const
Definition: fe_map.h:581
std::vector< std::vector< Real > > & get_d2phidzeta2_map()
Definition: fe_map.h:507
const std::vector< Real > & get_detady() const
Definition: fe_map.h:270
virtual void compute_null_map(const unsigned int dim, const std::vector< Real > &qw)
Definition: fe_map.C:1313
std::vector< std::vector< Real > > dpsidxi_map
Definition: fe_map.h:817
std::vector< Real > detady_map
Definition: fe_map.h:707
const std::vector< RealGradient > & get_d2xyzdzeta2() const
Definition: fe_map.h:203
const std::vector< std::vector< Real > > & get_dphidxi_map() const
Definition: fe_map.h:345
std::vector< std::vector< Point > > tangents
Definition: fe_map.h:849
void print_xyz(std::ostream &os) const
Definition: fe_map.C:1446
std::vector< std::vector< Real > > & get_dpsideta()
Definition: fe_map.h:415
std::vector< RealGradient > dxyzdeta_map
Definition: fe_map.h:630
void resize_quadrature_map_vectors(const unsigned int dim, unsigned int n_qp)
Definition: fe_map.C:1151
std::vector< std::vector< Real > > & get_psi()
Definition: fe_map.h:402
Real dxdxi_map(const unsigned int p) const
Definition: fe_map.h:549
const std::vector< RealGradient > & get_dxyzdeta() const
Definition: fe_map.h:172
void determine_calculations()
Definition: fe_map.h:527
const std::vector< RealGradient > & get_d2xyzdxidzeta() const
Definition: fe_map.h:221
Real dydxi_map(const unsigned int p) const
Definition: fe_map.h:557
virtual void compute_face_map(int dim, const std::vector< Real > &qw, const Elem *side)
Definition: fe_boundary.C:574
std::vector< Real > jac
Definition: fe_map.h:866
std::vector< Real > detadz_map
Definition: fe_map.h:713
std::vector< std::vector< Real > > d2phidzeta2_map
Definition: fe_map.h:804
std::vector< std::vector< Real > > & get_phi_map()
Definition: fe_map.h:443
void print_JxW(std::ostream &os) const
Definition: fe_map.C:1438
const std::vector< std::vector< Real > > & get_d2xidxyz2() const
Definition: fe_map.h:310
std::vector< std::vector< Real > > d2phidxi2_map
Definition: fe_map.h:779
Computes finite element mapping function values, gradients, etc.
Definition: fe_map.h:48
const std::vector< RealGradient > & get_d2xyzdxideta() const
Definition: fe_map.h:212
std::vector< Real > detadx_map
Definition: fe_map.h:701
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)
Definition: fe_map.C:412
std::vector< std::vector< Real > > dphideta_map
Definition: fe_map.h:767
std::vector< std::vector< Real > > & get_dphideta_map()
Definition: fe_map.h:457
void init_face_shape_functions(const std::vector< Point > &qp, const Elem *side)
Definition: fe_boundary.C:409
virtual ~FEMap()
Definition: fe_map.h:53
const std::vector< Real > & get_detadx() const
Definition: fe_map.h:262
std::vector< RealGradient > d2xyzdxidzeta_map
Definition: fe_map.h:662