https://mooseframework.inl.gov
ShaftConnectableUserObjectInterface.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 #include "MooseVariableScalar.h"
12 #include "UserObject.h"
13 
16 {
18  return params;
19 }
20 
22  const MooseObject * /*moose_object*/)
23  : _n_shaft_eq(1)
24 {
25  _omega_dof.resize(_n_shaft_eq);
28 }
29 
30 void
32 {
33  _torque = 0;
35 }
36 
37 void
39 {
40 }
41 
42 Real
44 {
45  return _torque;
46 }
47 
48 void
50  std::vector<dof_id_type> & dofs_j) const
51 {
52  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
53 
54  unsigned int n_dofs = 1 + _scalar_dofs.size();
55  for (unsigned int c = 0; c < _n_connections; c++)
57 
58  jacobian_block.resize(_n_shaft_eq, n_dofs);
59  dofs_j.resize(n_dofs);
60 
61  unsigned int k = 0;
62  jacobian_block(0, k) = _torque_jacobian_omega_var(0, 0);
63  dofs_j[k] = _omega_dof[0];
64  k++;
65 
66  // Store Jacobian entries w.r.t. scalar variables
67  for (unsigned int j = 0; j < _scalar_dofs.size(); j++, k++)
68  {
69  jacobian_block(0, k) = _torque_jacobian_scalar_vars(0, j);
70  dofs_j[k] = _scalar_dofs[j];
71  }
72  // Store Jacobian entries w.r.t. flow variables
73  for (unsigned int c = 0; c < _n_connections; c++)
74  {
75  for (unsigned int j = 0; j < _torque_jacobian_flow_channel_vars[c].n(); j++, k++)
76  {
77  jacobian_block(0, k) = _torque_jacobian_flow_channel_vars[c](0, j);
78  dofs_j[k] = _flow_channel_dofs[c][j];
79  }
80  }
81 }
82 
83 Real
85 {
86  return _moment_of_inertia;
87 }
88 
89 void
91  DenseMatrix<Real> & jacobian_block, std::vector<dof_id_type> & dofs_j) const
92 {
93  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
94 
95  unsigned int n_dofs = 1 + _scalar_dofs.size();
96  for (unsigned int c = 0; c < _n_connections; c++)
97  n_dofs += _moi_jacobian_flow_channel_vars[c].n();
98 
99  jacobian_block.resize(_n_shaft_eq, n_dofs);
100  dofs_j.resize(n_dofs);
101 
102  unsigned int k = 0;
103  jacobian_block(0, k) = _moi_jacobian_omega_var(0, 0);
104  dofs_j[k] = _omega_dof[0];
105  k++;
106 
107  // Store Jacobian entries w.r.t. scalar variables
108  for (unsigned int j = 0; j < _scalar_dofs.size(); j++, k++)
109  {
110  jacobian_block(0, k) = _moi_jacobian_scalar_vars(0, j);
111  dofs_j[k] = _scalar_dofs[j];
112  }
113  // Store Jacobian entries w.r.t. flow variables
114  for (unsigned int c = 0; c < _n_connections; c++)
115  {
116  for (unsigned int j = 0; j < _moi_jacobian_flow_channel_vars[c].n(); j++, k++)
117  {
118  jacobian_block(0, k) = _moi_jacobian_flow_channel_vars[c](0, j);
119  dofs_j[k] = _flow_channel_dofs[c][j];
120  }
121  }
122 }
123 
124 void
126  unsigned int n_flow_eq)
127 {
128  _n_connections = n_connections;
131 
132  _n_flow_eq = n_flow_eq;
133 }
134 
135 void
137  const std::vector<std::vector<std::vector<Real>>> & phi_face_values,
138  const std::vector<std::vector<dof_id_type>> & flow_channel_dofs)
139 {
140  _flow_channel_dofs = flow_channel_dofs;
141  _phi_face_values = phi_face_values;
142 }
143 
144 void
146 {
147  auto && dofs = omega_var->dofIndices();
148  mooseAssert(dofs.size() == 1,
149  "There should be exactly 1 coupled DoF index for the variable '" + omega_var->name() +
150  "'.");
151  _omega_dof = dofs;
152 }
153 
154 void
155 ShaftConnectableUserObjectInterface::setupJunctionData(std::vector<dof_id_type> & scalar_dofs)
156 {
157  _scalar_dofs = scalar_dofs;
158 }
159 
160 void
162  const DenseMatrix<Real> & jac, const unsigned int & c)
163 {
165  unsigned int jk = 0;
166  for (unsigned int j = 0; j < _n_flow_eq; j++)
167  {
168  for (unsigned int k = 0; k < _phi_face_values[c][j].size(); k++)
169  {
170  _torque_jacobian_flow_channel_vars[c](0, jk) = jac(0, j) * _phi_face_values[c][j][k];
171  jk++;
172  }
173  }
174 }
175 
176 void
178  const DenseMatrix<Real> & jac, const unsigned int & c)
179 {
181  unsigned int jk = 0;
182  for (unsigned int j = 0; j < _n_flow_eq; j++)
183  {
184  for (unsigned int k = 0; k < _phi_face_values[c][j].size(); k++)
185  {
186  _moi_jacobian_flow_channel_vars[c](0, jk) = jac(0, j) * _phi_face_values[c][j][k];
187  jk++;
188  }
189  }
190 }
191 
192 void
194 {
195 }
196 
197 void
199 {
200  const ShaftConnectableUserObjectInterface & sctc_uo =
201  dynamic_cast<const ShaftConnectableUserObjectInterface &>(uo);
202  _torque += sctc_uo._torque;
204 }
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)
const std::string & name() const override
virtual void setupJunctionData(std::vector< dof_id_type > &scalar_dofs)
Stores data associated with a junction component.
InputParameters emptyInputParameters()
virtual void computeTorqueScalarJacobianWRTFlowDofs(const DenseMatrix< Real > &jac, const unsigned int &c)
virtual void setupConnections(unsigned int n_connections, unsigned int n_flow_eq)
ShaftConnectableUserObjectInterface(const MooseObject *moose_object)
virtual void setConnectionData(const std::vector< std::vector< std::vector< Real >>> &phi_face_values, const std::vector< std::vector< dof_id_type >> &flow_channel_dofs)
Stores data computed by a volume-junction-like object associated with the conection.
unsigned int _n_flow_eq
Number of flow variables in connected flow channels.
virtual const std::vector< dof_id_type > & dofIndices() const
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)
virtual void setOmegaDofs(const MooseVariableScalar *omega_var)
virtual void computeMomentOfInertiaScalarJacobianWRTFlowDofs(const DenseMatrix< Real > &jac, const unsigned int &c)
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.
void resize(const unsigned int new_m, const unsigned int new_n)
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
std::vector< DenseMatrix< Real > > _torque_jacobian_flow_channel_vars
Cached scalar residual Jacobian matrices w.r.t.
static const std::string k
Definition: NS.h:130
virtual void getTorqueJacobianData(DenseMatrix< Real > &jacobian_block, std::vector< dof_id_type > &dofs_j) const
virtual void getMomentOfInertiaJacobianData(DenseMatrix< Real > &jacobian_block, std::vector< dof_id_type > &dofs_j) const
std::vector< DenseMatrix< Real > > _moi_jacobian_flow_channel_vars
unsigned int _n_connections
Number of flow channels the shaft connected component is attached to.