12#include <sundials/sundials_nvector.h>
13#include <sundials/sundials_core.hpp>
14#include <nvector/nvector_serial.h>
15#include <nvector/nvector_kokkos.hpp>
16#include <arkode/arkode_arkstep.h>
17#include <arkode/arkode_erkstep.h>
23namespace NeoN::sundials
31inline auto SUN_CONTEXT_DELETER = [](SUNContext* ctx)
44inline auto SUN_ARK_DELETER = [](
char* ark)
48 void* arkodMem =
reinterpret_cast<void*
>(ark);
49 ARKodeFree(&arkodMem);
59inline ARKODE_ERKTableID stringToERKTable(
const std::string& key)
61 if (key ==
"Forward-Euler")
return ARKODE_FORWARD_EULER_1_1;
64 NF_ERROR_EXIT(
"Currently unsupported until field time step-stage indexing resolved.");
65 return ARKODE_HEUN_EULER_2_1_2;
67 if (key ==
"Midpoint")
69 NF_ERROR_EXIT(
"Currently unsupported until field time step-stage indexing resolved.");
70 return ARKODE_EXPLICIT_MIDPOINT_EULER_2_1_2;
73 "Unsupported Runge-Kutta time integration method selectied: " + key +
".\n"
74 +
"Supported methods are: Forward-Euler, Heun, Midpoint."
76 return ARKODE_ERK_NONE;
87template<
typename SKVectorType,
typename ValueType>
90 auto view = ::sundials::kokkos::GetVec<SKVectorType>(vector)->View();
91 auto fieldView = field.
view();
93 field.
exec(), field.
range(), KOKKOS_LAMBDA(
const localIdx i) { view(i) = fieldView[i]; }
104template<
typename ValueType>
108 if (std::holds_alternative<NeoN::GPUExecutor>(field.
exec()))
110 fieldToSunNVectorImpl<::sundials::kokkos::Vector<Kokkos::DefaultExecutionSpace>>(
115 if (std::holds_alternative<NeoN::CPUExecutor>(field.
exec()))
117 fieldToSunNVectorImpl<::sundials::kokkos::Vector<Kokkos::DefaultHostExecutionSpace>>(
122 if (std::holds_alternative<NeoN::SerialExecutor>(field.
exec()))
124 fieldToSunNVectorImpl<::sundials::kokkos::Vector<Kokkos::Serial>>(field, vector);
138template<
typename SKVectorType,
typename ValueType>
141 auto view = ::sundials::kokkos::GetVec<SKVectorType>(vector)->View();
142 ValueType* fieldData = field.
data();
144 field.
exec(), field.
range(), KOKKOS_LAMBDA(
const localIdx i) { fieldData[i] = view(i); }
154template<
typename ValueType>
157 if (std::holds_alternative<NeoN::GPUExecutor>(field.
exec()))
159 sunNVectorToVectorImpl<::sundials::kokkos::Vector<Kokkos::DefaultExecutionSpace>>(
164 if (std::holds_alternative<NeoN::CPUExecutor>(field.
exec()))
166 sunNVectorToVectorImpl<::sundials::kokkos::Vector<Kokkos::DefaultHostExecutionSpace>>(
171 if (std::holds_alternative<NeoN::SerialExecutor>(field.
exec()))
173 sunNVectorToVectorImpl<::sundials::kokkos::Vector<Kokkos::Serial>>(vector, field);
193template<
typename SolutionVectorType>
194int explicitRKSolve([[maybe_unused]] sunrealtype t, N_Vector y, N_Vector ydot,
void* userData)
197 using ValueType =
typename SolutionVectorType::VectorValueType;
200 sunrealtype* yDotArray = N_VGetArrayPointer(ydot);
201 sunrealtype* yArray = N_VGetArrayPointer(y);
204 yDotArray !=
nullptr && yArray !=
nullptr && pdeExpre !=
nullptr,
205 "Failed to dereference pointers in sundails."
208 auto size =
static_cast<localIdx>(N_VGetLength(y));
211 if (std::holds_alternative<NeoN::GPUExecutor>(pdeExpre->
exec()))
215 NeoN::sundials::fieldToSunNVector(source, ydot);
229template<
typename Vector>
230void initNVector(
size_t size, std::shared_ptr<SUNContext> context, Vector& vec)
232 vec.initNVector(size, context);
241template<
typename Vector>
242const N_Vector& sunNVector(
const Vector& vec)
244 return vec.sunNVector();
253template<
typename Vector>
254N_Vector& sunNVector(Vector& vec)
256 return vec.sunNVector();
265template<
typename ValueType>
271 ~SKVectorSerial() =
default;
272 SKVectorSerial(
const SKVectorSerial& other)
273 : kvector_(other.kvector_), svector_(other.kvector_) {};
274 SKVectorSerial(SKVectorSerial&& other) noexcept
275 : kvector_(std::move(other.kvector_)), svector_(std::move(other.svector_)) {};
276 SKVectorSerial& operator=(
const SKVectorSerial& other) =
delete;
277 SKVectorSerial& operator=(SKVectorSerial&& other) =
delete;
280 using KVector = ::sundials::kokkos::Vector<Kokkos::Serial>;
281 void initNVector(
size_t size, std::shared_ptr<SUNContext> context)
283 kvector_ = KVector(size, *context);
286 const N_Vector& sunNVector()
const {
return svector_; };
287 N_Vector& sunNVector() {
return svector_; };
292 N_Vector svector_ {
nullptr};
300template<
typename ValueType>
301class SKVectorHostDefault
305 using KVector = ::sundials::kokkos::Vector<Kokkos::DefaultHostExecutionSpace>;
307 SKVectorHostDefault() =
default;
308 ~SKVectorHostDefault() =
default;
309 SKVectorHostDefault(
const SKVectorHostDefault& other)
310 : kvector_(other.kvector_), svector_(other.kvector_) {};
311 SKVectorHostDefault(SKVectorHostDefault&& other) noexcept
312 : kvector_(std::move(other.kvector_)), svector_(std::move(other.svector_)) {};
313 SKVectorHostDefault& operator=(
const SKVectorHostDefault& other) =
delete;
314 SKVectorHostDefault& operator=(SKVectorHostDefault&& other) =
delete;
316 void initNVector(
size_t size, std::shared_ptr<SUNContext> context)
318 kvector_ = KVector(size, *context);
321 const N_Vector& sunNVector()
const {
return svector_; };
322 N_Vector& sunNVector() {
return svector_; };
327 N_Vector svector_ {
nullptr};
335template<
typename ValueType>
340 using KVector = ::sundials::kokkos::Vector<Kokkos::DefaultExecutionSpace>;
342 SKVectorDefault() =
default;
343 ~SKVectorDefault() =
default;
344 SKVectorDefault(
const SKVectorDefault& other)
345 : kvector_(other.kvector_), svector_(other.kvector_) {};
346 SKVectorDefault(SKVectorDefault&& other) noexcept
347 : kvector_(std::move(other.kvector_)), svector_(std::move(other.svector_)) {};
348 SKVectorDefault& operator=(
const SKVectorDefault& other) =
delete;
349 SKVectorDefault& operator=(SKVectorDefault&& other) =
delete;
351 void initNVector(
size_t size, std::shared_ptr<SUNContext> context)
353 kvector_ = KVector(size, *context);
357 const N_Vector& sunNVector()
const {
return svector_; };
359 N_Vector& sunNVector() {
return svector_; };
364 N_Vector svector_ {
nullptr};
373template<
typename ValueType>
378 using SKVectorSerialV = SKVectorSerial<ValueType>;
379 using SKVectorHostDefaultV = SKVectorHostDefault<ValueType>;
380 using SKDefaultVectorV = SKVectorDefault<ValueType>;
381 using SKVectorVariant = std::variant<SKVectorSerialV, SKVectorHostDefaultV, SKDefaultVectorV>;
386 SKVector() { vector_.template emplace<SKVectorHostDefaultV>(); };
391 ~SKVector() =
default;
397 SKVector(
const SKVector&) =
default;
402 SKVector& operator=(
const SKVector&) =
delete;
408 SKVector(SKVector&&) noexcept = default;
413 SKVector& operator=(SKVector&&) noexcept = delete;
419 void setExecutor(const
NeoN::Executor& exec)
421 if (std::holds_alternative<NeoN::GPUExecutor>(exec))
423 vector_.template emplace<SKDefaultVectorV>();
426 if (std::holds_alternative<NeoN::CPUExecutor>(exec))
428 vector_.template emplace<SKVectorHostDefaultV>();
431 if (std::holds_alternative<NeoN::SerialExecutor>(exec))
433 vector_.template emplace<SKVectorSerialV>();
438 "Unsupported NeoN executor: "
439 << std::visit([](
const auto& e) {
return e.name(); }, exec) <<
"."
448 void initNVector(
size_t size, std::shared_ptr<SUNContext> context)
451 [size, &context](
auto& vec) { detail::initNVector(size, context, vec); }, vector_
459 const N_Vector& sunNVector()
const
462 [](
const auto& vec) ->
const N_Vector& {
return detail::sunNVector(vec); }, vector_
470 N_Vector& sunNVector()
472 return std::visit([](
auto& vec) -> N_Vector& {
return detail::sunNVector(vec); }, vector_);
479 const SKVectorVariant& variant()
const {
return vector_; }
485 SKVectorVariant& variant() {
return vector_; }
489 SKVectorVariant vector_;
A class to contain the data and executors for a field and define some basic operations.
ValueType * data()
Direct access to the underlying field data.
std::pair< localIdx, localIdx > range() const
Gets the range of the field.
const Executor & exec() const
Gets the executor associated with the field.
View< ValueType > view() &&=delete
Vector< ValueType > explicitOperation(localIdx nCells) const
const Executor & exec() const
#define NF_ERROR_EXIT(message)
Macro for printing an error message and aborting the program.
#define NF_ASSERT(condition, message)
Macro for asserting a condition and printing an error message if the condition is false.
SpatialOperator< scalar > source(fvcc::VolumeField< scalar > &coeff, fvcc::VolumeField< scalar > &phi)
void parallelFor(const Executor &exec, std::pair< localIdx, localIdx > range, Kernel kernel, std::string name="parallelFor")