MAST
|
#include <structural_element_base.h>
Public Member Functions | |
StructuralElementBase (MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p) | |
Constructor. More... | |
virtual | ~StructuralElementBase () |
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_perturbed_solution (const RealVectorX &vec, bool if_sens=false) |
stores vec as perturbed 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... | |
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_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_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... | |
const RealVectorX & | local_solution (bool if_sens=false) const |
const MAST::ElementPropertyCardBase & | elem_property () |
returns a constant reference to the finite element object More... | |
virtual bool | internal_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0 |
internal force contribution to system residual More... | |
virtual bool | linearized_internal_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac) |
internal force contribution to system residual of the linearized problem More... | |
virtual void | internal_residual_boundary_velocity (const MAST::FunctionBase &p, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0 |
calculates the term on side s: . More... | |
virtual bool | internal_residual_jac_dot_state_sensitivity (RealMatrixX &jac)=0 |
calculates d[J]/d{x} . More... | |
virtual bool | damping_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac) |
damping force contribution to system residual More... | |
virtual bool | inertial_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac) |
inertial force contribution to system residual More... | |
virtual bool | linearized_inertial_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac) |
inertial force contribution to system residual of linerized problem More... | |
bool | side_external_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc) |
side external force contribution to system residual. More... | |
bool | linearized_side_external_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc) |
side external force contribution to system residual. More... | |
bool | linearized_frequency_domain_side_external_residual (bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc) |
Calculates the external force due to frequency domain side external force contribution to system residual. More... | |
virtual bool | prestress_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0 |
prestress force contribution to system residual More... | |
bool | volume_external_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, 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, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > &bc, bool request_jacobian, RealVectorX &f, RealMatrixX &jac) |
boundary velocity contribution of volume external force. More... | |
bool | linearized_volume_external_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > &bc) |
volume external force contribution to system residual. More... | |
bool | linearized_frequency_domain_volume_external_residual (bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > &bc) |
Calculates the frequency domain volume external force contribution to system residual. More... | |
virtual bool | internal_residual_sensitivity (const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0 |
sensitivity of the internal force contribution to system residual More... | |
virtual bool | damping_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 More... | |
virtual bool | inertial_residual_sensitivity (const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac) |
sensitivity of the inertial force contribution to system residual More... | |
virtual void | inertial_residual_boundary_velocity (const MAST::FunctionBase &p, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac) |
sensitivity of the inertial force contribution to system residual More... | |
bool | side_external_residual_sensitivity (const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > &bc) |
sensitivity of the side external force contribution to system residual More... | |
virtual bool | prestress_residual_sensitivity (const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0 |
sensitivity of the prestress force contribution to system residual More... | |
bool | volume_external_residual_sensitivity (const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > &bc) |
sensitivity of the volume external force contribution to system residual More... | |
virtual bool | if_incompatible_modes () const =0 |
virtual unsigned int | incompatible_mode_size () const |
void | set_incompatible_mode_solution (RealVectorX &vec) |
sets the pointer to the incompatible mode solution vector. More... | |
virtual void | update_incompatible_mode_solution (const RealVectorX &dsol) |
updates the incompatible solution for this element. More... | |
virtual bool | calculate_stress (bool request_derivative, const MAST::FunctionBase *f, MAST::StressStrainOutputBase &output)=0 |
Calculates the stress tensor. More... | |
virtual void | calculate_stress_boundary_velocity (const MAST::FunctionBase &p, MAST::StressStrainOutputBase &output, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f)=0 |
Calculates the boundary velocity term contributions to the sensitivity of stress at the specified boundary of this element. More... | |
virtual void | calculate_stress_temperature_derivative (MAST::FEBase &fe_thermal, MAST::StressStrainOutputBase &output)=0 |
virtual void | thermal_residual_temperature_derivative (const MAST::FEBase &fe_thermal, RealMatrixX &m)=0 |
template<typename ValType > | |
void | transform_matrix_to_global_system (const ValType &local_mat, ValType &global_mat) const |
template<typename ValType > | |
void | transform_vector_to_local_system (const ValType &global_vec, ValType &local_vec) const |
template<typename ValType > | |
void | transform_vector_to_global_system (const ValType &local_vec, ValType &global_vec) const |
Public Member Functions inherited from MAST::ElementBase | |
ElementBase (MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem) | |
The default constructor. More... | |
virtual | ~ElementBase () |
Default virtual destructor. More... | |
MAST::SystemInitialization & | system_initialization () |
MAST::AssemblyBase & | assembly () |
MAST::NonlinearSystem & | system () |
const MAST::GeomElem & | elem () const |
const RealVectorX & | sol (bool if_sens=false) const |
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... | |
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... | |
Public Attributes | |
bool | follower_forces |
flag for follower forces More... | |
Protected Member Functions | |
virtual bool | surface_pressure_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0 |
Calculates the force vector and Jacobian due to surface pressure. More... | |
virtual bool | surface_pressure_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc) |
Calculates the force vector and Jacobian due to surface pressure applied on the entire element domain. More... | |
virtual void | surface_pressure_boundary_velocity (const MAST::FunctionBase &p, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc, bool request_jacobian, RealVectorX &f, RealMatrixX &jac) |
Calculates the force vector and Jacobian due to surface pressure applied on the entire element domain. More... | |
virtual bool | linearized_surface_pressure_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc) |
Calculates the force vector and Jacobian due to surface pressure applied on the entire element domain. More... | |
Real | piston_theory_cp (const unsigned int order, const Real vel_U, const Real gamma, const Real mach) |
Real | piston_theory_dcp_dv (const unsigned int order, const Real vel_U, const Real gamma, const Real mach) |
virtual bool | piston_theory_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)=0 |
Calculates the force vector and Jacobian due to piston-theory based surface pressure on the entire element domain. More... | |
virtual bool | piston_theory_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0 |
Calculates the force vector and Jacobian due to piston-theory based surface pressure on the side identified by side . More... | |
virtual bool | piston_theory_residual_sensitivity (const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)=0 |
Calculates the force vector and Jacobian due to piston-theory based surface pressure on the entire element domain. More... | |
virtual bool | piston_theory_residual_sensitivity (const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0 |
Calculates the force vector and Jacobian due to piston-theory based surface pressure on the side identified by side . More... | |
virtual bool | linearized_frequency_domain_surface_pressure_residual (bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0 |
Calculates the force vector and Jacobian due to small perturbation surface pressure. More... | |
virtual bool | linearized_frequency_domain_surface_pressure_residual (bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, MAST::BoundaryConditionBase &bc) |
Calculates the force vector and Jacobian due to surface pressure applied on the entire element domain. More... | |
virtual bool | linearized_frequency_domain_surface_pressure_residual_sensitivity (const MAST::FunctionBase &p, bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0 |
Calculates the sensitivity of force vector and Jacobian due to small is applicable for perturbation surface pressure. More... | |
virtual bool | linearized_frequency_domain_surface_pressure_residual_sensitivity (const MAST::FunctionBase &p, bool request_jacobian, ComplexVectorX &f, ComplexMatrixX &jac, MAST::BoundaryConditionBase &bc) |
Calculates the sensitivity of force vector and Jacobian due to surface pressure applied on the entire element domain. More... | |
virtual bool | surface_pressure_residual_sensitivity (const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, const unsigned int side, MAST::BoundaryConditionBase &bc)=0 |
Calculates the force vector and Jacobian due to surface pressure. More... | |
virtual bool | surface_pressure_residual_sensitivity (const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc) |
Calculates the force vector and Jacobian due to surface pressure. More... | |
virtual bool | thermal_residual (bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)=0 |
Calculates the force vector and Jacobian due to thermal stresses. More... | |
virtual bool | thermal_residual_sensitivity (const MAST::FunctionBase &p, bool request_jacobian, RealVectorX &f, RealMatrixX &jac, MAST::BoundaryConditionBase &bc)=0 |
Calculates the sensitivity of force vector and Jacobian due to thermal stresses. More... | |
virtual void | thermal_residual_boundary_velocity (const MAST::FunctionBase &p, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, MAST::BoundaryConditionBase &bc, bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0 |
Calculates the sensitivity of force vector and Jacobian due to thermal stresses. More... | |
Protected Attributes | |
const MAST::ElementPropertyCardBase & | _property |
element property More... | |
RealVectorX | _local_sol |
local solution More... | |
RealVectorX | _local_delta_sol |
local perturbed solution More... | |
RealVectorX | _local_sol_sens |
local solution sensitivity More... | |
RealVectorX | _local_delta_sol_sens |
local perturbed solution sensitivity More... | |
RealVectorX | _local_vel |
local velocity More... | |
RealVectorX | _local_delta_vel |
local perturbed velocity More... | |
RealVectorX | _local_vel_sens |
local velocity sensitivity More... | |
RealVectorX | _local_delta_vel_sens |
local perturbed velocity sensitivity More... | |
RealVectorX | _local_accel |
local acceleration More... | |
RealVectorX | _local_delta_accel |
local perturbed acceleration More... | |
RealVectorX | _local_accel_sens |
local acceleration sensitivity More... | |
RealVectorX | _local_delta_accel_sens |
local perturbed acceleration sensitivity More... | |
RealVectorX * | _incompatible_sol |
incompatible mode solution vector More... | |
Protected Attributes inherited from MAST::ElementBase | |
MAST::SystemInitialization & | _system |
SystemInitialization object associated with this element. More... | |
MAST::AssemblyBase & | _assembly |
Assembly object. More... | |
const MAST::GeomElem & | _elem |
geometric element for which the computations are performed More... | |
MAST::FunctionBase * | _active_sol_function |
pointer to the active solution mesh field function. More... | |
const Real & | _time |
time for which system is being assembled More... | |
RealVectorX | _sol |
local solution More... | |
RealVectorX | _sol_sens |
local solution sensitivity 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_sol |
local solution used for linearized analysis More... | |
RealVectorX | _delta_sol_sens |
local solution used for linearized analysis More... | |
RealVectorX | _vel |
local velocity More... | |
RealVectorX | _vel_sens |
local velocity More... | |
RealVectorX | _delta_vel |
local velocity More... | |
RealVectorX | _delta_vel_sens |
local velocity More... | |
RealVectorX | _accel |
local acceleration More... | |
RealVectorX | _accel_sens |
local acceleration More... | |
RealVectorX | _delta_accel |
local acceleration More... | |
RealVectorX | _delta_accel_sens |
local acceleration More... | |
Definition at line 42 of file structural_element_base.h.
MAST::StructuralElementBase::StructuralElementBase | ( | MAST::SystemInitialization & | sys, |
MAST::AssemblyBase & | assembly, | ||
const MAST::GeomElem & | elem, | ||
const MAST::ElementPropertyCardBase & | p | ||
) |
Constructor.
Definition at line 38 of file structural_element_base.cpp.
|
virtual |
Definition at line 51 of file structural_element_base.cpp.
|
pure virtual |
Calculates the stress tensor.
If derivative and sensitivity with respect to the parameter sesitivity_param
are calculated and provided if the respective flags are true.
Implemented in MAST::StructuralElement3D, MAST::StructuralElement1D, and MAST::StructuralElement2D.
|
pure virtual |
Calculates the boundary velocity term contributions to the sensitivity of stress at the specified boundary of this element.
Implemented in MAST::StructuralElement3D, MAST::StructuralElement1D, and MAST::StructuralElement2D.
|
pure virtual |
Implemented in MAST::StructuralElement1D, MAST::StructuralElement2D, and MAST::StructuralElement3D.
|
inlinevirtual |
damping force contribution to system residual
Definition at line 165 of file structural_element_base.h.
|
inlinevirtual |
sensitivity of the damping force contribution to system residual
Definition at line 309 of file structural_element_base.h.
|
inline |
returns a constant reference to the finite element object
Definition at line 120 of file structural_element_base.h.
|
pure virtual |
Implemented in MAST::StructuralElement1D, MAST::StructuralElement2D, and MAST::StructuralElement3D.
|
inlinevirtual |
Reimplemented in MAST::StructuralElement3D.
Definition at line 380 of file structural_element_base.h.
|
virtual |
inertial force contribution to system residual
Reimplemented in MAST::StructuralElement3D.
Definition at line 258 of file structural_element_base.cpp.
|
virtual |
sensitivity of the inertial force contribution to system residual
Definition at line 509 of file structural_element_base.cpp.
|
virtual |
sensitivity of the inertial force contribution to system residual
Definition at line 404 of file structural_element_base.cpp.
|
pure virtual |
internal force contribution to system residual
Implemented in MAST::StructuralElement2D, MAST::StructuralElement1D, and MAST::StructuralElement3D.
|
pure virtual |
calculates the term on side s:
.
Implemented in MAST::StructuralElement2D, MAST::StructuralElement1D, and MAST::StructuralElement3D.
|
pure virtual |
calculates d[J]/d{x} .
d{x}/dp
Implemented in MAST::StructuralElement1D, MAST::StructuralElement3D, and MAST::StructuralElement2D.
|
pure virtual |
sensitivity of the internal force contribution to system residual
Implemented in MAST::StructuralElement3D, MAST::StructuralElement2D, and MAST::StructuralElement1D.
bool MAST::StructuralElementBase::linearized_frequency_domain_side_external_residual | ( | bool | request_jacobian, |
ComplexVectorX & | f, | ||
ComplexMatrixX & | jac, | ||
std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > & | bc | ||
) |
Calculates the external force due to frequency domain side external force contribution to system residual.
This primarily handles the boundary conditions. If requested, the Jacobian due to x is returned in jac.
Definition at line 735 of file structural_element_base.cpp.
|
protectedpure virtual |
Calculates the force vector and Jacobian due to small perturbation surface pressure.
Implemented in MAST::StructuralElement1D, MAST::StructuralElement3D, and MAST::StructuralElement2D.
|
protectedvirtual |
Calculates the force vector and Jacobian due to surface pressure applied on the entire element domain.
This is applicable for only 1D and 2D elements. The implementation can be used as the evaluation of , or the contribution of the off-diagonal Jacobian times fluid solution perturbation.
Definition at line 1452 of file structural_element_base.cpp.
|
protectedpure virtual |
Calculates the sensitivity of force vector and Jacobian due to small is applicable for perturbation surface pressure.
Implemented in MAST::StructuralElement1D, MAST::StructuralElement3D, and MAST::StructuralElement2D.
|
inlineprotectedvirtual |
Calculates the sensitivity of force vector and Jacobian due to surface pressure applied on the entire element domain.
This is applicable for only 1D and 2D elements.
Definition at line 633 of file structural_element_base.h.
bool MAST::StructuralElementBase::linearized_frequency_domain_volume_external_residual | ( | bool | request_jacobian, |
ComplexVectorX & | f, | ||
ComplexMatrixX & | jac, | ||
std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > & | bc | ||
) |
Calculates the frequency domain volume external force contribution to system residual.
If requested, the Jacobian due to x is returned in jac.
Definition at line 905 of file structural_element_base.cpp.
|
virtual |
inertial force contribution to system residual of linerized problem
Definition at line 361 of file structural_element_base.cpp.
|
virtual |
internal force contribution to system residual of the linearized problem
Definition at line 225 of file structural_element_base.cpp.
bool MAST::StructuralElementBase::linearized_side_external_residual | ( | bool | request_jacobian, |
RealVectorX & | f, | ||
RealMatrixX & | jac_xdot, | ||
RealMatrixX & | jac, | ||
std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > & | bc | ||
) |
side external force contribution to system residual.
This primarily handles the boundary conditions. If requested, the Jacobians of the residual due to xdot will be returned in jac_xdot
and the Jacobian due to x is returned in jac.
Definition at line 690 of file structural_element_base.cpp.
|
protectedvirtual |
Calculates the force vector and Jacobian due to surface pressure applied on the entire element domain.
This is applicable for only 1D and 2D elements.
Definition at line 1374 of file structural_element_base.cpp.
bool MAST::StructuralElementBase::linearized_volume_external_residual | ( | bool | request_jacobian, |
RealVectorX & | f, | ||
RealMatrixX & | jac_xdot, | ||
RealMatrixX & | jac, | ||
std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > & | bc | ||
) |
volume external force contribution to system residual.
If requested, the Jacobians of the residual due to xdot will be returned in jac_xdot
and the Jacobian due to x is returned in jac.
Definition at line 844 of file structural_element_base.cpp.
|
inline |
if_sens
is true) in the local element coordinate system Definition at line 109 of file structural_element_base.h.
|
protected |
Definition at line 1675 of file structural_element_base.cpp.
|
protected |
Definition at line 1706 of file structural_element_base.cpp.
|
protectedpure virtual |
Calculates the force vector and Jacobian due to piston-theory based surface pressure on the entire element domain.
This is applicable for only 1D and 2D elements. The order of the boundary condition and direction of fluid flow are obtained from the BoundaryConditionBase object.
Implemented in MAST::StructuralElement1D, MAST::StructuralElement2D, and MAST::StructuralElement3D.
|
protectedpure virtual |
Calculates the force vector and Jacobian due to piston-theory based surface pressure on the side identified by side
.
The order of the boundary condition and direction of fluid flow are obtained from the BoundaryConditionBase object.
Implemented in MAST::StructuralElement1D, MAST::StructuralElement3D, and MAST::StructuralElement2D.
|
protectedpure virtual |
Calculates the force vector and Jacobian due to piston-theory based surface pressure on the entire element domain.
This is applicable for only 1D and 2D elements. The order of the boundary condition and direction of fluid flow are obtained from the BoundaryConditionBase object.
Implemented in MAST::StructuralElement1D, MAST::StructuralElement2D, and MAST::StructuralElement3D.
|
protectedpure virtual |
Calculates the force vector and Jacobian due to piston-theory based surface pressure on the side identified by side
.
The order of the boundary condition and direction of fluid flow are obtained from the BoundaryConditionBase object.
Implemented in MAST::StructuralElement1D, MAST::StructuralElement3D, and MAST::StructuralElement2D.
|
pure virtual |
prestress force contribution to system residual
Implemented in MAST::StructuralElement1D, MAST::StructuralElement2D, and MAST::StructuralElement3D.
|
pure virtual |
sensitivity of the prestress force contribution to system residual
Implemented in MAST::StructuralElement1D, MAST::StructuralElement2D, and MAST::StructuralElement3D.
|
virtual |
stores vec
as acceleration for element level calculations, or its sensitivity if if_sens
is true.
Reimplemented from MAST::ElementBase.
Definition at line 171 of file structural_element_base.cpp.
|
inline |
sets the pointer to the incompatible mode solution vector.
This is stored as a pointer, and any modifications will overwrite the original vector
Definition at line 393 of file structural_element_base.h.
|
virtual |
stores vec
as perturbed acceleration for element level calculations, or its sensitivity if if_sens
is true.
Reimplemented from MAST::ElementBase.
Definition at line 198 of file structural_element_base.cpp.
|
virtual |
stores vec
as perturbed solution for element level calculations, or its sensitivity if if_sens
is true.
Reimplemented from MAST::ElementBase.
Definition at line 87 of file structural_element_base.cpp.
|
virtual |
stores vec
as perturbed velocity for element level calculations, or its sensitivity if if_sens
is true.
Reimplemented from MAST::ElementBase.
Definition at line 144 of file structural_element_base.cpp.
|
virtual |
stores vec
as solution for element level calculations, or its sensitivity if if_sens
is true.
Reimplemented from MAST::ElementBase.
Definition at line 58 of file structural_element_base.cpp.
|
virtual |
stores vec
as velocity for element level calculations, or its sensitivity if if_sens
is true.
Reimplemented from MAST::ElementBase.
Definition at line 117 of file structural_element_base.cpp.
bool MAST::StructuralElementBase::side_external_residual | ( | bool | request_jacobian, |
RealVectorX & | f, | ||
RealMatrixX & | jac_xdot, | ||
RealMatrixX & | jac, | ||
std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > & | bc | ||
) |
side external force contribution to system residual.
This primarily handles the boundary conditions. If requested, the Jacobians of the residual due to xdot will be returned in jac_xdot
and the Jacobian due to x is returned in jac.
Definition at line 631 of file structural_element_base.cpp.
bool MAST::StructuralElementBase::side_external_residual_sensitivity | ( | const MAST::FunctionBase & | p, |
bool | request_jacobian, | ||
RealVectorX & | f, | ||
RealMatrixX & | jac_xdot, | ||
RealMatrixX & | jac, | ||
std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase * > & | bc | ||
) |
sensitivity of the side external force contribution to system residual
Definition at line 954 of file structural_element_base.cpp.
|
protectedvirtual |
Calculates the force vector and Jacobian due to surface pressure applied on the entire element domain.
This is applicable for only 1D and 2D elements.
Definition at line 1281 of file structural_element_base.cpp.
|
protectedpure virtual |
Calculates the force vector and Jacobian due to surface pressure.
Implemented in MAST::StructuralElement1D, MAST::StructuralElement2D, and MAST::StructuralElement3D.
|
protectedvirtual |
Calculates the force vector and Jacobian due to surface pressure applied on the entire element domain.
This is applicable for only 1D and 2D elements.
Definition at line 1132 of file structural_element_base.cpp.
|
protectedpure virtual |
Calculates the force vector and Jacobian due to surface pressure.
Implemented in MAST::StructuralElement1D, MAST::StructuralElement2D, and MAST::StructuralElement3D.
|
protectedvirtual |
Calculates the force vector and Jacobian due to surface pressure.
this should be implemented for each element type
Definition at line 1207 of file structural_element_base.cpp.
|
protectedpure virtual |
Calculates the force vector and Jacobian due to thermal stresses.
this should be implemented for each element type
Implemented in MAST::StructuralElement1D, MAST::StructuralElement2D, and MAST::StructuralElement3D.
|
protectedpure virtual |
Calculates the sensitivity of force vector and Jacobian due to thermal stresses.
this should be implemented for each element type
Implemented in MAST::StructuralElement1D, MAST::StructuralElement2D, and MAST::StructuralElement3D.
|
protectedpure virtual |
Calculates the sensitivity of force vector and Jacobian due to thermal stresses.
this should be implemented for each element type
Implemented in MAST::StructuralElement1D, MAST::StructuralElement2D, and MAST::StructuralElement3D.
|
pure virtual |
Implemented in MAST::StructuralElement1D, MAST::StructuralElement2D, and MAST::StructuralElement3D.
template void MAST::StructuralElementBase::transform_matrix_to_global_system< ComplexMatrixX > | ( | const ValType & | local_mat, |
ValType & | global_mat | ||
) | const |
Definition at line 1531 of file structural_element_base.cpp.
template void MAST::StructuralElementBase::transform_vector_to_global_system< ComplexVectorX > | ( | const ValType & | local_vec, |
ValType & | global_vec | ||
) | const |
Definition at line 1607 of file structural_element_base.cpp.
template void MAST::StructuralElementBase::transform_vector_to_local_system< ComplexVectorX > | ( | const ValType & | global_vec, |
ValType & | local_vec | ||
) | const |
Definition at line 1571 of file structural_element_base.cpp.
|
inlinevirtual |
updates the incompatible solution for this element.
dsol
is the update to the element solution for the current nonlinear step.
Reimplemented in MAST::StructuralElement3D.
Definition at line 403 of file structural_element_base.h.
bool MAST::StructuralElementBase::volume_external_residual | ( | bool | request_jacobian, |
RealVectorX & | f, | ||
RealMatrixX & | jac_xdot, | ||
RealMatrixX & | jac, | ||
std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > & | bc | ||
) |
volume external force contribution to system residual.
If requested, the Jacobians of the residual due to xdot will be returned in jac_xdot
and the Jacobian due to x is returned in jac.
Definition at line 787 of file structural_element_base.cpp.
void MAST::StructuralElementBase::volume_external_residual_boundary_velocity | ( | const MAST::FunctionBase & | p, |
const unsigned int | s, | ||
const MAST::FieldFunction< RealVectorX > & | vel_f, | ||
std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > & | bc, | ||
bool | request_jacobian, | ||
RealVectorX & | f, | ||
RealMatrixX & | jac | ||
) |
boundary velocity contribution of volume external force.
Definition at line 1076 of file structural_element_base.cpp.
bool MAST::StructuralElementBase::volume_external_residual_sensitivity | ( | const MAST::FunctionBase & | p, |
bool | request_jacobian, | ||
RealVectorX & | f, | ||
RealMatrixX & | jac_xdot, | ||
RealMatrixX & | jac, | ||
std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > & | bc | ||
) |
sensitivity of the volume external force contribution to system residual
Definition at line 1017 of file structural_element_base.cpp.
|
protected |
incompatible mode solution vector
Definition at line 779 of file structural_element_base.h.
|
protected |
local acceleration
Definition at line 755 of file structural_element_base.h.
|
protected |
local acceleration sensitivity
Definition at line 767 of file structural_element_base.h.
|
protected |
local perturbed acceleration
Definition at line 761 of file structural_element_base.h.
|
protected |
local perturbed acceleration sensitivity
Definition at line 773 of file structural_element_base.h.
|
protected |
local perturbed solution
Definition at line 713 of file structural_element_base.h.
|
protected |
local perturbed solution sensitivity
Definition at line 725 of file structural_element_base.h.
|
protected |
local perturbed velocity
Definition at line 737 of file structural_element_base.h.
|
protected |
local perturbed velocity sensitivity
Definition at line 749 of file structural_element_base.h.
|
protected |
local solution
Definition at line 707 of file structural_element_base.h.
|
protected |
local solution sensitivity
Definition at line 719 of file structural_element_base.h.
|
protected |
local velocity
Definition at line 731 of file structural_element_base.h.
|
protected |
local velocity sensitivity
Definition at line 743 of file structural_element_base.h.
|
protected |
element property
Definition at line 701 of file structural_element_base.h.
bool MAST::StructuralElementBase::follower_forces |
flag for follower forces
Definition at line 440 of file structural_element_base.h.