Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Fix nightly GPU failures after kokkos 4.5.0 #1100

Merged
merged 1 commit into from
Dec 14, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
31 changes: 0 additions & 31 deletions src/PHAL_Utilities_Def.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -260,35 +260,4 @@ void scale (ArrayT& a, const T& val) {
loop(sl, a);
}

template< class T , class ... P >
inline
typename std::enable_if<
!Kokkos::is_dynrankview_fad<Kokkos::DynRankView<T,P...>>::value,
typename Kokkos::DynRankView<T,P...>::non_const_type >::type
create_copy( const std::string& name,
const Kokkos::DynRankView<T,P...> & src )
{
using dst_type = typename Kokkos::DynRankView<T,P...>::non_const_type;
auto layout = Kokkos::Impl::reconstructLayout(src.layout(), src.rank());
return dst_type( name , layout );
}

template< class T , class ... P >
inline
typename std::enable_if<
Kokkos::is_dynrankview_fad<Kokkos::DynRankView<T,P...>>::value,
typename Kokkos::DynRankView<T,P...>::non_const_type >::type
create_copy( const std::string& /* name */,
const Kokkos::DynRankView<T,P...> & src )
{
using Dst = typename Kokkos::DynRankView<T,P...>::non_const_type;
auto sm = src.impl_map();
auto sl = sm.layout();
auto fad_rank = src.rank();
sl.dimension[fad_rank] = sm.dimension_scalar();
auto real_rank = fad_rank + 1;
auto ml = Kokkos::Impl::reconstructLayout(sl, real_rank);
return Dst ( src.label(), ml );
}

} // namespace PHAL
1 change: 1 addition & 0 deletions src/corePDEs/evaluators/PHAL_HeatEqResid.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -60,6 +60,7 @@ class HeatEqResid : public PHX::EvaluatorWithBaseImpl<Traits>,
bool haverhoCp;
unsigned int numQPs, numDims, numNodes, worksetSize;
Kokkos::DynRankView<ScalarT, PHX::Device> flux;
Kokkos::DynRankView<ScalarT, PHX::Device> neg_source;
Kokkos::DynRankView<ScalarT, PHX::Device> aterm;
Kokkos::DynRankView<ScalarT, PHX::Device> convection;
};
Expand Down
9 changes: 4 additions & 5 deletions src/corePDEs/evaluators/PHAL_HeatEqResid_Def.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -107,9 +107,10 @@ postRegistrationSetup(typename Traits::SetupData /* d */,
this->utils.setFieldData(TResidual,fm);

// Allocate workspace
flux = Kokkos::createDynRankView(Temperature.get_view(), "XXX", worksetSize, numQPs, numDims);
if (haveAbsorption) aterm = Kokkos::createDynRankView(Temperature.get_view(), "XXX", worksetSize, numQPs);
if (haveConvection) convection = Kokkos::createDynRankView(Temperature.get_view(), "XXX", worksetSize, numQPs);
flux = Kokkos::createDynRankView(Temperature.get_view(), "flux", worksetSize, numQPs, numDims);
if (haveSource) neg_source = Kokkos::createDynRankView(Temperature.get_view(), "neg_source", worksetSize, numQPs);
if (haveAbsorption) aterm = Kokkos::createDynRankView(Temperature.get_view(), "aterm", worksetSize, numQPs);
if (haveConvection) convection = Kokkos::createDynRankView(Temperature.get_view(), "convection", worksetSize, numQPs);
}

//**********************************************************************
Expand All @@ -124,8 +125,6 @@ evaluateFields(typename Traits::EvalData workset)
FST::integrate(TResidual.get_view(), flux, wGradBF.get_view(), false); // "false" overwrites

