LCOV - code coverage report
Current view: top level - include/kokkos/bcs - KokkosIntegratedBC.h (source / functions) Hit Total Coverage
Test: idaholab/moose framework: 6f668f Lines: 82 87 94.3 %
Date: 2025-09-22 20:01:15 Functions: 63 102 61.8 %
Legend: Lines: hit not hit

          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 "KokkosIntegratedBCBase.h"
      13             : 
      14             : namespace Moose
      15             : {
      16             : namespace Kokkos
      17             : {
      18             : 
      19             : /**
      20             :  * The base class for a user to derive their own Kokkos integrated boundary conditions.
      21             :  *
      22             :  * The polymorphic design of the original MOOSE is reproduced statically by leveraging the Curiously
      23             :  * Recurring Template Pattern (CRTP), a programming idiom that involves a class template inheriting
      24             :  * from a template instantiation of itself. When the user derives their Kokkos object from this
      25             :  * class, the inheritance structure will look like:
      26             :  *
      27             :  * class UserIntegratedBC final : public Moose::Kokkos::IntegratedBC<UserIntegratedBC>
      28             :  *
      29             :  * It is important to note that the template argument should point to the last derived class.
      30             :  * Therefore, if the user wants to define a derived class that can be further inherited, the derived
      31             :  * class should be a class template as well. Otherwise, it is recommended to mark the derived class
      32             :  * as final to prevent its inheritence by mistake.
      33             :  *
      34             :  * The user is expected to define computeQpResidual(), computeQpJacobian(), and
      35             :  * computeQpOffDiagJacobian() as inlined public methods in their derived class (not virtual
      36             :  * override). The signature of computeQpResidual() expected to be defined in the derived class is as
      37             :  * follows:
      38             :  *
      39             :  * @param i The element-local DOF index
      40             :  * @param qp The local quadrature point index
      41             :  * @param datum The ResidualDatum object of the current thread
      42             :  * @returns The residual contribution
      43             :  *
      44             :  * KOKKOS_FUNCTION Real computeQpResidual(const unsigned int i,
      45             :  *                                        const unsigned int qp,
      46             :  *                                        ResidualDatum & datum) const;
      47             :  *
      48             :  * The signatures of computeQpJacobian() and computeOffDiagQpJacobian() can be found in the code
      49             :  * below, and their definition in the derived class is optional. If they are defined in the derived
      50             :  * class, they will hide the default definitions in the base class.
      51             :  */
      52             : template <typename Derived>
      53             : class IntegratedBC : public IntegratedBCBase
      54             : {
      55             : public:
      56             :   static InputParameters validParams();
      57             : 
      58             :   /**
      59             :    * Constructor
      60             :    */
      61             :   IntegratedBC(const InputParameters & parameters);
      62             : 
      63             :   /**
      64             :    * Dispatch residual calculation
      65             :    */
      66             :   virtual void computeResidual() override;
      67             :   /**
      68             :    * Dispatch diagonal and off-diagonal Jacobian calculation
      69             :    */
      70             :   virtual void computeJacobian() override;
      71             : 
      72             :   /**
      73             :    * Default methods to prevent compile errors even when these methods were not defined in the
      74             :    * derived class
      75             :    */
      76             :   ///@{
      77             :   /**
      78             :    * Compute diagonal Jacobian contribution on a quadrature point
      79             :    * @param i The test function DOF index
      80             :    * @param j The trial function DOF index
      81             :    * @param qp The local quadrature point index
      82             :    * @param datum The ResidualDatum object of the current thread
      83             :    * @returns The diagonal Jacobian contribution
      84             :    */
      85           0 :   KOKKOS_FUNCTION Real computeQpJacobian(const unsigned int /* i */,
      86             :                                          const unsigned int /* j */,
      87             :                                          const unsigned int /* qp */,
      88             :                                          ResidualDatum & /* datum */) const
      89             :   {
      90           0 :     return 0;
      91             :   }
      92             :   /**
      93             :    * Compute off-diagonal Jacobian contribution on a quadrature point
      94             :    * @param i The test function DOF index
      95             :    * @param j The trial function DOF index
      96             :    * @param jvar The variable number for column
      97             :    * @param qp The local quadrature point index
      98             :    * @param datum The ResidualDatum object of the current thread
      99             :    * @returns The off-diagonal Jacobian contribution
     100             :    */
     101           0 :   KOKKOS_FUNCTION Real computeQpOffDiagJacobian(const unsigned int /* i */,
     102             :                                                 const unsigned int /* j */,
     103             :                                                 const unsigned int /* jvar */,
     104             :                                                 const unsigned int /* qp */,
     105             :                                                 ResidualDatum & /* datum */) const
     106             :   {
     107           0 :     return 0;
     108             :   }
     109             :   ///@}
     110             : 
     111             :   /**
     112             :    * The parallel computation entry functions called by Kokkos
     113             :    */
     114             :   ///@{
     115             :   KOKKOS_FUNCTION void operator()(ResidualLoop, const ThreadID tid) const;
     116             :   KOKKOS_FUNCTION void operator()(JacobianLoop, const ThreadID tid) const;
     117             :   KOKKOS_FUNCTION void operator()(OffDiagJacobianLoop, const ThreadID tid) const;
     118             :   ///@}
     119             : 
     120             :   /**
     121             :    * The parallel computation bodies that can be customized in the derived class by defining
     122             :    * them in the derived class with the same signature.
     123             :    * Make sure to define them as inlined public methods if to be defined in the derived class.
     124             :    */
     125             :   ///@{
     126             :   /**
     127             :    * Compute residual
     128             :    * @param bc The boundary condition object of the final derived type
     129             :    * @param datum The ResidualDatum object of the current thread
     130             :    */
     131             :   KOKKOS_FUNCTION void computeResidualInternal(const Derived * bc, ResidualDatum & datum) const;
     132             :   /**
     133             :    * Compute diagonal Jacobian
     134             :    * @param bc The boundary condition object of the final derived type
     135             :    * @param datum The ResidualDatum object of the current thread
     136             :    */
     137             :   KOKKOS_FUNCTION void computeJacobianInternal(const Derived * bc, ResidualDatum & datum) const;
     138             :   /**
     139             :    * Compute off-diagonal Jacobian
     140             :    * @param bc The boundary condition object of the final derived type
     141             :    * @param datum The ResidualDatum object of the current thread
     142             :    */
     143             :   KOKKOS_FUNCTION void computeOffDiagJacobianInternal(const Derived * bc,
     144             :                                                       ResidualDatum & datum) const;
     145             :   ///@}
     146             : 
     147             : protected:
     148             :   /**
     149             :    * Current test function
     150             :    */
     151             :   const VariableTestValue _test;
     152             :   /**
     153             :    * Gradient of the current test function
     154             :    */
     155             :   const VariableTestGradient _grad_test;
     156             :   /**
     157             :    * Current shape function
     158             :    */
     159             :   const VariablePhiValue _phi;
     160             :   /**
     161             :    * Gradient of the current shape function
     162             :    */
     163             :   const VariablePhiGradient _grad_phi;
     164             :   /**
     165             :    * Current solution at quadrature points
     166             :    */
     167             :   const VariableValue _u;
     168             :   /**
     169             :    * Gradient of the current solution at quadrature points
     170             :    */
     171             :   const VariableGradient _grad_u;
     172             : 
     173             : protected:
     174             :   /**
     175             :    * Get whether computeQpJacobian() was not defined in the derived class
     176             :    * @returns Whether computeQpJacobian() was not defined in the derived class
     177             :    */
     178         730 :   virtual bool defaultJacobian() const
     179             :   {
     180         730 :     return &Derived::computeQpJacobian == &IntegratedBC::computeQpJacobian;
     181             :   }
     182             :   /**
     183             :    * Get whether computeQpOffDiagJacobian() was not defined in the derived class
     184             :    * @returns Whether computeQpOffDiagJacobian() was not defined in the derived class
     185             :    */
     186         730 :   virtual bool defaultOffDiagJacobian() const
     187             :   {
     188         730 :     return &Derived::computeQpOffDiagJacobian == &IntegratedBC::computeQpOffDiagJacobian;
     189             :   }
     190             : };
     191             : 
     192             : template <typename Derived>
     193             : InputParameters
     194       55058 : IntegratedBC<Derived>::validParams()
     195             : {
     196       55058 :   InputParameters params = IntegratedBCBase::validParams();
     197       55058 :   return params;
     198             : }
     199             : 
     200             : template <typename Derived>
     201         201 : IntegratedBC<Derived>::IntegratedBC(const InputParameters & parameters)
     202             :   : IntegratedBCBase(parameters, Moose::VarFieldType::VAR_FIELD_STANDARD),
     203             :     _test(),
     204             :     _grad_test(),
     205             :     _phi(),
     206             :     _grad_phi(),
     207         154 :     _u(_var),
     208         154 :     _grad_u(_var)
     209             : {
     210         201 :   addMooseVariableDependency(&_var);
     211         201 : }
     212             : 
     213             : template <typename Derived>
     214             : void
     215        8336 : IntegratedBC<Derived>::computeResidual()
     216             : {
     217        8336 :   ::Kokkos::RangePolicy<ResidualLoop, ExecSpace, ::Kokkos::IndexType<ThreadID>> policy(
     218             :       0, numKokkosBoundarySides());
     219        8336 :   ::Kokkos::parallel_for(policy, *static_cast<Derived *>(this));
     220        8336 :   ::Kokkos::fence();
     221        8336 : }
     222             : 
     223             : template <typename Derived>
     224             : void
     225         730 : IntegratedBC<Derived>::computeJacobian()
     226             : {
     227         730 :   if (!defaultJacobian())
     228             :   {
     229         222 :     ::Kokkos::RangePolicy<JacobianLoop, ExecSpace, ::Kokkos::IndexType<ThreadID>> policy(
     230             :         0, numKokkosBoundarySides());
     231         222 :     ::Kokkos::parallel_for(policy, *static_cast<Derived *>(this));
     232         222 :     ::Kokkos::fence();
     233         222 :   }
     234             : 
     235         730 :   if (!defaultOffDiagJacobian())
     236             :   {
     237          74 :     auto & sys = kokkosSystem(_kokkos_var.sys());
     238             : 
     239         148 :     _thread.resize({sys.getCoupling(_kokkos_var.var()).size(), numKokkosBoundarySides()});
     240             : 
     241          74 :     ::Kokkos::RangePolicy<OffDiagJacobianLoop, ExecSpace, ::Kokkos::IndexType<ThreadID>> policy(
     242             :         0, _thread.size());
     243          74 :     ::Kokkos::parallel_for(policy, *static_cast<Derived *>(this));
     244          74 :     ::Kokkos::fence();
     245          74 :   }
     246         730 : }
     247             : 
     248             : template <typename Derived>
     249             : KOKKOS_FUNCTION void
     250       43230 : IntegratedBC<Derived>::operator()(ResidualLoop, const ThreadID tid) const
     251             : {
     252       43230 :   auto bc = static_cast<const Derived *>(this);
     253       43230 :   auto [elem, side] = kokkosBoundaryElementSideID(tid);
     254             : 
     255       86460 :   ResidualDatum datum(
     256       43230 :       elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, _kokkos_var.var());
     257             : 
     258       43230 :   bc->computeResidualInternal(bc, datum);
     259       43230 : }
     260             : 
     261             : template <typename Derived>
     262             : KOKKOS_FUNCTION void
     263        1340 : IntegratedBC<Derived>::operator()(JacobianLoop, const ThreadID tid) const
     264             : {
     265        1340 :   auto bc = static_cast<const Derived *>(this);
     266        1340 :   auto [elem, side] = kokkosBoundaryElementSideID(tid);
     267             : 
     268        2680 :   ResidualDatum datum(
     269        1340 :       elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, _kokkos_var.var());
     270             : 
     271        1340 :   bc->computeJacobianInternal(bc, datum);
     272        1340 : }
     273             : 
     274             : template <typename Derived>
     275             : KOKKOS_FUNCTION void
     276          20 : IntegratedBC<Derived>::operator()(OffDiagJacobianLoop, const ThreadID tid) const
     277             : {
     278          20 :   auto bc = static_cast<const Derived *>(this);
     279          20 :   auto [elem, side] = kokkosBoundaryElementSideID(_thread(tid, 1));
     280             : 
     281          20 :   auto & sys = kokkosSystem(_kokkos_var.sys());
     282          20 :   auto jvar = sys.getCoupling(_kokkos_var.var())[_thread(tid, 0)];
     283             : 
     284          20 :   if (!sys.isVariableActive(jvar, kokkosMesh().getElementInfo(elem).subdomain))
     285           0 :     return;
     286             : 
     287          20 :   ResidualDatum datum(elem, side, kokkosAssembly(), kokkosSystems(), _kokkos_var, jvar);
     288             : 
     289          20 :   bc->computeOffDiagJacobianInternal(bc, datum);
     290             : }
     291             : 
     292             : template <typename Derived>
     293             : KOKKOS_FUNCTION void
     294       43230 : IntegratedBC<Derived>::computeResidualInternal(const Derived * bc, ResidualDatum & datum) const
     295             : {
     296       43230 :   ResidualObject::computeResidualInternal(
     297             :       datum,
     298       86460 :       [&](Real * local_re, const unsigned int ib, const unsigned int ie)
     299             :       {
     300      129184 :         for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
     301             :         {
     302       85954 :           datum.reinit();
     303             : 
     304      428758 :           for (unsigned int i = ib; i < ie; ++i)
     305      342804 :             local_re[i] += datum.JxW(qp) * bc->computeQpResidual(i, qp, datum);
     306             :         }
     307             :       });
     308       43230 : }
     309             : 
     310             : template <typename Derived>
     311             : KOKKOS_FUNCTION void
     312        1340 : IntegratedBC<Derived>::computeJacobianInternal(const Derived * bc, ResidualDatum & datum) const
     313             : {
     314        1340 :   ResidualObject::computeJacobianInternal(
     315             :       datum,
     316        2680 :       [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
     317             :       {
     318        4020 :         for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
     319             :         {
     320        2680 :           datum.reinit();
     321             : 
     322       45560 :           for (unsigned int ij = ijb; ij < ije; ++ij)
     323             :           {
     324       42880 :             unsigned int i = ij % datum.n_jdofs();
     325       42880 :             unsigned int j = ij / datum.n_jdofs();
     326             : 
     327       42880 :             local_ke[ij] += datum.JxW(qp) * bc->computeQpJacobian(i, j, qp, datum);
     328             :           }
     329             :         }
     330             :       });
     331        1340 : }
     332             : 
     333             : template <typename Derived>
     334             : KOKKOS_FUNCTION void
     335          20 : IntegratedBC<Derived>::computeOffDiagJacobianInternal(const Derived * bc,
     336             :                                                       ResidualDatum & datum) const
     337             : {
     338          20 :   ResidualObject::computeJacobianInternal(
     339             :       datum,
     340          40 :       [&](Real * local_ke, const unsigned int ijb, const unsigned int ije)
     341             :       {
     342          40 :         for (unsigned int qp = 0; qp < datum.n_qps(); ++qp)
     343             :         {
     344          20 :           datum.reinit();
     345             : 
     346         100 :           for (unsigned int ij = ijb; ij < ije; ++ij)
     347             :           {
     348          80 :             unsigned int i = ij % datum.n_jdofs();
     349          80 :             unsigned int j = ij / datum.n_jdofs();
     350             : 
     351          80 :             local_ke[ij] +=
     352          80 :                 datum.JxW(qp) * bc->computeQpOffDiagJacobian(i, j, datum.jvar(), qp, datum);
     353             :           }
     354             :         }
     355             :       });
     356          20 : }
     357             : 
     358             : } // namespace Kokkos
     359             : } // namespace Moose
     360             : 
     361             : #define usingKokkosIntegratedBCMembers(T)                                                          \
     362             :   usingKokkosIntegratedBCBaseMembers;                                                              \
     363             :                                                                                                    \
     364             : protected:                                                                                         \
     365             :   using Moose::Kokkos::IntegratedBC<T>::_test;                                                     \
     366             :   using Moose::Kokkos::IntegratedBC<T>::_grad_test;                                                \
     367             :   using Moose::Kokkos::IntegratedBC<T>::_phi;                                                      \
     368             :   using Moose::Kokkos::IntegratedBC<T>::_grad_phi;                                                 \
     369             :   using Moose::Kokkos::IntegratedBC<T>::_u;                                                        \
     370             :   using Moose::Kokkos::IntegratedBC<T>::_grad_u;                                                   \
     371             :                                                                                                    \
     372             : public:                                                                                            \
     373             :   using Moose::Kokkos::IntegratedBC<T>::operator()

Generated by: LCOV version 1.14