Line data Source code
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 : 10 : #include "ShaftConnectableUserObjectInterface.h" 11 : #include "MooseVariableScalar.h" 12 : #include "UserObject.h" 13 : 14 : InputParameters 15 10 : ShaftConnectableUserObjectInterface::validParams() 16 : { 17 10 : InputParameters params = emptyInputParameters(); 18 10 : return params; 19 : } 20 : 21 5 : ShaftConnectableUserObjectInterface::ShaftConnectableUserObjectInterface( 22 5 : const MooseObject * /*moose_object*/) 23 5 : : _n_shaft_eq(1) 24 : { 25 5 : _omega_dof.resize(_n_shaft_eq); 26 5 : _torque_jacobian_omega_var.resize(1, _n_shaft_eq); 27 5 : _moi_jacobian_omega_var.resize(1, _n_shaft_eq); 28 5 : } 29 : 30 : void 31 0 : ShaftConnectableUserObjectInterface::initialize() 32 : { 33 0 : _torque = 0; 34 0 : _moment_of_inertia = 0; 35 0 : } 36 : 37 : void 38 0 : ShaftConnectableUserObjectInterface::execute() 39 : { 40 0 : } 41 : 42 : Real 43 0 : ShaftConnectableUserObjectInterface::getTorque() const 44 : { 45 0 : return _torque; 46 : } 47 : 48 : void 49 0 : ShaftConnectableUserObjectInterface::getTorqueJacobianData(DenseMatrix<Real> & jacobian_block, 50 : std::vector<dof_id_type> & dofs_j) const 51 : { 52 : Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 53 : 54 0 : unsigned int n_dofs = 1 + _scalar_dofs.size(); 55 0 : for (unsigned int c = 0; c < _n_connections; c++) 56 0 : n_dofs += _torque_jacobian_flow_channel_vars[c].n(); 57 : 58 0 : jacobian_block.resize(_n_shaft_eq, n_dofs); 59 0 : dofs_j.resize(n_dofs); 60 : 61 : unsigned int k = 0; 62 0 : jacobian_block(0, k) = _torque_jacobian_omega_var(0, 0); 63 0 : dofs_j[k] = _omega_dof[0]; 64 : k++; 65 : 66 : // Store Jacobian entries w.r.t. scalar variables 67 0 : for (unsigned int j = 0; j < _scalar_dofs.size(); j++, k++) 68 : { 69 0 : jacobian_block(0, k) = _torque_jacobian_scalar_vars(0, j); 70 0 : dofs_j[k] = _scalar_dofs[j]; 71 : } 72 : // Store Jacobian entries w.r.t. flow variables 73 0 : for (unsigned int c = 0; c < _n_connections; c++) 74 : { 75 0 : for (unsigned int j = 0; j < _torque_jacobian_flow_channel_vars[c].n(); j++, k++) 76 : { 77 0 : jacobian_block(0, k) = _torque_jacobian_flow_channel_vars[c](0, j); 78 0 : dofs_j[k] = _flow_channel_dofs[c][j]; 79 : } 80 : } 81 0 : } 82 : 83 : Real 84 0 : ShaftConnectableUserObjectInterface::getMomentOfInertia() const 85 : { 86 0 : return _moment_of_inertia; 87 : } 88 : 89 : void 90 0 : ShaftConnectableUserObjectInterface::getMomentOfInertiaJacobianData( 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 0 : unsigned int n_dofs = 1 + _scalar_dofs.size(); 96 0 : for (unsigned int c = 0; c < _n_connections; c++) 97 0 : n_dofs += _moi_jacobian_flow_channel_vars[c].n(); 98 : 99 0 : jacobian_block.resize(_n_shaft_eq, n_dofs); 100 0 : dofs_j.resize(n_dofs); 101 : 102 : unsigned int k = 0; 103 0 : jacobian_block(0, k) = _moi_jacobian_omega_var(0, 0); 104 0 : dofs_j[k] = _omega_dof[0]; 105 : k++; 106 : 107 : // Store Jacobian entries w.r.t. scalar variables 108 0 : for (unsigned int j = 0; j < _scalar_dofs.size(); j++, k++) 109 : { 110 0 : jacobian_block(0, k) = _moi_jacobian_scalar_vars(0, j); 111 0 : dofs_j[k] = _scalar_dofs[j]; 112 : } 113 : // Store Jacobian entries w.r.t. flow variables 114 0 : for (unsigned int c = 0; c < _n_connections; c++) 115 : { 116 0 : for (unsigned int j = 0; j < _moi_jacobian_flow_channel_vars[c].n(); j++, k++) 117 : { 118 0 : jacobian_block(0, k) = _moi_jacobian_flow_channel_vars[c](0, j); 119 0 : dofs_j[k] = _flow_channel_dofs[c][j]; 120 : } 121 : } 122 0 : } 123 : 124 : void 125 0 : ShaftConnectableUserObjectInterface::setupConnections(unsigned int n_connections, 126 : unsigned int n_flow_eq) 127 : { 128 0 : _n_connections = n_connections; 129 0 : _torque_jacobian_flow_channel_vars.resize(_n_connections); 130 0 : _moi_jacobian_flow_channel_vars.resize(_n_connections); 131 : 132 0 : _n_flow_eq = n_flow_eq; 133 0 : } 134 : 135 : void 136 0 : ShaftConnectableUserObjectInterface::setConnectionData( 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 0 : _flow_channel_dofs = flow_channel_dofs; 141 0 : _phi_face_values = phi_face_values; 142 0 : } 143 : 144 : void 145 0 : ShaftConnectableUserObjectInterface::setOmegaDofs(const MooseVariableScalar * omega_var) 146 : { 147 0 : 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 0 : _omega_dof = dofs; 152 0 : } 153 : 154 : void 155 0 : ShaftConnectableUserObjectInterface::setupJunctionData(std::vector<dof_id_type> & scalar_dofs) 156 : { 157 0 : _scalar_dofs = scalar_dofs; 158 0 : } 159 : 160 : void 161 0 : ShaftConnectableUserObjectInterface::computeTorqueScalarJacobianWRTFlowDofs( 162 : const DenseMatrix<Real> & jac, const unsigned int & c) 163 : { 164 0 : _torque_jacobian_flow_channel_vars[c].resize(1, _flow_channel_dofs[c].size()); 165 : unsigned int jk = 0; 166 0 : for (unsigned int j = 0; j < _n_flow_eq; j++) 167 : { 168 0 : for (unsigned int k = 0; k < _phi_face_values[c][j].size(); k++) 169 : { 170 0 : _torque_jacobian_flow_channel_vars[c](0, jk) = jac(0, j) * _phi_face_values[c][j][k]; 171 0 : jk++; 172 : } 173 : } 174 0 : } 175 : 176 : void 177 0 : ShaftConnectableUserObjectInterface::computeMomentOfInertiaScalarJacobianWRTFlowDofs( 178 : const DenseMatrix<Real> & jac, const unsigned int & c) 179 : { 180 0 : _moi_jacobian_flow_channel_vars[c].resize(1, _flow_channel_dofs[c].size()); 181 : unsigned int jk = 0; 182 0 : for (unsigned int j = 0; j < _n_flow_eq; j++) 183 : { 184 0 : for (unsigned int k = 0; k < _phi_face_values[c][j].size(); k++) 185 : { 186 0 : _moi_jacobian_flow_channel_vars[c](0, jk) = jac(0, j) * _phi_face_values[c][j][k]; 187 0 : jk++; 188 : } 189 : } 190 0 : } 191 : 192 : void 193 0 : ShaftConnectableUserObjectInterface::finalize() 194 : { 195 0 : } 196 : 197 : void 198 0 : ShaftConnectableUserObjectInterface::threadJoin(const UserObject & uo) 199 : { 200 : const ShaftConnectableUserObjectInterface & sctc_uo = 201 0 : dynamic_cast<const ShaftConnectableUserObjectInterface &>(uo); 202 0 : _torque += sctc_uo._torque; 203 0 : _moment_of_inertia += sctc_uo._moment_of_inertia; 204 0 : }