MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (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-2020 Manav Bhatia and MAST authors
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 
49  const MAST::GeomElem& elem,
50  const MAST::FlightCondition& f);
51 
53 
54 
58  virtual bool
59  internal_residual (bool request_jacobian,
60  RealVectorX& f,
61  RealMatrixX& jac);
62 
63 
68  virtual bool
69  linearized_internal_residual (bool request_jacobian,
70  RealVectorX& f,
71  RealMatrixX& jac);
72 
73 
77  virtual bool
78  velocity_residual (bool request_jacobian,
79  RealVectorX& f,
80  RealMatrixX& jac_xdot,
81  RealMatrixX& jac);
82 
87  virtual bool
88  linearized_velocity_residual (bool request_jacobian,
89  RealVectorX& f,
90  RealMatrixX& jac_xdot,
91  RealMatrixX& jac);
92 
96  bool
97  side_external_residual (bool request_jacobian,
98  RealVectorX& f,
99  RealMatrixX& jac,
100  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
101 
102 
106  bool
107  linearized_side_external_residual (bool request_jacobian,
108  RealVectorX& f,
109  RealMatrixX& jac,
110  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
111 
115  virtual bool
117  bool request_jacobian,
118  RealVectorX& f,
119  RealMatrixX& jac);
123  virtual bool
125  bool request_jacobian,
126  RealVectorX& f,
127  RealMatrixX& jac_xdot,
128  RealMatrixX& jac);
129 
133  bool
135  bool request_jacobian,
136  RealVectorX& f,
137  RealMatrixX& jac,
138  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
139 
143  virtual void
144  side_integrated_force(const unsigned int s,
145  RealVectorX& f,
146  RealMatrixX* dfdX = nullptr);
147 
148 
149  virtual void
151  const unsigned int s,
152  RealVectorX& f);
153 
154 
158  virtual bool symmetry_surface_residual(bool request_jacobian,
159  RealVectorX& f,
160  RealMatrixX& jac,
161  const unsigned int s,
163 
168  bool request_jacobian,
169  RealVectorX& f,
170  RealMatrixX& jac,
171  const unsigned int s,
173 
177  virtual bool slip_wall_surface_residual(bool request_jacobian,
178  RealVectorX& f,
179  RealMatrixX& jac,
180  const unsigned int s,
182 
183 
187  virtual bool
188  linearized_slip_wall_surface_residual(bool request_jacobian,
189  RealVectorX& f,
190  RealMatrixX& jac,
191  const unsigned int s,
193 
194 
199  bool request_jacobian,
200  RealVectorX& f,
201  RealMatrixX& jac,
202  const unsigned int s,
204 
208  virtual bool noslip_wall_surface_residual(bool request_jacobian,
209  RealVectorX& f,
210  RealMatrixX& jac,
211  const unsigned int s,
213 
214 
219  bool request_jacobian,
220  RealVectorX& f,
221  RealMatrixX& jac,
222  const unsigned int s,
224 
228  virtual bool far_field_surface_residual(bool request_jacobian,
229  RealVectorX& f,
230  RealMatrixX& jac,
231  const unsigned int s,
233 
238  bool request_jacobian,
239  RealVectorX& f,
240  RealMatrixX& jac,
241  const unsigned int s,
243 
244  protected:
245 
246 
250  void _calculate_surface_integrated_load(bool request_derivative,
251  const MAST::FunctionBase* p,
252  const unsigned int s,
254 
255 
258  void _initialize_fem_interpolation_operator(const unsigned int qp,
259  const unsigned int dim,
260  const MAST::FEBase& fe,
262 
263 
273  void _initialize_fem_gradient_operator(const unsigned int qp,
274  const unsigned int dim,
275  const MAST::FEBase& fe,
276  std::vector<MAST::FEMOperatorMatrix>& dBmat);
277 
281  void
282  _initialize_fem_second_derivative_operator(const unsigned int qp,
283  const unsigned int dim,
284  const MAST::FEBase& fe,
285  std::vector<std::vector<MAST::FEMOperatorMatrix>>& d2Bmat);
286 
287  };
288 }
289 
290 
291 #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)
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.
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 velocity_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial 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)
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
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
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
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
const MAST::GeomElem & elem() const
Definition: elem_base.h:109
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)
ConservativeFluidElementBase(MAST::SystemInitialization &sys, const MAST::GeomElem &elem, const MAST::FlightCondition &f)
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