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