12 #include "libmesh/dense_matrix.h" 28 virtual void execute();
29 virtual void finalize();
30 virtual void threadJoin(
const UserObject & uo);
32 virtual Real getTorque()
const;
34 std::vector<dof_id_type> & dofs_j)
const;
35 virtual Real getMomentOfInertia()
const;
37 std::vector<dof_id_type> & dofs_j)
const;
40 virtual void setupConnections(
unsigned int n_connections,
unsigned int n_flow_eq);
45 setConnectionData(
const std::vector<std::vector<std::vector<Real>>> & phi_face_values,
46 const std::vector<std::vector<dof_id_type>> & flow_channel_dofs);
51 virtual void setupJunctionData(std::vector<dof_id_type> & scalar_dofs);
53 virtual void computeMomentOfInertiaScalarJacobianWRTFlowDofs(
const DenseMatrix<Real> & jac,
54 const unsigned int &
c);
55 virtual void computeTorqueScalarJacobianWRTFlowDofs(
const DenseMatrix<Real> & jac,
56 const unsigned int &
c);
std::vector< std::vector< std::vector< Real > > > _phi_face_values
Side shape function value (i.e. side from the flow channels)
DenseMatrix< Real > _moi_jacobian_omega_var
Jacobian entries of moment of inertia wrt to omega variable (from shaft)
DenseMatrix< Real > _torque_jacobian_scalar_vars
Jacobian entries of torque wrt to scalar variables (from junction)
Real _torque
Total torque.
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
void initialize(EquationSystems &es, const std::string &system_name)
InputParameters validParams()
unsigned int _n_flow_eq
Number of flow variables in connected flow channels.
Interface class for user objects that are connected to a shaft.
unsigned int _n_shaft_eq
Number of equation in the shaft component.
std::vector< dof_id_type > _omega_dof
Degrees of freedom for omega variable (from shaft)
DenseMatrix< Real > _moi_jacobian_scalar_vars
Jacobian entries of moment of inertia wrt to omega scalar variables (from junction) ...
std::vector< dof_id_type > _scalar_dofs
Degrees of freedom for scalar variables (from junction)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
DenseMatrix< Real > _torque_jacobian_omega_var
Jacobian entries of torque wrt to omega variable (from shaft)
std::vector< std::vector< dof_id_type > > _flow_channel_dofs
Degrees of freedom for flow channel variables, for each connection.
std::vector< DenseMatrix< Real > > _torque_jacobian_flow_channel_vars
Cached scalar residual Jacobian matrices w.r.t.
std::vector< DenseMatrix< Real > > _moi_jacobian_flow_channel_vars
unsigned int _n_connections
Number of flow channels the shaft connected component is attached to.
Real _moment_of_inertia
Moment of inertia.