www.mooseframework.org
EigenKernel.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
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 "EigenKernel.h"
11 
12 // MOOSE includes
13 #include "Assembly.h"
14 #include "EigenExecutionerBase.h"
15 #include "Executioner.h"
16 #include "MooseApp.h"
17 #include "MooseEigenSystem.h"
18 #include "MooseVariableFE.h"
19 
20 #include "libmesh/quadrature.h"
21 
24 {
26  params.addParam<bool>(
27  "eigen", true, "Use for eigenvalue problem (true) or source problem (false)");
28  params.addParam<PostprocessorName>(
29  "eigen_postprocessor", 1.0, "The name of the postprocessor that provides the eigenvalue.");
30  params.registerBase("EigenKernel");
31  return params;
32 }
33 
35  : Kernel(parameters),
36  _eigen(getParam<bool>("eigen")),
37  _eigen_sys(
38  dynamic_cast<MooseEigenSystem *>(&_fe_problem.getNonlinearSystemBase(_sys.number()))),
39  _eigenvalue(NULL)
40 {
41  // The name to the postprocessor storing the eigenvalue
42  std::string eigen_pp_name;
43 
44  // If the "eigen_postprocessor" is given, use it. The isParamValid does not work here because of
45  // the default value, which
46  // you don't want to use if an EigenExecutioner exists.
47  if (!isDefaultPostprocessorValue("eigen_postprocessor"))
48  eigen_pp_name = getPostprocessorName("eigen_postprocessor");
49 
50  // Attempt to extract the eigenvalue postprocessor from the Executioner
51  else
52  {
54  if (exec)
55  eigen_pp_name = exec->getParam<PostprocessorName>("bx_norm");
56  }
57 
58  // If the postprocessor name was not provided and an EigenExecutionerBase is not being used,
59  // use the default value from the "eigen_postprocessor" parameter
60  if (eigen_pp_name.empty())
61  _eigenvalue = &getPostprocessorValue("eigen_postprocessor");
62 
63  // If the name does exist, then use the postprocessor value
64  else
65  {
66  if (_is_implicit)
67  _eigenvalue = &getPostprocessorValueByName(eigen_pp_name);
68  else
69  {
71  if (exec)
72  _eigenvalue = &exec->eigenvalueOld();
73  else
75  }
76  }
77 }
78 
79 void
81 {
83 
85 
86  mooseAssert(*_eigenvalue != 0.0, "Can't divide by zero eigenvalue in EigenKernel!");
87  Real one_over_eigen = 1.0 / *_eigenvalue;
88  for (_i = 0; _i < _test.size(); _i++)
89  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
90  _local_re(_i) += _JxW[_qp] * _coord[_qp] * one_over_eigen * computeQpResidual();
91 
93  if (_has_save_in)
94  {
95  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
96  for (const auto & var : _save_in)
97  var->sys().solution().add_vector(_local_re, var->dofIndices());
98  }
99 }
100 
101 void
103 {
104  if (!_is_implicit)
105  return;
106 
108 
110  mooseAssert(*_eigenvalue != 0.0, "Can't divide by zero eigenvalue in EigenKernel!");
111  Real one_over_eigen = 1.0 / *_eigenvalue;
112  for (_i = 0; _i < _test.size(); _i++)
113  for (_j = 0; _j < _phi.size(); _j++)
114  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
115  _local_ke(_i, _j) += _JxW[_qp] * _coord[_qp] * one_over_eigen * computeQpJacobian();
116 
118 
120  {
122  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
123  for (const auto & var : _diag_save_in)
124  var->sys().solution().add_vector(diag, var->dofIndices());
125  }
126 }
127 
128 void
129 EigenKernel::computeOffDiagJacobian(const unsigned int jvar_num)
130 {
131  if (!_is_implicit)
132  return;
133 
134  if (jvar_num == _var.number())
135  computeJacobian();
136  else
137  {
138  const auto & jvar = getVariable(jvar_num);
139 
140  prepareMatrixTag(_assembly, _var.number(), jvar_num);
141 
142  // This (undisplaced) jvar could potentially yield the wrong phi size if this object is acting
143  // on the displaced mesh
144  auto phi_size = jvar.dofIndices().size();
145  mooseAssert(
146  phi_size * jvar.count() == _local_ke.n(),
147  "The size of the phi container does not match the number of local Jacobian columns");
148 
149  if (_local_ke.m() != _test.size())
150  return;
151 
152  precalculateOffDiagJacobian(jvar_num);
153  mooseAssert(*_eigenvalue != 0.0, "Can't divide by zero eigenvalue in EigenKernel!");
154  Real one_over_eigen = 1.0 / *_eigenvalue;
155  if (jvar.count() == 1)
156  {
157  for (_i = 0; _i < _test.size(); _i++)
158  for (_j = 0; _j < phi_size; _j++)
159  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
160  _local_ke(_i, _j) +=
161  _JxW[_qp] * _coord[_qp] * one_over_eigen * computeQpOffDiagJacobian(jvar.number());
162  }
163  else
164  {
165  unsigned int n = phi_size;
166  for (_i = 0; _i < _test.size(); _i++)
167  for (_j = 0; _j < n; _j++)
168  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
169  {
170  RealEigenVector v = _JxW[_qp] * _coord[_qp] * one_over_eigen *
171  computeQpOffDiagJacobianArray(static_cast<ArrayMooseVariable &>(
172  const_cast<MooseVariableFieldBase &>(jvar)));
173  for (unsigned int k = 0; k < v.size(); ++k)
174  _local_ke(_i, _j + k * n) += v(k);
175  }
176  }
177 
179  }
180 }
181 
182 bool
184 {
185  bool flag = MooseObject::enabled();
186  if (_eigen)
187  {
188  if (_is_implicit)
189  return flag && (!_eigen_sys->activeOnOld());
190  else
191  return flag && _eigen_sys->activeOnOld();
192  }
193  else
194  return flag;
195 }
MooseEigenSystem * _eigen_sys
EigenKernel always lives in EigenSystem.
Definition: EigenKernel.h:41
static InputParameters validParams()
Definition: Kernel.C:23
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
virtual Real computeQpResidual()=0
Compute this Kernel&#39;s contribution to the residual at the current quadrature point.
virtual RealEigenVector computeQpOffDiagJacobianArray(const ArrayMooseVariable &jvar)
For coupling array variables.
Definition: Kernel.h:66
MooseVariable & _var
This is a regular kernel so we cast to a regular MooseVariable.
Definition: Kernel.h:72
std::vector< MooseVariableFEBase * > _save_in
Definition: KernelBase.h:64
const MooseArray< Real > & _JxW
The current quadrature point weight value.
Definition: KernelBase.h:51
virtual void precalculateOffDiagJacobian(unsigned int)
const Real * _eigenvalue
A pointer to the eigenvalue that is stored in a postprocessor This is a pointer so that the method fo...
Definition: EigenKernel.h:47
unsigned int number() const
Get variable number coming from libMesh.
virtual Real computeQpOffDiagJacobian(unsigned int)
For coupling standard variables.
Definition: Kernel.h:56
bool activeOnOld()
Return if eigen kernels should be on old solution.
const Real & eigenvalueOld()
The old eigenvalue used by inverse power iterations.
virtual void precalculateResidual()
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
bool computingScalingJacobian() const
Whether we are computing an initial Jacobian for automatic variable scaling.
Definition: SystemBase.C:1481
unsigned int m() const
const PostprocessorValue & getPostprocessorValueOldByName(const PostprocessorName &name) const
const PostprocessorValue & getPostprocessorValue(const std::string &param_name, const unsigned int index=0) const
doco-normal-methods-begin Retrieve the value of a Postprocessor or one of it&#39;s old or older values ...
bool _has_diag_save_in
The aux variables to save the diagonal Jacobian contributions to.
Definition: KernelBase.h:68
DenseVector< Real > getJacobianDiagonal(DenseMatrix< Number > &ke)
Definition: Assembly.h:1840
DenseMatrix< Number > _local_ke
Holds local Jacobian entries as they are accumulated by this Kernel.
static InputParameters validParams()
Definition: EigenKernel.C:23
const MooseVariableFieldBase & getVariable(unsigned int jvar_num) const
Retrieve the variable object from our system associated with jvar_num.
SystemBase & _sys
Reference to the EquationSystem object.
void registerBase(const std::string &value)
This method must be called from every base "Moose System" to create linkage with the Action System...
virtual void computeResidual() override
Compute this Kernel&#39;s contribution to the residual.
Definition: EigenKernel.C:80
unsigned int size() const
The number of elements that can currently be stored in the array.
Definition: MooseArray.h:256
const VariableTestValue & _test
the current test function
Definition: Kernel.h:75
virtual void precalculateJacobian()
std::vector< MooseVariableFEBase * > _diag_save_in
Definition: KernelBase.h:69
This class provides reusable routines for eigenvalue executioners.
const QBase *const & _qrule
active quadrature rule
Definition: KernelBase.h:48
bool _has_save_in
The aux variables to save the residual contributions to.
Definition: KernelBase.h:63
virtual Real computeQpJacobian()
Compute this Kernel&#39;s contribution to the Jacobian at the current quadrature point.
Definition: Kernel.h:51
const T & getParam(const std::string &name) const
Retrieve a parameter for the object.
virtual bool enabled() const override
Return the enabled status of the object.
Definition: EigenKernel.C:183
unsigned int _i
current index for the test function
Definition: KernelBase.h:57
virtual void computeOffDiagJacobian(unsigned int jvar) override
Computes d-residual / d-jvar... storing the result in Ke.
Definition: EigenKernel.C:129
MooseApp & _app
The MOOSE application this is associated with.
Definition: MooseBase.h:69
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
bool isDefaultPostprocessorValue(const std::string &param_name, const unsigned int index=0) const
Determine whether or not the Postprocessor is a default value.
Executioner * getExecutioner() const
Retrieve the Executioner for this App.
Definition: MooseApp.C:1482
const MooseArray< Real > & _coord
The scaling factor to convert from cartesian to another coordinate system (e.g rz, spherical, etc.)
Definition: KernelBase.h:54
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
virtual const PostprocessorValue & getPostprocessorValueByName(const PostprocessorName &name) const
Retrieve the value of the Postprocessor.
unsigned int _j
current index for the shape function
Definition: KernelBase.h:60
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
Definition: Kernel.h:15
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an option parameter and a documentation string to the InputParameters object...
virtual void computeJacobian() override
Compute this Kernel&#39;s contribution to the diagonal Jacobian entries.
Definition: EigenKernel.C:102
Eigen::Matrix< Real, Eigen::Dynamic, 1 > RealEigenVector
Definition: MooseTypes.h:138
unsigned int n() const
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual according to active tags.
void prepareMatrixTag(Assembly &assembly, unsigned int ivar, unsigned int jvar)
Prepare data for computing element jacobian according to the active tags.
const VariablePhiValue & _phi
the current shape functions
Definition: Kernel.h:81
EigenKernel(const InputParameters &parameters)
Definition: EigenKernel.C:34
bool _is_implicit
If the object is using implicit or explicit form.
virtual bool enabled() const
Return the enabled status of the object.
Definition: MooseObject.h:49
const PostprocessorName & getPostprocessorName(const std::string &param_name, const unsigned int index=0) const
Get the name of a postprocessor.
bool _eigen
flag for as an eigen kernel or a normal kernel
Definition: EigenKernel.h:38
unsigned int _qp
The current quadrature point index.
Definition: KernelBase.h:42