20 #ifndef __mast__structural_element_base__ 21 #define __mast__structural_element_base__ 35 class ElementPropertyCardBase;
36 class BoundaryConditionBase;
37 class FEMOperatorMatrix;
38 class StressStrainOutputBase;
39 template <
typename ValType>
class FieldFunction;
62 bool if_sens =
false);
70 bool if_sens =
false);
78 bool if_sens =
false);
85 bool if_sens =
false);
93 bool if_sens =
false);
101 bool if_sens =
false);
148 const unsigned int s,
150 bool request_jacobian,
204 std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
214 (
bool request_jacobian,
218 std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
230 (
bool request_jacobian,
233 std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
253 std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
261 const unsigned int s,
263 std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc,
264 bool request_jacobian,
277 (
bool request_jacobian,
281 std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
291 (
bool request_jacobian,
294 std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
302 bool request_jacobian,
310 bool request_jacobian,
321 bool request_jacobian,
332 const unsigned int s,
334 bool request_jacobian,
344 bool request_jacobian,
348 std::multimap<libMesh::boundary_id_type, MAST::BoundaryConditionBase*>& bc);
354 bool request_jacobian,
362 bool request_jacobian,
366 std::multimap<libMesh::subdomain_id_type, MAST::BoundaryConditionBase*>& bc);
426 const unsigned int s,
443 template <
typename ValType>
445 ValType& global_mat)
const;
447 template <
typename ValType>
449 ValType& local_vec)
const;
451 template <
typename ValType>
453 ValType& global_vec)
const;
464 const unsigned int side,
486 const unsigned int s,
489 bool request_jacobian,
549 const unsigned int side,
562 bool request_jacobian,
575 bool request_jacobian,
579 const unsigned int side,
589 (
bool request_jacobian,
592 const unsigned int side,
606 (
bool request_jacobian,
619 bool request_jacobian,
622 const unsigned int side,
634 bool request_jacobian,
648 bool request_jacobian,
651 const unsigned int side,
660 bool request_jacobian,
681 bool request_jacobian,
691 const unsigned int s,
694 bool request_jacobian,
787 std::unique_ptr<MAST::StructuralElementBase>
796 #endif // __mast__structural_element_base__
RealVectorX _local_accel_sens
local acceleration sensitivity
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: .
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...
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.
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.
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.
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.
virtual void update_incompatible_mode_solution(const RealVectorX &dsol)
updates the incompatible solution for this element.
Data structure provides the mechanism to store stress and strain output from a structural analysis...
RealVectorX _local_delta_vel_sens
local perturbed velocity sensitivity
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
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 el...
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
std::unique_ptr< MAST::StructuralElementBase > build_structural_element(MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
builds the structural element for the specified element type
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.
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...
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...
const RealVectorX & local_solution(bool if_sens=false) const
const MAST::GeomElem & elem() const
Matrix< Complex, Dynamic, 1 > ComplexVectorX
void transform_vector_to_local_system(const ValType &global_vec, ValType &local_vec) const
RealVectorX _local_delta_sol_sens
local perturbed solution sensitivity
virtual bool inertial_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xddot, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual
RealVectorX _local_delta_accel_sens
local perturbed acceleration sensitivity
virtual bool damping_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
damping force contribution to system residual
RealVectorX _local_accel
local acceleration
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.
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
RealVectorX _local_sol_sens
local solution sensitivity
MAST::AssemblyBase & assembly()
virtual bool if_incompatible_modes() const =0
virtual void thermal_residual_temperature_derivative(const MAST::FEBase &fe_thermal, RealMatrixX &m)=0
Real piston_theory_cp(const unsigned int order, const Real vel_U, const Real gamma, const Real mach)
RealVectorX _local_delta_vel
local perturbed velocity
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...
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
void transform_matrix_to_global_system(const ValType &local_mat, ValType &global_mat) const
RealVectorX _local_delta_sol
local perturbed solution
bool follower_forces
flag for follower forces
const MAST::ElementPropertyCardBase & elem_property()
returns a constant reference to the finite element object
Matrix< Real, Dynamic, Dynamic > RealMatrixX
Matrix< Complex, Dynamic, Dynamic > ComplexMatrixX
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
Real piston_theory_dcp_dv(const unsigned int order, const Real vel_U, const Real gamma, const Real mach)
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.
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...
virtual bool calculate_stress(bool request_derivative, const MAST::FunctionBase *f, MAST::StressStrainOutputBase &output)=0
Calculates the stress tensor.
Matrix< Real, Dynamic, 1 > RealVectorX
virtual ~StructuralElementBase()
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 s...
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 el...
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.
const MAST::ElementPropertyCardBase & _property
element property
virtual unsigned int incompatible_mode_size() const
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 bou...
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
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 tru...
RealVectorX _local_vel_sens
local velocity sensitivity
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...
virtual void calculate_stress_temperature_derivative(MAST::FEBase &fe_thermal, MAST::StressStrainOutputBase &output)=0
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 tru...
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
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 resi...
StructuralElementBase(MAST::SystemInitialization &sys, MAST::AssemblyBase &assembly, const MAST::GeomElem &elem, const MAST::ElementPropertyCardBase &p)
Constructor.
virtual bool prestress_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
prestress force contribution to system residual
virtual bool internal_residual_jac_dot_state_sensitivity(RealMatrixX &jac)=0
calculates d[J]/d{x} .
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
virtual bool linearized_internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual of the linearized problem
RealVectorX * _incompatible_sol
incompatible mode solution vector
RealVectorX _local_delta_accel
local perturbed acceleration
RealVectorX _local_sol
local solution
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.
RealVectorX _local_vel
local velocity
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.
void set_incompatible_mode_solution(RealVectorX &vec)
sets the pointer to the incompatible mode solution vector.
virtual bool internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)=0
internal force contribution to system residual
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
void transform_vector_to_global_system(const ValType &local_vec, ValType &global_vec) const
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.
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.
This is the base class for elements that implement calculation of finite element quantities over the ...