MAST
conservative_fluid_element_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 #ifndef __mast__conservative_fluid_element_base__
21 #define __mast__conservative_fluid_element_base__
22 
23 
24 // MAST includes
25 #include "base/elem_base.h"
26 #include "fluid/fluid_elem_base.h"
27 
28 
29 namespace MAST {
30 
31  // Forward declerations
32  class ElementPropertyCardBase;
33  class BoundaryConditionBase;
34  class FEMOperatorMatrix;
35  class OutputAssemblyElemOperations;
36 
37 
43  public MAST::FluidElemBase,
44  public MAST::ElementBase {
45 
46  public:
47 
50  const MAST::GeomElem& elem,
51  const MAST::FlightCondition& f);
52 
54 
55 
59  virtual bool
60  internal_residual (bool request_jacobian,
61  RealVectorX& f,
62  RealMatrixX& jac);
63 
64 
69  virtual bool
70  linearized_internal_residual (bool request_jacobian,
71  RealVectorX& f,
72  RealMatrixX& jac);
73 
74 
78  virtual bool
79  velocity_residual (bool request_jacobian,
80  RealVectorX& f,
81  RealMatrixX& jac_xdot,
82  RealMatrixX& jac);
83 
88  virtual bool
89  linearized_velocity_residual (bool request_jacobian,
90  RealVectorX& f,
91  RealMatrixX& jac_xdot,
92  RealMatrixX& jac);
93 
97  bool
98  side_external_residual (bool request_jacobian,
99  RealVectorX& f,
100  RealMatrixX& jac,
101  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
102 
103 
107  bool
108  linearized_side_external_residual (bool request_jacobian,
109  RealVectorX& f,
110  RealMatrixX& jac,
111  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
112 
116  virtual bool
118  bool request_jacobian,
119  RealVectorX& f,
120  RealMatrixX& jac);
124  virtual bool
126  bool request_jacobian,
127  RealVectorX& f,
128  RealMatrixX& jac_xdot,
129  RealMatrixX& jac);
130 
134  bool
136  bool request_jacobian,
137  RealVectorX& f,
138  RealMatrixX& jac,
139  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
140 
144  virtual void
145  side_integrated_force(const unsigned int s,
146  RealVectorX& f,
147  RealMatrixX* dfdX = nullptr);
148 
149 
150  virtual void
152  const unsigned int s,
153  RealVectorX& f);
154 
155 
159  virtual bool symmetry_surface_residual(bool request_jacobian,
160  RealVectorX& f,
161  RealMatrixX& jac,
162  const unsigned int s,
164 
169  bool request_jacobian,
170  RealVectorX& f,
171  RealMatrixX& jac,
172  const unsigned int s,
174 
178  virtual bool slip_wall_surface_residual(bool request_jacobian,
179  RealVectorX& f,
180  RealMatrixX& jac,
181  const unsigned int s,
183 
184 
188  virtual bool
189  linearized_slip_wall_surface_residual(bool request_jacobian,
190  RealVectorX& f,
191  RealMatrixX& jac,
192  const unsigned int s,
194 
195 
200  bool request_jacobian,
201  RealVectorX& f,
202  RealMatrixX& jac,
203  const unsigned int s,
205 
209  virtual bool noslip_wall_surface_residual(bool request_jacobian,
210  RealVectorX& f,
211  RealMatrixX& jac,
212  const unsigned int s,
214 
215 
220  bool request_jacobian,
221  RealVectorX& f,
222  RealMatrixX& jac,
223  const unsigned int s,
225 
229  virtual bool far_field_surface_residual(bool request_jacobian,
230  RealVectorX& f,
231  RealMatrixX& jac,
232  const unsigned int s,
234 
239  bool request_jacobian,
240  RealVectorX& f,
241  RealMatrixX& jac,
242  const unsigned int s,
244 
245  protected:
246 
247 
251  void _calculate_surface_integrated_load(bool request_derivative,
252  const MAST::FunctionBase* p,
253  const unsigned int s,
255 
256 
259  void _initialize_fem_interpolation_operator(const unsigned int qp,
260  const unsigned int dim,
261  const MAST::FEBase& fe,
263 
264 
274  void _initialize_fem_gradient_operator(const unsigned int qp,
275  const unsigned int dim,
276  const MAST::FEBase& fe,
277  std::vector<MAST::FEMOperatorMatrix>& dBmat);
278 
282  void
283  _initialize_fem_second_derivative_operator(const unsigned int qp,
284  const unsigned int dim,
285  const MAST::FEBase& fe,
286  std::vector<std::vector<MAST::FEMOperatorMatrix>>& d2Bmat);
287 
288  };
289 }
290 
291 
292 #endif // __mast__conservative_fluid_element_base__
void _initialize_fem_gradient_operator(const unsigned int qp, const unsigned int dim, const MAST::FEBase &fe, std::vector< MAST::FEMOperatorMatrix > &dBmat)
For mass = true, the FEM operator matrix is initilized to the weak form of the Laplacian ...
This class provides the necessary functionality for spatial discretization of the conservative fluid ...
virtual bool linearized_velocity_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual of the linearized problem
This class provides the necessary functions to evaluate the flux vectors and their Jacobians for both...
virtual bool noslip_wall_surface_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
virtual void side_integrated_force_sensitivity(const MAST::FunctionBase &p, const unsigned int s, RealVectorX &f)
const MAST::GeomElem & elem() const
Definition: elem_base.h:117
This provides the base class for definitin of element level contribution of output quantity in an ana...
virtual bool linearized_internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual of the linearized problem.
virtual bool velocity_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual
bool side_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc)
side external force contribution to system residual
ConservativeFluidElementBase(MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem, const MAST::FlightCondition &f)
MAST::AssemblyBase & assembly()
Definition: elem_base.h:103
bool side_external_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc)
sensitivity of the side external force contribution to system residual
virtual bool slip_wall_surface_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
virtual bool slip_wall_surface_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
virtual void side_integrated_force(const unsigned int s, RealVectorX &f, RealMatrixX *dfdX=nullptr)
surface integrated force
virtual bool symmetry_surface_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
virtual bool linearized_slip_wall_surface_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
virtual bool velocity_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
sensitivity of the damping force contribution to system residual
bool linearized_side_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc)
side external force contribution to system residual
Matrix< Real, Dynamic, Dynamic > RealMatrixX
void _calculate_surface_integrated_load(bool request_derivative, const MAST::FunctionBase *p, const unsigned int s, MAST::OutputAssemblyElemOperations &output)
calculates the surface integrated force vector
Matrix< Real, Dynamic, 1 > RealVectorX
virtual bool noslip_wall_surface_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
virtual bool internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
void _initialize_fem_interpolation_operator(const unsigned int qp, const unsigned int dim, const MAST::FEBase &fe, MAST::FEMOperatorMatrix &Bmat)
virtual bool internal_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
sensitivity of the internal force contribution to system residual
virtual bool far_field_surface_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
void _initialize_fem_second_derivative_operator(const unsigned int qp, const unsigned int dim, const MAST::FEBase &fe, std::vector< std::vector< MAST::FEMOperatorMatrix >> &d2Bmat)
d2Bmat[i][j] is the derivative d2B/dxi dxj
virtual bool symmetry_surface_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
virtual bool far_field_surface_residual_sensitivity(const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
This is the base class for elements that implement calculation of finite element quantities over the ...
Definition: elem_base.h:72