Line data Source code
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 "KokkosTimeKernel.h" 13 : 14 : class KokkosTimeDerivative : public Moose::Kokkos::TimeKernel 15 : { 16 : public: 17 : static InputParameters validParams(); 18 : 19 : KokkosTimeDerivative(const InputParameters & parameters); 20 : 21 : template <typename Derived> 22 : KOKKOS_FUNCTION void computeJacobianInternal(const Derived & kernel, ResidualDatum & datum) const; 23 : 24 : KOKKOS_FUNCTION Real computeQpResidual(const unsigned int i, 25 : const unsigned int qp, 26 : ResidualDatum & datum) const; 27 : KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int i, 28 : const unsigned int j, 29 : const unsigned int qp, 30 : ResidualDatum & datum) const; 31 : 32 : protected: 33 : const bool _lumping; 34 : }; 35 : 36 : template <typename Derived> 37 : KOKKOS_FUNCTION void 38 192050 : KokkosTimeDerivative::computeJacobianInternal(const Derived & kernel, ResidualDatum & datum) const 39 : { 40 : using Moose::Kokkos::MAX_CACHED_DOF; 41 : 42 : Real local_ke[MAX_CACHED_DOF]; 43 : 44 192050 : unsigned int num_batches = datum.n_idofs() * datum.n_jdofs() / MAX_CACHED_DOF; 45 : 46 192050 : if ((datum.n_idofs() * datum.n_jdofs()) % MAX_CACHED_DOF) 47 192050 : ++num_batches; 48 : 49 407908 : for (unsigned int batch = 0; batch < num_batches; ++batch) 50 : { 51 215858 : unsigned int ijb = batch * MAX_CACHED_DOF; 52 215858 : unsigned int ije = ::Kokkos::min(ijb + MAX_CACHED_DOF, datum.n_idofs() * datum.n_jdofs()); 53 : 54 2964610 : for (unsigned int ij = ijb; ij < ije; ++ij) 55 2748752 : local_ke[ij - ijb] = 0; 56 : 57 1072898 : for (unsigned int qp = 0; qp < datum.n_qps(); ++qp) 58 : { 59 857040 : datum.reinit(); 60 : 61 14302512 : for (unsigned int ij = ijb; ij < ije; ++ij) 62 : { 63 13445472 : unsigned int i = ij % datum.n_jdofs(); 64 13445472 : unsigned int j = ij / datum.n_jdofs(); 65 : 66 13445472 : local_ke[ij - ijb] += datum.JxW(qp) * kernel.computeQpJacobian(i, j, qp, datum); 67 : } 68 : } 69 : 70 2964610 : for (unsigned int ij = ijb; ij < ije; ++ij) 71 : { 72 2748752 : unsigned int i = ij % datum.n_jdofs(); 73 2748752 : unsigned int j = ij / datum.n_jdofs(); 74 : 75 5497504 : accumulateTaggedElementalMatrix( 76 5497504 : local_ke[ij - ijb], datum.elem().id, i, _lumping ? i : j, datum.jvar()); 77 : } 78 : } 79 192050 : } 80 : 81 : KOKKOS_FUNCTION inline Real 82 17933936 : KokkosTimeDerivative::computeQpResidual(const unsigned int i, 83 : const unsigned int qp, 84 : ResidualDatum & datum) const 85 : { 86 17933936 : return _test(datum, i, qp) * _u_dot(datum, qp); 87 : } 88 : 89 : KOKKOS_FUNCTION inline Real 90 13445472 : KokkosTimeDerivative::computeQpJacobian(const unsigned int i, 91 : const unsigned int j, 92 : const unsigned int qp, 93 : ResidualDatum & datum) const 94 : { 95 13445472 : return _test(datum, i, qp) * _phi(datum, j, qp) * _du_dot_du; 96 : }