if (haveSource) {
auto neg_source = PHAL::create_copy("neg_source", Source.get_view());

for (size_t i =0; i< Source.extent(0); i++)
for (size_t j =0; j< Source.extent(1); j++)
neg_source(i,j) = Source(i,j) * -1.0;
Expand Down
3 changes: 3 additions & 0 deletions src/demoPDEs/evaluators/PHAL_HelmholtzResid.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -54,7 +54,10 @@ class HelmholtzResid : public PHX::EvaluatorWithBaseImpl<Traits>,

bool haveSource;

unsigned int worksetSize, numQPs;
ScalarT ksqr;
Kokkos::DynRankView<ScalarT, PHX::Device> U_ksqr;
Kokkos::DynRankView<ScalarT, PHX::Device> V_ksqr;

// Output:
PHX::MDField<ScalarT,Cell,Node> UResidual;
Expand Down
13 changes: 10 additions & 3 deletions src/demoPDEs/evaluators/PHAL_HelmholtzResid_Def.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -61,6 +61,12 @@ HelmholtzResid(const Teuchos::ParameterList& p) :
p.get< Teuchos::RCP<ParamLib> >("Parameter Library");
this->registerSacadoParameter("Ksqr", paramLib);

Teuchos::RCP<PHX::DataLayout> qp_scalar_dl =
p.get< Teuchos::RCP<PHX::DataLayout> >("QP Scalar Data Layout");
std::vector<PHX::DataLayout::size_type> dims;
qp_scalar_dl->dimensions(dims);
worksetSize = dims[0];
numQPs = dims[1];
}

//**********************************************************************
Expand All @@ -81,6 +87,10 @@ postRegistrationSetup(typename Traits::SetupData /* d */,
}
this->utils.setFieldData(UResidual,fm);
this->utils.setFieldData(VResidual,fm);

// Allocate workspace
U_ksqr = Kokkos::createDynRankView(U.get_view(), "U_ksqr", worksetSize, numQPs);
V_ksqr = Kokkos::createDynRankView(U.get_view(), "V_ksqr", worksetSize, numQPs);
}

//**********************************************************************
Expand All @@ -102,9 +112,6 @@ evaluateFields(typename Traits::EvalData /* workset */)
FST::integrate(VResidual.get_view(), VSource.get_view(), wBF.get_view(), true);
}

auto U_ksqr = create_copy("U_ksqr", U.get_view());
auto V_ksqr = create_copy("V_ksqr", V.get_view());

RST::scale(U_ksqr, U.get_view(), ksqr);
RST::scale(V_ksqr, V.get_view(), ksqr);

Expand Down
2 changes: 2 additions & 0 deletions src/demoPDEs/evaluators/PHAL_PoissonResid.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -44,6 +44,8 @@ class PoissonResid : public PHX::EvaluatorWithBaseImpl<Traits>,
PHX::MDField<const ScalarT> Source;

bool haveSource;
unsigned int worksetSize, numQPs;
Kokkos::DynRankView<ScalarT, PHX::Device> neg_source;

// Output:
PHX::MDField<ScalarT> PhiResidual;
Expand Down
8 changes: 6 additions & 2 deletions src/demoPDEs/evaluators/PHAL_PoissonResid_Def.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -23,7 +23,9 @@ PoissonResid(const Teuchos::ParameterList& p,
Source (p.get<std::string> ("Source Name"), dl->qp_scalar ),
haveSource (p.get<bool>("Have Source")),
PhiResidual (p.get<std::string> ("Residual Name"), dl->node_scalar ),
PhiFlux (p.get<std::string> ("Flux QP Variable Name"), dl->qp_gradient)
PhiFlux (p.get<std::string> ("Flux QP Variable Name"), dl->qp_gradient),
worksetSize (dl->qp_scalar->extent(0)),
numQPs (dl->qp_scalar->extent(1))
{
this->addDependentField(wBF);
//this->addDependentField(Potential);
Expand Down Expand Up @@ -53,6 +55,9 @@ postRegistrationSetup(typename Traits::SetupData /* d */,
if (haveSource) this->utils.setFieldData(Source,fm);

this->utils.setFieldData(PhiResidual,fm);

// Allocate workspace
if (haveSource) neg_source = Kokkos::createDynRankView(Source.get_view(), "neg_source", worksetSize, numQPs);
}

//**********************************************************************
Expand All @@ -67,7 +72,6 @@ evaluateFields(typename Traits::EvalData /* workset */)

FST::integrate(PhiResidual.get_view(), PhiFlux.get_view(), wGradBF.get_view(), false); // "false" overwrites

auto neg_source = PHAL::create_copy("neg_Source", Source.get_view());
if (haveSource) {
for (unsigned int i=0; i<Source.extent(0); i++)
for (unsigned int j=0; j<Source.extent(1); j++)
Expand Down
1 change: 1 addition & 0 deletions tests/demoPDEs/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -4,6 +4,7 @@
## in the file "license.txt" in the top-level Albany directory //
##*****************************************************************//

add_subdirectory(Helmholtz2D)
add_subdirectory(AdvDiff)
add_subdirectory(ReactDiffSystem)
add_subdirectory(ODE)
Expand Down