www.mooseframework.org
ODEKernel.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 "ODEKernel.h"
11 
12 // MOOSE includes
13 #include "Assembly.h"
14 #include "MooseVariableScalar.h"
15 #include "SystemBase.h"
16 
17 template <>
20 {
22  return params;
23 }
24 
25 ODEKernel::ODEKernel(const InputParameters & parameters) : ScalarKernel(parameters) {}
26 
27 void
29 {
30 }
31 
32 void
34 {
36 
37  for (_i = 0; _i < _var.order(); _i++)
39 
41 }
42 
43 void
45 {
47 
48  for (_i = 0; _i < _var.order(); _i++)
49  for (_j = 0; _j < _var.order(); _j++)
51 
53 
54  // compute off-diagonal jacobians wrt scalar variables
55  const std::vector<MooseVariableScalar *> & scalar_vars = _sys.getScalarVariables(_tid);
56  for (const auto & var : scalar_vars)
57  computeOffDiagJacobian(var->number());
58 }
59 
60 void
62 {
63  if (_sys.isScalarVariable(jvar))
64  {
66 
68  for (_i = 0; _i < _var.order(); _i++)
69  for (_j = 0; _j < var_j.order(); _j++)
70  {
71  if (jvar != _var.number())
73  }
74 
76  }
77 }
78 
79 Real
81 {
82  return 0.;
83 }
84 
85 Real
87 {
88  return 0.;
89 }
virtual void computeJacobian() override
Definition: ODEKernel.C:44
InputParameters validParams< ODEKernel >()
Definition: ODEKernel.C:19
InputParameters validParams< ScalarKernel >()
Definition: ScalarKernel.C:20
const std::vector< MooseVariableScalar * > & getScalarVariables(THREAD_ID tid)
Definition: SystemBase.h:614
void accumulateTaggedLocalResidual()
Local residual blocks will be appended by adding the current local kernel residual.
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...
Assembly & _assembly
Definition: ScalarKernel.h:73
DenseMatrix< Number > _local_ke
Holds residual entries as they are accumulated by this Kernel.
virtual void computeOffDiagJacobian(unsigned int jvar) override
Definition: ODEKernel.C:61
ODEKernel(const InputParameters &parameters)
Definition: ODEKernel.C:25
MooseVariableScalar & _var
Scalar variable.
Definition: ScalarKernel.h:75
THREAD_ID _tid
Definition: ScalarKernel.h:71
void accumulateTaggedLocalMatrix()
Local Jacobian blocks will be appended by adding the current local kernel Jacobian.
SystemBase & _sys
Definition: ScalarKernel.h:69
Order order() const
Get the order of this variable Note: Order enum can be implicitly converted to unsigned int...
virtual void computeResidual() override
Definition: ODEKernel.C:33
DenseVector< Number > _local_re
Holds residual entries as they are accumulated by this Kernel.
Class for scalar variables (they are different).
virtual Real computeQpResidual()=0
virtual bool isScalarVariable(unsigned int var_name) const
Definition: SystemBase.C:686
void prepareVectorTag(Assembly &assembly, unsigned int ivar)
Prepare data for computing element residual the according to active tags.
unsigned int _i
Definition: ScalarKernel.h:78
void prepareMatrixTag(Assembly &assembly, unsigned int ivar, unsigned int jvar)
Prepare data for computing element jacobian according to the ative tags.
virtual void reinit() override
Definition: ODEKernel.C:28
virtual MooseVariableScalar & getScalarVariable(THREAD_ID tid, const std::string &var_name)
Gets a reference to a scalar variable with specified number.
Definition: SystemBase.C:149
virtual Real computeQpOffDiagJacobian(unsigned int jvar)
Definition: ODEKernel.C:86
virtual Real computeQpJacobian()
Definition: ODEKernel.C:80
unsigned int _j
Definition: ScalarKernel.h:78