This element implements the Galerkin discretization of the heat conduction problem
with the flux provided on the boundary with Neumann boundary conditions.
The discrete form is represented as
, where
, and
, where is the boundary normal surface flux on and is the volumetric heat generation.
Definition at line 53 of file heat_conduction_elem_base.h.
#include <heat_conduction_elem_base.h>
Public Member Functions | |
HeatConductionElementBase (MAST::SystemInitialization &sys, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p) | |
Constructor. More... | |
virtual | ~HeatConductionElementBase () |
const MAST::ElementPropertyCardBase & | elem_property () |
returns a constant reference to the finite element object More... | |
virtual void | internal_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac) |
internal force contribution to system residual More... | |
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 More... | |
virtual void | internal_residual_sensitivity (const MAST::FunctionBase &p, RealVectorX &f) |
sensitivity of the internal force contribution to system residual More... | |
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 More... | |
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 More... | |
virtual void | velocity_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac) |
inertial force contribution to system residual More... | |
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 More... | |
virtual void | velocity_residual_sensitivity (const MAST::FunctionBase &p, RealVectorX &f) |
sensitivity of the damping force contribution to system residual More... | |
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 More... | |
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. More... | |
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 More... | |
Public Member Functions inherited from MAST::ElementBase | |
ElementBase (MAST::SystemInitialization &sys, const MAST::GeomElem &elem) | |
The default constructor. More... | |
virtual | ~ElementBase () |
Default virtual destructor. More... | |
void | attach_active_solution_function (MAST::FunctionBase &f) |
Attaches the function that represents the system solution. More... | |
void | detach_active_solution_function () |
Detaches the function object that may have been attached to the element. More... | |
const MAST::GeomElem & | elem () const |
virtual void | set_acceleration (const RealVectorX &vec, bool if_sens=false) |
stores vec as acceleration for element level calculations, or its sensitivity if if_sens is true. More... | |
virtual void | set_complex_solution (const ComplexVectorX &vec, bool if_sens=false) |
This provides the complex solution (or its sensitivity if if_sens is true.) for frequecy-domain analysis. More... | |
virtual void | set_perturbed_acceleration (const RealVectorX &vec, bool if_sens=false) |
stores vec as perturbed acceleration for element level calculations, or its sensitivity if if_sens is true. More... | |
virtual void | set_perturbed_solution (const RealVectorX &vec, bool if_sens=false) |
This provides the perturbed solution (or its sensitivity if if_sens is true.) for linearized analysis. More... | |
virtual void | set_perturbed_velocity (const RealVectorX &vec, bool if_sens=false) |
stores vec as perturbed velocity for element level calculations, or its sensitivity if if_sens is true. More... | |
virtual void | set_solution (const RealVectorX &vec, bool if_sens=false) |
stores vec as solution for element level calculations, or its sensitivity if if_sens is true. More... | |
virtual void | set_velocity (const RealVectorX &vec, bool if_sens=false) |
stores vec as velocity for element level calculations, or its sensitivity if if_sens is true. More... | |
const RealVectorX & | sol (bool if_sens=false) const |
MAST::NonlinearSystem & | system () |
MAST::SystemInitialization & | system_initialization () |
Protected Member Functions | |
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
. More... | |
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
. More... | |
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. More... | |
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. More... | |
virtual void | surface_convection_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc) |
Calculates the residual vector and Jacobian due to surface convection on the element domain. More... | |
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. More... | |
virtual void | surface_convection_residual_sensitivity (const MAST::FunctionBase &p, RealVectorX &f, MAST::BoundaryConditionBase &bc) |
Calculates the residual vector sensitivity and Jacobian due to surface convection on element domain. More... | |
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. More... | |
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 . More... | |
virtual void | surface_flux_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc) |
Calculates the residual vector sensitivity due to surface flux on element volumetric domain. More... | |
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 . More... | |
virtual void | surface_flux_residual_sensitivity (const MAST::FunctionBase &p, RealVectorX &f, MAST::BoundaryConditionBase &bc) |
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain. More... | |
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. More... | |
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. More... | |
virtual void | surface_radiation_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc) |
Calculates the residual vector and Jacobian due to surface radiation flux on element domain. More... | |
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. More... | |
virtual void | surface_radiation_residual_sensitivity (const MAST::FunctionBase &p, RealVectorX &f, MAST::BoundaryConditionBase &bc) |
Calculates the residual vector sensitivity and Jacobian due to surface radiation flux on element domain. More... | |
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. More... | |
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. More... | |
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. More... | |
Protected Attributes | |
const MAST::ElementPropertyCardBase & | _property |
element property More... | |
Protected Attributes inherited from MAST::ElementBase | |
RealVectorX | _accel |
local acceleration More... | |
RealVectorX | _accel_sens |
local acceleration More... | |
MAST::FunctionBase * | _active_sol_function |
pointer to the active solution mesh field function. More... | |
ComplexVectorX | _complex_sol |
local solution used for frequency domain analysis More... | |
ComplexVectorX | _complex_sol_sens |
local solution used for frequency domain analysis More... | |
RealVectorX | _delta_accel |
local acceleration More... | |
RealVectorX | _delta_accel_sens |
local acceleration More... | |
RealVectorX | _delta_sol |
local solution used for linearized analysis More... | |
RealVectorX | _delta_sol_sens |
local solution used for linearized analysis More... | |
RealVectorX | _delta_vel |
local velocity More... | |
RealVectorX | _delta_vel_sens |
local velocity More... | |
const MAST::GeomElem & | _elem |
geometric element for which the computations are performed More... | |
RealVectorX | _sol |
local solution More... | |
RealVectorX | _sol_sens |
local solution sensitivity More... | |
MAST::SystemInitialization & | _system |
SystemInitialization object associated with this element. More... | |
const Real & | _time |
time for which system is being assembled More... | |
RealVectorX | _vel |
local velocity More... | |
RealVectorX | _vel_sens |
local velocity More... | |
MAST::HeatConductionElementBase::HeatConductionElementBase | ( | MAST::SystemInitialization & | sys, |
const MAST::GeomElem & | elem, | ||
const MAST::ElementPropertyCardBase & | p | ||
) |
Constructor.
Definition at line 36 of file heat_conduction_elem_base.cpp.
|
virtual |
Definition at line 46 of file heat_conduction_elem_base.cpp.
|
protected |
For mass
= true, the FEM operator matrix is initilized to the weak form of the Laplacian
.
Definition at line 1880 of file heat_conduction_elem_base.cpp.
|
protected |
When mass
= false, initializes the FEM operator matrix to the shape functions as
.
Definition at line 1857 of file heat_conduction_elem_base.cpp.
|
inline |
returns a constant reference to the finite element object
Definition at line 72 of file heat_conduction_elem_base.h.
|
virtual |
internal force contribution to system residual
Definition at line 54 of file heat_conduction_elem_base.cpp.
|
virtual |
sensitivity of the internal force contribution to system residual
Definition at line 622 of file heat_conduction_elem_base.cpp.
|
virtual |
sensitivity of the internal force contribution to system residual
Definition at line 554 of file heat_conduction_elem_base.cpp.
void MAST::HeatConductionElementBase::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
Definition at line 250 of file heat_conduction_elem_base.cpp.
void MAST::HeatConductionElementBase::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
Definition at line 369 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain.
This is used only for 1D or 2D elements.
Definition at line 1333 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector and Jacobian due to surface convection.
Definition at line 1101 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector and Jacobian due to surface convection on the element domain.
This is relevant only for 1D and 2D elements.
Definition at line 1162 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector sensitivity and Jacobian due to surface convection.
Definition at line 1218 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector sensitivity and Jacobian due to surface convection on element domain.
This is relevant only for 1D and 2D elements.
Definition at line 1282 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain.
This is used only for 1D or 2D elements.
Definition at line 1044 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector and Jacobian due to surface flux on element side s
.
Definition at line 877 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector sensitivity due to surface flux on element volumetric domain.
This is used only for 2D or 3D elements.
Definition at line 921 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector sensitivity and Jacobian due to surface flux on element side s
.
Definition at line 958 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain.
This is used only for 1D or 2D elements.
Definition at line 1009 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector and Jacobian due to surface flux on element volumetric domain.
This is used only for 1D or 2D elements.
Definition at line 1636 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector and Jacobian due to surface radiation flux on side s.
Definition at line 1404 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector and Jacobian due to surface radiation flux on element domain.
Definition at line 1467 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector sensitivity and Jacobian due to surface radiation flux on element side.
Definition at line 1525 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector sensitivity and Jacobian due to surface radiation flux on element domain.
Definition at line 1589 of file heat_conduction_elem_base.cpp.
|
virtual |
inertial force contribution to system residual
Definition at line 156 of file heat_conduction_elem_base.cpp.
|
virtual |
sensitivity of the internal force contribution to system residual
Definition at line 783 of file heat_conduction_elem_base.cpp.
|
virtual |
sensitivity of the damping force contribution to system residual
Definition at line 709 of file heat_conduction_elem_base.cpp.
void MAST::HeatConductionElementBase::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
Definition at line 312 of file heat_conduction_elem_base.cpp.
void MAST::HeatConductionElementBase::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.
Definition at line 486 of file heat_conduction_elem_base.cpp.
void MAST::HeatConductionElementBase::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
Definition at line 429 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector and Jacobian due to volume heat source on element volumetric domain.
Definition at line 1796 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector and Jacobian due to volume heat source.
Definition at line 1711 of file heat_conduction_elem_base.cpp.
|
protectedvirtual |
Calculates the residual vector and Jacobian due to volume heat source.
Definition at line 1751 of file heat_conduction_elem_base.cpp.
|
protected |
element property
Definition at line 387 of file heat_conduction_elem_base.h.