MAST
heat_conduction_elem_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__heat_conduction_elem_base__
22 #define __mast__heat_conduction_elem_base__
23 
24 // MAST includes
25 #include "base/elem_base.h"
26 
27 namespace MAST {
28 
29  // Forward declerations
30  class ElementPropertyCardBase;
31  class BoundaryConditionBase;
32  class FEMOperatorMatrix;
33  template <typename ValType> class FieldFunction;
34 
35 
54  public MAST::ElementBase
55  {
56  public:
62  const MAST::GeomElem& elem,
64 
65 
67 
68 
74  return _property;
75  }
76 
77 
81  virtual void
82  internal_residual (bool request_jacobian,
83  RealVectorX& f,
84  RealMatrixX& jac);
85 
86 
90  virtual void
91  velocity_residual (bool request_jacobian,
92  RealVectorX& f,
93  RealMatrixX& jac_xdot,
94  RealMatrixX& jac);
95 
99  void
100  side_external_residual (bool request_jacobian,
101  RealVectorX& f,
102  RealMatrixX& jac,
103  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
104 
108  void
109  volume_external_residual (bool request_jacobian,
110  RealVectorX& f,
111  RealMatrixX& jac,
112  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
113 
117  virtual void
119  RealVectorX& f);
120 
124  virtual void
126  RealVectorX& f,
127  const unsigned int s,
128  const MAST::FieldFunction<RealVectorX>& vel_f);
129 
133  virtual void
135  RealVectorX& f);
136 
140  virtual void
142  RealVectorX& f,
143  const unsigned int s,
144  const MAST::FieldFunction<RealVectorX>& vel_f);
145 
149  void
151  RealVectorX& f,
152  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
153 
157  void
159  RealVectorX& f,
160  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
161 
165  void
167  RealVectorX& f,
168  const unsigned int s,
170  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
171 
172  protected:
173 
174 
175 
180  virtual void surface_flux_residual(bool request_jacobian,
181  RealVectorX& f,
182  RealMatrixX& jac,
183  const unsigned int s,
185 
191  virtual void surface_flux_residual(bool request_jacobian,
192  RealVectorX& f,
193  RealMatrixX& jac,
195 
201  RealVectorX& f,
202  const unsigned int s,
204 
211  RealVectorX& f,
213 
220  RealVectorX& f,
221  const unsigned int s,
224 
225 
230  virtual void surface_convection_residual(bool request_jacobian,
231  RealVectorX& f,
232  RealMatrixX& jac,
233  const unsigned int s,
235 
236 
242  virtual void surface_convection_residual(bool request_jacobian,
243  RealVectorX& f,
244  RealMatrixX& jac,
246 
247 
252  virtual void
254  RealVectorX& f,
255  const unsigned int s,
257 
263  virtual void
265  RealVectorX& f,
267 
273  virtual void
275  RealVectorX& f,
276  const unsigned int s,
279 
284  virtual void surface_radiation_residual(bool request_jacobian,
285  RealVectorX& f,
286  RealMatrixX& jac,
287  const unsigned int s,
289 
290 
295  virtual void surface_radiation_residual(bool request_jacobian,
296  RealVectorX& f,
297  RealMatrixX& jac,
299 
300 
305  virtual void
307  RealVectorX& f,
308  const unsigned int s,
310 
311 
316  virtual void
318  RealVectorX& f,
320 
326  virtual void
328  RealVectorX& f,
329  const unsigned int s,
332 
337  virtual void volume_heat_source_residual(bool request_jacobian,
338  RealVectorX& f,
339  RealMatrixX& jac,
341 
342 
347  virtual void
349  RealVectorX& f,
351 
356  virtual void
358  RealVectorX& f,
359  const unsigned int s,
362 
368  void _initialize_mass_fem_operator(const unsigned int qp,
369  const MAST::FEBase& fe,
371 
372 
380  void _initialize_fem_gradient_operator(const unsigned int qp,
381  const unsigned int dim,
382  const MAST::FEBase& fe,
383  std::vector<MAST::FEMOperatorMatrix>& dBmat);
384 
389  };
390 
391 }
392 
393 
394 #endif // __mast__heat_conduction_elem_base__
virtual void surface_convection_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector sensitivity and Jacobian due to surface convection.
void _initialize_mass_fem_operator(const unsigned int qp, const MAST::FEBase &fe, MAST::FEMOperatorMatrix &Bmat)
When mass = false, initializes the FEM operator matrix to the shape functions as .
void 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
HeatConductionElementBase(MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
Constructor.
virtual void volume_heat_source_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to volume heat source.
virtual void internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual
virtual void surface_flux_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector sensitivity and Jacobian due to surface flux on element side s...
virtual void surface_flux_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface flux on element side s.
void side_external_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc)
sensitivity of the side external force contribution to system residual
const MAST::GeomElem & elem() const
Definition: elem_base.h:117
virtual void volume_heat_source_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to volume heat source on element volumetric domain...
void volume_external_residual_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > &bc)
boundary velocity contribution of volume external force.
MAST::AssemblyBase & assembly()
Definition: elem_base.h:103
virtual void internal_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f)
sensitivity of the internal force contribution to system residual
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 ...
virtual void surface_radiation_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain...
Matrix< Real, Dynamic, Dynamic > RealMatrixX
virtual void surface_convection_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface convection.
Matrix< Real, Dynamic, 1 > RealVectorX
virtual void surface_convection_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain...
virtual void volume_heat_source_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to volume heat source.
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
const MAST::ElementPropertyCardBase & elem_property()
returns a constant reference to the finite element object
void volume_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > &bc)
volume external force contribution to system residual
virtual void surface_flux_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain...
This element implements the Galerkin discretization of the heat conduction problem with the flux pro...
virtual void velocity_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual
virtual void surface_radiation_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector sensitivity and Jacobian due to surface radiation flux on element side...
virtual void velocity_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f)
sensitivity of the damping force contribution to system residual
virtual void internal_residual_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f)
sensitivity of the internal force contribution to system residual
virtual void surface_radiation_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int s, MAST::BoundaryConditionBase &bc)
Calculates the residual vector and Jacobian due to surface radiation flux on side s...
const MAST::ElementPropertyCardBase & _property
element property
void volume_external_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > &bc)
sensitivity of the volume external force contribution to system residual
virtual void velocity_residual_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f)
sensitivity of the internal force contribution to system residual
This is the base class for elements that implement calculation of finite element quantities over the ...
Definition: elem_base.h:72