https://mooseframework.inl.gov
ShaftTimeDerivativeScalarKernel.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 
12 #include "UserObject.h"
13 #include "MooseVariableScalar.h"
14 #include "Assembly.h"
15 
17 
20 {
22  params.addRequiredParam<std::vector<UserObjectName>>("uo_names",
23  "Names of shaft-connectable user objects");
24  params.addClassDescription("Adds a time derivative term to the shaft ODE");
25  params.set<MultiMooseEnum>("vector_tags") = "time";
26  params.set<MultiMooseEnum>("matrix_tags") = "system time";
27 
28  return params;
29 }
30 
32  : ScalarKernel(parameters),
33  _u_dot(_var.uDot()),
34  _du_dot_du(_var.duDotDu()),
35  _uo_names(getParam<std::vector<UserObjectName>>("uo_names")),
36  _n_components(_uo_names.size())
37 {
39  for (unsigned int i = 0; i < _n_components; ++i)
40  {
42  &getUserObjectByName<ShaftConnectableUserObjectInterface>(_uo_names[i]);
43  }
44 }
45 
46 void
48 {
49 }
50 
51 void
53 {
55 
56  Real sum_inertias = 0;
57  for (unsigned int i = 0; i < _n_components; ++i)
58  sum_inertias += _shaft_connected_uos[i]->getMomentOfInertia();
59 
60  _local_re(0) += sum_inertias * _u_dot[0];
61 
63 }
64 
65 void
67 {
69 
70  Real sum_inertias = 0;
71  for (unsigned int i = 0; i < _n_components; ++i)
72  sum_inertias += _shaft_connected_uos[i]->getMomentOfInertia();
73 
74  // FIXME: Add dI_domega * domega_dt
75  _local_ke(0, 0) += sum_inertias * _du_dot_du[0];
76 
78 
79  for (unsigned int i = 0; i < _n_components; ++i)
80  {
81  DenseMatrix<Real> jacobian_block;
82  std::vector<dof_id_type> dofs_j;
83  _shaft_connected_uos[i]->getMomentOfInertiaJacobianData(jacobian_block, dofs_j);
84  jacobian_block.scale(_u_dot[0]);
85  addJacobian(_assembly, jacobian_block, _var.dofIndices(), dofs_j, _var.scalingFactor());
86  }
87 }
void accumulateTaggedLocalResidual()
unsigned int number() const
std::vector< const ShaftConnectableUserObjectInterface * > _shaft_connected_uos
List of shaft connected user objects.
T & set(const std::string &name, bool quiet_mode=false)
const VariableValue & _du_dot_du
Derivative of u_dot wrt u.
MooseVariableScalar & _var
DenseMatrix< Number > _local_ke
void addRequiredParam(const std::string &name, const std::string &doc_string)
const VariableValue & _u_dot
Time derivative of u.
Time derivative for angular speed of shaft.
void addJacobian(Assembly &assembly, const Residuals &residuals, const Indices &dof_indices, Real scaling_factor)
ShaftTimeDerivativeScalarKernel(const InputParameters &parameters)
virtual const std::vector< dof_id_type > & dofIndices() const
void accumulateTaggedLocalMatrix()
unsigned int _n_components
Number of shaft connected user objects.
Assembly & _assembly
void scale(const Real factor)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
DenseVector< Number > _local_re
void addClassDescription(const std::string &doc_string)
registerMooseObject("ThermalHydraulicsApp", ShaftTimeDerivativeScalarKernel)
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
void prepareMatrixTag(Assembly &assembly, unsigned int ivar, unsigned int jvar)
static InputParameters validParams()
void scalingFactor(const std::vector< Real > &factor)
const std::vector< UserObjectName > & _uo_names
List of names of shaft connected user objects.