https://mooseframework.inl.gov
ODEKernel.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 
10 #include "ODEKernel.h"
11 
12 // MOOSE includes
13 #include "Assembly.h"
14 #include "MooseVariableScalar.h"
15 #include "SystemBase.h"
16 
19 {
21  return params;
22 }
23 
24 ODEKernel::ODEKernel(const InputParameters & parameters) : ScalarKernel(parameters) {}
25 
26 void
28 {
29 }
30 
31 void
33 {
35 
36  for (_i = 0; _i < _var.order(); _i++)
38 
40 }
41 
42 void
44 {
46 
47  for (_i = 0; _i < _var.order(); _i++)
48  for (_j = 0; _j < _var.order(); _j++)
50 
52 
53  // compute off-diagonal jacobians wrt scalar variables
54  const std::vector<MooseVariableScalar *> & scalar_vars = _sys.getScalarVariables(_tid);
55  for (const auto & var : scalar_vars)
56  computeOffDiagJacobianScalar(var->number());
57 }
58 
59 void
61 {
63 
65  for (_i = 0; _i < _var.order(); _i++)
66  for (_j = 0; _j < var_j.order(); _j++)
67  {
68  if (jvar != _var.number())
70  }
71 
73 }
74 
75 Real
77 {
78  return 0.;
79 }
virtual void computeJacobian() override
Compute this object&#39;s contribution to the diagonal Jacobian entries.
Definition: ODEKernel.C:43
const std::vector< MooseVariableScalar * > & getScalarVariables(THREAD_ID tid)
Definition: SystemBase.h:757
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
virtual Real computeQpOffDiagJacobianScalar(unsigned int jvar)
Definition: ODEKernel.h:31
unsigned int number() const
Get variable number coming from libMesh.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
MooseVariableScalar & _var
Scalar variable.
THREAD_ID _tid
The thread ID for this kernel.
DenseMatrix< Number > _local_ke
Holds local Jacobian entries as they are accumulated by this Kernel.
SystemBase & _sys
Reference to the EquationSystem object.
ODEKernel(const InputParameters &parameters)
Definition: ODEKernel.C:24
virtual void computeOffDiagJacobianScalar(unsigned int jvar) override
Computes jacobian block with respect to a scalar variable.
Definition: ODEKernel.C:60
virtual Real computeQpResidual()
Definition: ScalarKernel.h:24
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
virtual Real computeQpJacobian() override
Definition: ODEKernel.C:76
Assembly & _assembly
Reference to this Kernel&#39;s assembly object.
libMesh::Order order() const
Get the order of this variable Note: Order enum can be implicitly converted to unsigned int...
virtual MooseVariableScalar & getScalarVariable(THREAD_ID tid, const std::string &var_name) const
Gets a reference to a scalar variable with specified number.
Definition: SystemBase.C:144
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void computeResidual() override
Compute this object&#39;s contribution to the residual.
Definition: ODEKernel.C:32
DenseVector< Number > _local_re
Holds local residual entries as they are accumulated by this Kernel.
Class for scalar variables (they are different).
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.
virtual void reinit() override
Reinitialization method called before each call to computeResidual()
Definition: ODEKernel.C:27
static InputParameters validParams()
Definition: ODEKernel.C:18
static InputParameters validParams()
Definition: ScalarKernel.C:14