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 "KokkosNodalBC.h" 13 : 14 : namespace Moose 15 : { 16 : namespace Kokkos 17 : { 18 : 19 : /** 20 : * The base Kokkos boundary condition of a Dirichlet type 21 : */ 22 : template <typename Derived> 23 : class DirichletBCBase : public NodalBC<Derived> 24 : { 25 : usingKokkosNodalBCMembers(Derived); 26 : 27 : public: 28 : static InputParameters validParams(); 29 : 30 : /** 31 : * Constructor 32 : */ 33 : DirichletBCBase(const InputParameters & parameters); 34 : 35 : /** 36 : * Get whether the value is to be preset 37 : * @returns Whether the value is to be preset 38 : */ 39 1189 : virtual bool preset() const override { return _preset; } 40 : 41 : /** 42 : * Dispatch solution vector preset 43 : * @param tag The tag associated with the solution vector to be preset 44 : */ 45 : virtual void presetSolution(TagID tag) override; 46 : 47 : /** 48 : * The preset function called by Kokkos 49 : */ 50 : KOKKOS_FUNCTION void operator()(const ThreadID tid) const; 51 : 52 : /** 53 : * Compute residual contribution on a node 54 : * @param node The contiguous node ID 55 : * @returns The residual contribution 56 : */ 57 : KOKKOS_FUNCTION Real computeQpResidual(const ContiguousNodeID node) const; 58 : 59 : private: 60 : /** 61 : * Flag whether the value is to be preset 62 : */ 63 : const bool _preset; 64 : /** 65 : * Tag associated with the solution vector to be preset 66 : */ 67 : TagID _solution_tag; 68 : }; 69 : 70 : template <typename Derived> 71 : InputParameters 72 20608 : DirichletBCBase<Derived>::validParams() 73 : { 74 20608 : InputParameters params = NodalBC<Derived>::validParams(); 75 41782 : params.addParam<bool>( 76 40084 : "preset", true, "Whether or not to preset the BC (apply the value before the solve begins)."); 77 20608 : return params; 78 0 : } 79 : 80 : template <typename Derived> 81 1472 : DirichletBCBase<Derived>::DirichletBCBase(const InputParameters & parameters) 82 1812 : : NodalBC<Derived>(parameters), _preset(this->template getParam<bool>("preset")) 83 : { 84 1189 : } 85 : 86 : template <typename Derived> 87 : void 88 3362 : DirichletBCBase<Derived>::presetSolution(TagID tag) 89 : { 90 3362 : _solution_tag = tag; 91 : 92 3362 : ::Kokkos::parallel_for(::Kokkos::RangePolicy<ExecSpace, ::Kokkos::IndexType<ThreadID>>( 93 : 0, this->numKokkosBoundaryNodes()), 94 : *static_cast<Derived *>(this)); 95 3362 : } 96 : 97 : template <typename Derived> 98 : KOKKOS_FUNCTION void 99 18075 : DirichletBCBase<Derived>::operator()(const ThreadID tid) const 100 : { 101 18075 : auto bc = static_cast<const Derived *>(this); 102 18075 : auto node = kokkosBoundaryNodeID(tid); 103 18075 : auto & sys = kokkosSystem(_kokkos_var.sys()); 104 18075 : auto dof = sys.getNodeLocalDofIndex(node, _kokkos_var.var()); 105 : 106 18075 : if (dof == libMesh::DofObject::invalid_id) 107 0 : return; 108 : 109 18075 : sys.getVectorDofValue(dof, _solution_tag) = bc->computeValue(node); 110 : } 111 : 112 : template <typename Derived> 113 : KOKKOS_FUNCTION Real 114 269929 : DirichletBCBase<Derived>::computeQpResidual(const ContiguousNodeID node) const 115 : { 116 269929 : auto bc = static_cast<const Derived *>(this); 117 : 118 269929 : return _u(node) - bc->computeValue(node); 119 : } 120 : 121 : } // namespace Kokkos 122 : } // namespace Moose 123 : 124 : #define usingKokkosDirichletBCBaseMembers(T) \ 125 : usingKokkosNodalBCMembers(T); \ 126 : \ 127 : public: \ 128 : using Moose::Kokkos::DirichletBCBase<T>::operator()