MAST
fe_base.h
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2019 Manav Bhatia
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 
21 #ifndef __mast_fe_base_h__
22 #define __mast_fe_base_h__
23 
24 // MAST includes
25 #include "base/mast_data_types.h"
26 
27 // libMesh includes
28 #include "libmesh/elem.h"
29 #include "libmesh/fe_base.h"
30 #include "libmesh/quadrature.h"
31 
32 
33 namespace MAST {
34 
35  // forward declerations
36  class SystemInitialization;
37  class GeomElem;
38 
39  class FEBase {
40 
41  public:
42 
44 
45  virtual ~FEBase();
46 
51  void set_extra_quadrature_order(int n);
52 
57 
69  virtual void init(const MAST::GeomElem& elem,
70  bool init_grads,
71  const std::vector<libMesh::Point>* pts = nullptr);
72 
73 
81  virtual void init_for_side(const MAST::GeomElem& elem,
82  unsigned int s,
83  bool if_calculate_dphi);
84 
85 
86  libMesh::FEType
87  get_fe_type() const;
88 
89  virtual const std::vector<Real>&
90  get_JxW() const;
91 
96  virtual const std::vector<libMesh::Point>&
97  get_xyz() const;
98 
99  virtual unsigned int
100  n_shape_functions() const;
101 
102  virtual const std::vector<std::vector<Real> >&
103  get_phi() const;
104 
105  virtual const std::vector<std::vector<libMesh::RealVectorValue> >&
106  get_dphi() const;
107 
108  virtual const std::vector<std::vector<libMesh::RealTensorValue>>&
109  get_d2phi() const;
110 
111  virtual const std::vector<Real>&
112  get_dxidx() const;
113 
114  virtual const std::vector<Real>&
115  get_dxidy() const;
116 
117  virtual const std::vector<Real>&
118  get_dxidz() const;
119 
120  virtual const std::vector<Real>&
121  get_detadx() const;
122 
123  virtual const std::vector<Real>&
124  get_detady() const;
125 
126  virtual const std::vector<Real>&
127  get_detadz() const;
128 
129  virtual const std::vector<Real>&
130  get_dzetadx() const;
131 
132  virtual const std::vector<Real>&
133  get_dzetady() const;
134 
135  virtual const std::vector<Real>&
136  get_dzetadz() const;
137 
138  virtual const std::vector<libMesh::RealVectorValue>&
139  get_dxyzdxi() const;
140 
141  virtual const std::vector<libMesh::RealVectorValue>&
142  get_dxyzdeta() const;
143 
144  virtual const std::vector<libMesh::RealVectorValue>&
145  get_dxyzdzeta() const;
146 
147  virtual const std::vector<std::vector<Real> >&
148  get_dphidxi() const;
149 
150  virtual const std::vector<std::vector<Real> >&
151  get_dphideta() const;
152 
153  virtual const std::vector<std::vector<Real> >&
154  get_dphidzeta() const;
155 
160  virtual const std::vector<libMesh::Point>&
162 
167  virtual const std::vector<libMesh::Point>&
169 
170 
171  virtual const std::vector<libMesh::Point>&
172  get_qpoints() const;
173 
174  virtual const libMesh::QBase&
175  get_qrule() const;
176 
177  protected:
178 
184  libMesh::FEBase* _fe;
185  libMesh::QBase* _qrule;
186  std::vector<libMesh::Point> _qpoints;
187  std::vector<libMesh::Point> _global_xyz;
188  std::vector<libMesh::Point> _local_normals;
189  std::vector<libMesh::Point> _global_normals;
190  };
191 }
192 
193 
194 #endif // __mast_fe_base_h__
virtual const std::vector< Real > & get_detadz() const
Definition: fe_base.cpp:311
virtual ~FEBase()
Definition: fe_base.cpp:39
virtual const std::vector< Real > & get_dzetadx() const
Definition: fe_base.cpp:319
std::vector< libMesh::Point > _global_normals
Definition: fe_base.h:189
virtual const std::vector< libMesh::RealVectorValue > & get_dxyzdzeta() const
Definition: fe_base.cpp:359
libMesh::QBase * _qrule
Definition: fe_base.h:185
virtual const std::vector< Real > & get_JxW() const
Definition: fe_base.cpp:220
virtual const std::vector< Real > & get_dxidy() const
Definition: fe_base.cpp:279
virtual const std::vector< Real > & get_dxidz() const
Definition: fe_base.cpp:287
virtual const libMesh::QBase & get_qrule() const
Definition: fe_base.cpp:418
FEBase(const MAST::SystemInitialization &sys)
Definition: fe_base.cpp:27
virtual const std::vector< Real > & get_detadx() const
Definition: fe_base.cpp:295
virtual const std::vector< libMesh::Point > & get_xyz() const
physical location of the quadrature point in the global coordinate system for the reference element ...
Definition: fe_base.cpp:228
bool _init_second_order_derivatives
Definition: fe_base.h:181
virtual const std::vector< std::vector< Real > > & get_dphidzeta() const
Definition: fe_base.cpp:383
virtual const std::vector< libMesh::Point > & get_qpoints() const
Definition: fe_base.cpp:407
virtual const std::vector< libMesh::RealVectorValue > & get_dxyzdeta() const
Definition: fe_base.cpp:351
virtual const std::vector< Real > & get_dzetady() const
Definition: fe_base.cpp:327
std::vector< libMesh::Point > _local_normals
Definition: fe_base.h:188
std::vector< libMesh::Point > _qpoints
Definition: fe_base.h:186
bool _initialized
Definition: fe_base.h:182
virtual void init_for_side(const MAST::GeomElem &elem, unsigned int s, bool if_calculate_dphi)
Initializes the quadrature and finite element for element side integration.
Definition: fe_base.cpp:137
virtual const std::vector< Real > & get_dzetadz() const
Definition: fe_base.cpp:335
const MAST::SystemInitialization & _sys
Definition: fe_base.h:179
void set_evaluate_second_order_derivatives(bool f)
sets the flag for evaluation of second order derivative
Definition: fe_base.cpp:57
virtual const std::vector< std::vector< Real > > & get_dphideta() const
Definition: fe_base.cpp:375
libMesh::FEBase * _fe
Definition: fe_base.h:184
void set_extra_quadrature_order(int n)
this is used, in addition to libMesh::System::extra_quadrature_order to set the quadrature rule...
Definition: fe_base.cpp:50
virtual const std::vector< std::vector< Real > > & get_dphidxi() const
Definition: fe_base.cpp:367
virtual unsigned int n_shape_functions() const
Definition: fe_base.cpp:239
virtual const std::vector< libMesh::Point > & get_normals_for_reference_coordinate() const
normals defined in the global coordinate system for the reference element
Definition: fe_base.cpp:391
std::vector< libMesh::Point > _global_xyz
Definition: fe_base.h:187
virtual const std::vector< Real > & get_detady() const
Definition: fe_base.cpp:303
virtual const std::vector< Real > & get_dxidx() const
Definition: fe_base.cpp:271
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
virtual void init(const MAST::GeomElem &elem, bool init_grads, const std::vector< libMesh::Point > *pts=nullptr)
Initializes the quadrature and finite element for element volume integration.
Definition: fe_base.cpp:64
virtual const std::vector< std::vector< libMesh::RealVectorValue > > & get_dphi() const
Definition: fe_base.cpp:255
virtual const std::vector< libMesh::RealVectorValue > & get_dxyzdxi() const
Definition: fe_base.cpp:343
virtual const std::vector< std::vector< libMesh::RealTensorValue > > & get_d2phi() const
Definition: fe_base.cpp:263
virtual const std::vector< std::vector< Real > > & get_phi() const
Definition: fe_base.cpp:247
libMesh::FEType get_fe_type() const
Definition: fe_base.cpp:212
virtual const std::vector< libMesh::Point > & get_normals_for_local_coordinate() const
normals defined in the coordinate system for the local reference element.
Definition: fe_base.cpp:399
const MAST::GeomElem * _elem
Definition: fe_base.h:183
unsigned int _extra_quadrature_order
Definition: fe_base.h:180