13#include <sundials/sundials_nvector.h>
14#include <sundials/sundials_core.hpp>
15#include <nvector/nvector_serial.h>
16#include <nvector/nvector_kokkos.hpp>
17#include <arkode/arkode_arkstep.h>
18#include <arkode/arkode_erkstep.h>
24namespace NeoN::sundials
32inline auto SUN_CONTEXT_DELETER = [](SUNContext* ctx)
45inline auto SUN_ARK_DELETER = [](
char* ark)
49 void* arkodMem =
reinterpret_cast<void*
>(ark);
50 ARKodeFree(&arkodMem);
60inline ARKODE_ERKTableID stringToERKTable(
const std::string& key)
62 if (key ==
"Forward-Euler")
return ARKODE_FORWARD_EULER_1_1;
65 NF_ERROR_EXIT(
"Currently unsupported until field time step-stage indexing resolved.");
66 return ARKODE_HEUN_EULER_2_1_2;
68 if (key ==
"Midpoint")
70 NF_ERROR_EXIT(
"Currently unsupported until field time step-stage indexing resolved.");
71 return ARKODE_EXPLICIT_MIDPOINT_EULER_2_1_2;
74 "Unsupported Runge-Kutta time integration method selectied: " + key +
".\n"
75 +
"Supported methods are: Forward-Euler, Heun, Midpoint."
77 return ARKODE_ERK_NONE;
88template<
typename SKVectorType,
typename ValueType>
91 auto view = ::sundials::kokkos::GetVec<SKVectorType>(vector)->View();
92 auto fieldView = field.
view();
96 KOKKOS_LAMBDA(
const localIdx i) { view(i) = fieldView[i]; },
108template<
typename ValueType>
112 if (std::holds_alternative<NeoN::GPUExecutor>(field.
exec()))
114 fieldToSunNVectorImpl<::sundials::kokkos::Vector<Kokkos::DefaultExecutionSpace>>(
119 if (std::holds_alternative<NeoN::CPUExecutor>(field.
exec()))
121 fieldToSunNVectorImpl<::sundials::kokkos::Vector<Kokkos::DefaultHostExecutionSpace>>(
126 if (std::holds_alternative<NeoN::SerialExecutor>(field.
exec()))
128 fieldToSunNVectorImpl<::sundials::kokkos::Vector<Kokkos::Serial>>(field, vector);
142template<
typename SKVectorType,
typename ValueType>
145 auto view = ::sundials::kokkos::GetVec<SKVectorType>(vector)->View();
146 ValueType* fieldData = field.
data();
150 KOKKOS_LAMBDA(
const localIdx i) { fieldData[i] = view(i); },
151 "sunNVectorToVectorImpl"
161template<
typename ValueType>
164 if (std::holds_alternative<NeoN::GPUExecutor>(field.
exec()))
166 sunNVectorToVectorImpl<::sundials::kokkos::Vector<Kokkos::DefaultExecutionSpace>>(
171 if (std::holds_alternative<NeoN::CPUExecutor>(field.
exec()))
173 sunNVectorToVectorImpl<::sundials::kokkos::Vector<Kokkos::DefaultHostExecutionSpace>>(
178 if (std::holds_alternative<NeoN::SerialExecutor>(field.
exec()))
180 sunNVectorToVectorImpl<::sundials::kokkos::Vector<Kokkos::Serial>>(vector, field);
200template<
typename SolutionVectorType>
201int explicitRKSolve([[maybe_unused]] sunrealtype t, N_Vector y, N_Vector ydot,
void* userData)
204 using ValueType =
typename SolutionVectorType::VectorValueType;
207 sunrealtype* yDotArray = N_VGetArrayPointer(ydot);
208 sunrealtype* yArray = N_VGetArrayPointer(y);
211 yDotArray !=
nullptr && yArray !=
nullptr && pdeExpre !=
nullptr,
212 "Failed to dereference pointers in sundails."
215 auto size =
static_cast<localIdx>(N_VGetLength(y));
219 NeoN::sundials::fieldToSunNVector(source, ydot);
233template<
typename Vector>
234void initNVector(
size_t size, std::shared_ptr<SUNContext> context, Vector& vec)
236 vec.initNVector(size, context);
245template<
typename Vector>
246const N_Vector& sunNVector(
const Vector& vec)
248 return vec.sunNVector();
257template<
typename Vector>
258N_Vector& sunNVector(Vector& vec)
260 return vec.sunNVector();
269template<
typename ValueType>
275 ~SKVectorSerial() =
default;
276 SKVectorSerial(
const SKVectorSerial& other)
277 : kvector_(other.kvector_), svector_(other.kvector_) {};
278 SKVectorSerial(SKVectorSerial&& other) noexcept
279 : kvector_(std::move(other.kvector_)), svector_(std::move(other.svector_)) {};
280 SKVectorSerial& operator=(
const SKVectorSerial& other) =
delete;
281 SKVectorSerial& operator=(SKVectorSerial&& other) =
delete;
284 using KVector = ::sundials::kokkos::Vector<Kokkos::Serial>;
285 void initNVector(
size_t size, std::shared_ptr<SUNContext> context)
287 kvector_ = KVector(size, *context);
290 const N_Vector& sunNVector()
const {
return svector_; };
291 N_Vector& sunNVector() {
return svector_; };
296 N_Vector svector_ {
nullptr};
304template<
typename ValueType>
305class SKVectorHostDefault
309 using KVector = ::sundials::kokkos::Vector<Kokkos::DefaultHostExecutionSpace>;
311 SKVectorHostDefault() =
default;
312 ~SKVectorHostDefault() =
default;
313 SKVectorHostDefault(
const SKVectorHostDefault& other)
314 : kvector_(other.kvector_), svector_(other.kvector_) {};
315 SKVectorHostDefault(SKVectorHostDefault&& other) noexcept
316 : kvector_(std::move(other.kvector_)), svector_(std::move(other.svector_)) {};
317 SKVectorHostDefault& operator=(
const SKVectorHostDefault& other) =
delete;
318 SKVectorHostDefault& operator=(SKVectorHostDefault&& other) =
delete;
320 void initNVector(
size_t size, std::shared_ptr<SUNContext> context)
322 kvector_ = KVector(size, *context);
325 const N_Vector& sunNVector()
const {
return svector_; };
326 N_Vector& sunNVector() {
return svector_; };
331 N_Vector svector_ {
nullptr};
339template<
typename ValueType>
344 using KVector = ::sundials::kokkos::Vector<Kokkos::DefaultExecutionSpace>;
346 SKVectorDefault() =
default;
347 ~SKVectorDefault() =
default;
348 SKVectorDefault(
const SKVectorDefault& other)
349 : kvector_(other.kvector_), svector_(other.kvector_) {};
350 SKVectorDefault(SKVectorDefault&& other) noexcept
351 : kvector_(std::move(other.kvector_)), svector_(std::move(other.svector_)) {};
352 SKVectorDefault& operator=(
const SKVectorDefault& other) =
delete;
353 SKVectorDefault& operator=(SKVectorDefault&& other) =
delete;
355 void initNVector(
size_t size, std::shared_ptr<SUNContext> context)
357 kvector_ = KVector(size, *context);
361 const N_Vector& sunNVector()
const {
return svector_; };
363 N_Vector& sunNVector() {
return svector_; };
368 N_Vector svector_ {
nullptr};
377template<
typename ValueType>
382 using SKVectorSerialV = SKVectorSerial<ValueType>;
383 using SKVectorHostDefaultV = SKVectorHostDefault<ValueType>;
384 using SKDefaultVectorV = SKVectorDefault<ValueType>;
385 using SKVectorVariant = std::variant<SKVectorSerialV, SKVectorHostDefaultV, SKDefaultVectorV>;
390 SKVector() { vector_.template emplace<SKVectorHostDefaultV>(); };
395 ~SKVector() =
default;
401 SKVector(
const SKVector&) =
default;
406 SKVector& operator=(
const SKVector&) =
delete;
412 SKVector(SKVector&&) noexcept = default;
417 SKVector& operator=(SKVector&&) noexcept = delete;
423 void setExecutor(const
NeoN::Executor& exec)
425 if (std::holds_alternative<NeoN::GPUExecutor>(exec))
427 vector_.template emplace<SKDefaultVectorV>();
430 if (std::holds_alternative<NeoN::CPUExecutor>(exec))
432 vector_.template emplace<SKVectorHostDefaultV>();
435 if (std::holds_alternative<NeoN::SerialExecutor>(exec))
437 vector_.template emplace<SKVectorSerialV>();
442 "Unsupported NeoN executor: "
443 << std::visit([](
const auto& e) {
return e.name(); }, exec) <<
"."
452 void initNVector(
size_t size, std::shared_ptr<SUNContext> context)
455 [size, &context](
auto& vec) { detail::initNVector(size, context, vec); }, vector_
463 const N_Vector& sunNVector()
const
466 [](
const auto& vec) ->
const N_Vector& {
return detail::sunNVector(vec); }, vector_
474 N_Vector& sunNVector()
476 return std::visit([](
auto& vec) -> N_Vector& {
return detail::sunNVector(vec); }, vector_);
483 const SKVectorVariant& variant()
const {
return vector_; }
489 SKVectorVariant& variant() {
return vector_; }
493 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 fence(const Executor &exec)
void parallelFor(const ExecutorType &exec, std::pair< localIdx, localIdx > range, Kernel kernel, std::string name)