MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (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-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 
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:
61  const MAST::GeomElem& elem,
63 
64 
66 
67 
73  return _property;
74  }
75 
76 
80  virtual void
81  internal_residual (bool request_jacobian,
82  RealVectorX& f,
83  RealMatrixX& jac);
84 
85 
89  virtual void
90  velocity_residual (bool request_jacobian,
91  RealVectorX& f,
92  RealMatrixX& jac_xdot,
93  RealMatrixX& jac);
94 
98  void
99  side_external_residual (bool request_jacobian,
100  RealVectorX& f,
101  RealMatrixX& jac,
102  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
103 
107  void
108  volume_external_residual (bool request_jacobian,
109  RealVectorX& f,
110  RealMatrixX& jac,
111  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
112 
116  virtual void
118  RealVectorX& f);
119 
123  virtual void
125  RealVectorX& f,
126  const unsigned int s,
127  const MAST::FieldFunction<RealVectorX>& vel_f);
128 
132  virtual void
134  RealVectorX& f);
135 
139  virtual void
141  RealVectorX& f,
142  const unsigned int s,
143  const MAST::FieldFunction<RealVectorX>& vel_f);
144 
148  void
150  RealVectorX& f,
151  std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
152 
156  void
158  RealVectorX& f,
159  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
160 
164  void
166  RealVectorX& f,
167  const unsigned int s,
169  std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
170 
171  protected:
172 
173 
174 
179  virtual void surface_flux_residual(bool request_jacobian,
180  RealVectorX& f,
181  RealMatrixX& jac,
182  const unsigned int s,
184 
190  virtual void surface_flux_residual(bool request_jacobian,
191  RealVectorX& f,
192  RealMatrixX& jac,
194 
200  RealVectorX& f,
201  const unsigned int s,
203 
210  RealVectorX& f,
212 
219  RealVectorX& f,
220  const unsigned int s,
223 
224 
229  virtual void surface_convection_residual(bool request_jacobian,
230  RealVectorX& f,
231  RealMatrixX& jac,
232  const unsigned int s,
234 
235 
241  virtual void surface_convection_residual(bool request_jacobian,
242  RealVectorX& f,
243  RealMatrixX& jac,
245 
246 
251  virtual void
253  RealVectorX& f,
254  const unsigned int s,
256 
262  virtual void
264  RealVectorX& f,
266 
272  virtual void
274  RealVectorX& f,
275  const unsigned int s,
278 
283  virtual void surface_radiation_residual(bool request_jacobian,
284  RealVectorX& f,
285  RealMatrixX& jac,
286  const unsigned int s,
288 
289 
294  virtual void surface_radiation_residual(bool request_jacobian,
295  RealVectorX& f,
296  RealMatrixX& jac,
298 
299 
304  virtual void
306  RealVectorX& f,
307  const unsigned int s,
309 
310 
315  virtual void
317  RealVectorX& f,
319 
325  virtual void
327  RealVectorX& f,
328  const unsigned int s,
331 
336  virtual void volume_heat_source_residual(bool request_jacobian,
337  RealVectorX& f,
338  RealMatrixX& jac,
340 
341 
346  virtual void
348  RealVectorX& f,
350 
355  virtual void
357  RealVectorX& f,
358  const unsigned int s,
361 
367  void _initialize_mass_fem_operator(const unsigned int qp,
368  const MAST::FEBase& fe,
370 
371 
379  void _initialize_fem_gradient_operator(const unsigned int qp,
380  const unsigned int dim,
381  const MAST::FEBase& fe,
382  std::vector<MAST::FEMOperatorMatrix>& dBmat);
383 
388  };
389 
390 }
391 
392 
393 #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 .
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.
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...
virtual void internal_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f)
sensitivity of the internal force contribution to system residual
HeatConductionElementBase(MAST::SystemInitialization &sys, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
Constructor.
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 ...
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
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.
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
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
const MAST::GeomElem & elem() const
Definition: elem_base.h:109
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
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
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.
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
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