Line data Source code
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 : #pragma once 11 : 12 : #include "KokkosKernel.h" 13 : 14 : namespace Moose::Kokkos 15 : { 16 : 17 : /** 18 : * The base class for a user to derive their own Kokkos kernels where the residual is of the form 19 : * 20 : * $(\dots, \psi_i)$ 21 : * 22 : * i.e. the test function $(\psi_i)$ can be factored out for optimization. 23 : * 24 : * The user should still define computeQpResidual(), computeQpJacobian(), and 25 : * computeQpOffDiagJacobian(), but their signatures are different from the base class. The signature 26 : * of computeQpResidual() expected to be defined in the derived class is as follows: 27 : * 28 : * @tparam Derived The object type 29 : * @param qp The local quadrature point index 30 : * @param datum The AssemblyDatum object of the current thread 31 : * @returns The component of the residual contribution that will be multiplied by the test function 32 : * 33 : * template <typename Derived> 34 : * KOKKOS_FUNCTION Real computeQpResidual(const unsigned int qp, 35 : * AssemblyDatum & datum) const; 36 : * 37 : * The signature of computeQpJacobian() and computeQpOffDiagJacobian() can be found in the code 38 : * below. 39 : */ 40 : class KernelValue : public Kernel 41 : { 42 : public: 43 : static InputParameters validParams(); 44 : 45 : /** 46 : * Constructor 47 : */ 48 : KernelValue(const InputParameters & parameters); 49 : 50 : /** 51 : * Default methods to prevent compile errors even when these methods were not defined in the 52 : * derived class 53 : */ 54 : ///@{ 55 : /** 56 : * Compute diagonal Jacobian contribution on a quadrature point 57 : * @tparam Derived The object type 58 : * @param j The trial function DOF index 59 : * @param qp The local quadrature point index 60 : * @param datum The AssemblyDatum object of the current thread 61 : * @returns The component of the Jacobian contribution that will be multiplied by the test 62 : * function 63 : */ 64 : template <typename Derived> 65 0 : KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int /* j */, 66 : const unsigned int /* qp */, 67 : AssemblyDatum & /* datum */) const 68 : { 69 0 : ::Kokkos::abort("Default computeQpJacobian() should never be called. Make sure you properly " 70 : "redefined this method in your class without typos."); 71 : 72 : return 0; 73 : } 74 : /** 75 : * Compute off-diagonal Jacobian contribution on a quadrature point 76 : * @tparam Derived The object type 77 : * @param j The trial function DOF index 78 : * @param jvar The variable number for column 79 : * @param qp The local quadrature point index 80 : * @param datum The AssemblyDatum object of the current thread 81 : * @returns The component of the off-diagonal Jacobian contribution that will be multiplied by the 82 : * test function 83 : */ 84 : template <typename Derived> 85 0 : KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int /* j */, 86 : const unsigned int /* jvar */, 87 : const unsigned int /* qp */, 88 : AssemblyDatum & /* datum */) const 89 : { 90 0 : ::Kokkos::abort( 91 : "Default computeQpOffDiagJacobian() should never be called. Make sure you properly " 92 : "redefined this method in your class without typos."); 93 : 94 : return 0; 95 : } 96 : ///@} 97 : 98 : /** 99 : * Functions used to check if users have overriden the hook methods, whose calculations can be 100 : * skipped when not overriden 101 : * @returns The function pointer of the default hook method 102 : */ 103 : ///@{ 104 : template <typename Derived> 105 242080 : static auto defaultJacobian() 106 : { 107 242080 : return &KernelValue::computeQpJacobian<Derived>; 108 : } 109 : template <typename Derived> 110 242080 : static auto defaultOffDiagJacobian() 111 : { 112 242080 : return &KernelValue::computeQpOffDiagJacobian<Derived>; 113 : } 114 : ///@} 115 : 116 : /** 117 : * The parallel computation bodies that hide the base class methods to optimize for factoring 118 : * out the test function 119 : */ 120 : ///@{ 121 : template <typename Derived> 122 : KOKKOS_FUNCTION void computeResidualInternal(const Derived & kernel, AssemblyDatum & datum) const; 123 : template <typename Derived> 124 : KOKKOS_FUNCTION void computeJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const; 125 : template <typename Derived> 126 : KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived & kernel, 127 : AssemblyDatum & datum) const; 128 : ///@} 129 : }; 130 : 131 : template <typename Derived> 132 : KOKKOS_FUNCTION void 133 598588 : KernelValue::computeResidualInternal(const Derived & kernel, AssemblyDatum & datum) const 134 : { 135 598588 : ResidualObject::computeResidualInternal( 136 : datum, 137 1197176 : [&](Real * local_re, const unsigned int ib, const unsigned int ie) 138 : { 139 2243180 : for (unsigned int qp = 0; qp < datum.n_qps(); ++qp) 140 : { 141 1644592 : Real value = datum.JxW(qp) * kernel.template computeQpResidual<Derived>(qp, datum); 142 : 143 2914304 : for (unsigned int i = ib; i < ie; ++i) 144 1269712 : local_re[i] += value * _test(datum, i, qp); 145 : } 146 : }); 147 598588 : } 148 : 149 : template <typename Derived> 150 : KOKKOS_FUNCTION void 151 96 : KernelValue::computeJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const 152 : { 153 96 : ResidualObject::computeJacobianInternal( 154 : datum, 155 192 : [&](Real * local_ke, const unsigned int ib, const unsigned int ie, const unsigned int j) 156 : { 157 480 : for (unsigned int qp = 0; qp < datum.n_qps(); ++qp) 158 : { 159 384 : Real value = datum.JxW(qp) * kernel.template computeQpJacobian<Derived>(j, qp, datum); 160 : 161 1920 : for (unsigned int i = ib; i < ie; ++i) 162 1536 : local_ke[i] += value * _test(datum, i, qp); 163 : } 164 : }); 165 96 : } 166 : 167 : template <typename Derived> 168 : KOKKOS_FUNCTION void 169 5000 : KernelValue::computeOffDiagJacobianInternal(const Derived & kernel, AssemblyDatum & datum) const 170 : { 171 5000 : ResidualObject::computeJacobianInternal( 172 : datum, 173 10000 : [&](Real * local_ke, const unsigned int ib, const unsigned int ie, const unsigned int j) 174 : { 175 25000 : for (unsigned int qp = 0; qp < datum.n_qps(); ++qp) 176 : { 177 20000 : Real value = datum.JxW(qp) * kernel.template computeQpOffDiagJacobian<Derived>( 178 : j, datum.jvar(), qp, datum); 179 : 180 100000 : for (unsigned int i = ib; i < ie; ++i) 181 80000 : local_ke[i] += value * _test(datum, i, qp); 182 : } 183 : }); 184 5000 : } 185 : 186 : } // namespace Moose::Kokkos