https://mooseframework.inl.gov
KokkosVariableTimeIntegrationAux.h
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 #pragma once
11 
12 #include "KokkosAuxKernel.h"
13 
20 {
21 public:
23 
25 
26  KOKKOS_FUNCTION Real computeValue(const unsigned int qp, AssemblyDatum & datum) const;
27 
28 protected:
29  KOKKOS_FUNCTION Real getIntegralValue(const unsigned int qp, AssemblyDatum & datum) const;
30 
32  const Real _coef;
33  const unsigned int _order;
35 
40 };
41 
42 KOKKOS_FUNCTION inline Real
44 {
45  Real integral = getIntegralValue(qp, datum);
46 
47  if (_order == 3)
48  return _u_older(datum, qp) + _coef * integral;
49 
50  return _u_old(datum, qp) + _coef * integral;
51 }
52 
53 KOKKOS_FUNCTION inline Real
55  AssemblyDatum & datum) const
56 {
57  Real integral_value = 0.0;
58 
59  for (unsigned int i = 0; i < _order; ++i)
60  integral_value += _integration_coef[i] * _coupled_vars[i](datum, qp) * _dt;
61 
69  if (_order == 3 && _dt != _dt_old)
70  {
71  Real x0 = 0.0;
72  Real x1 = _dt_old;
73  Real x2 = _dt + _dt_old;
74  Real y0 = _coupled_vars[2](datum, qp);
75  Real y1 = _coupled_vars[1](datum, qp);
76  Real y2 = _coupled_vars[0](datum, qp);
77  Real term1 = (x2 - x0) * (y0 + (x2 - x0) * (y1 - y0) / (2.0 * (x1 - x0)));
78  Real term2 = (2.0 * x2 * x2 - x0 * x2 - x0 * x0 + 3.0 * x0 * x1 - 3.0 * x1 * x2) / 6.0;
79  Real term3 = (y2 - y1) / (x2 - x1) - (y1 - y0) / (x1 - x0);
80  integral_value = term1 + term2 * term3;
81  }
82 
83  return integral_value;
84 }
Moose::Kokkos::Array< Moose::Kokkos::VariableValue > _coupled_vars
const Moose::Kokkos::VariableValue _u_old
The old variable value (zero if order == 3)
The base class for a user to derive their own Kokkos auxiliary kernels.
const InputParameters & parameters() const
Get the parameters of the object.
Definition: MooseBase.h:131
KokkosVariableTimeIntegrationAux(const InputParameters &parameters)
const Moose::Kokkos::VariableValue _u_older
The older variable value (zero if order != 3)
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
KOKKOS_FUNCTION Real computeValue(const unsigned int qp, AssemblyDatum &datum) const
An AuxKernel that can be used to integrate a field variable in time using a variety of different inte...
KOKKOS_FUNCTION Real getIntegralValue(const unsigned int qp, AssemblyDatum &datum) const
Scalar< Real > _dt_old
Size of the old time step.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
The Kokkos object that holds thread-private data in the parallel operations of Kokkos kernels...
Definition: KokkosDatum.h:244
The Kokkos wrapper classes for MOOSE-like variable value access.
Scalar< Real > _dt
Time step size.
static InputParameters validParams()