commit 95b21ade145802442b2b9b99d73fd479641a5a85
parent 5191a48df06ea47104ca67db19fa82058d09c20a
Author: d.levin256@gmail.com <d.levin256@gmail.com>
Date: Mon, 18 Nov 2024 14:40:02 +0100
Add multidimensional DFT to C API
Diffstat:
4 files changed, 683 insertions(+), 377 deletions(-)
diff --git a/include/kfr/capi.h b/include/kfr/capi.h
@@ -185,6 +185,31 @@ typedef enum KFR_DFT_PACK_FORMAT
KFR_API_SPEC KFR_DFT_PLAN_F32* kfr_dft_create_plan_f32(size_t size);
/**
+ * @brief Create a 2D complex DFT plan (Single precision).
+ * @param size1 Size of the first dimension.
+ * @param size2 Size of the second dimension.
+ * @return Pointer to the created 2D DFT plan. Use `kfr_dft_delete_plan_f**` to free.
+ */
+KFR_API_SPEC KFR_DFT_PLAN_F32* kfr_dft_create_2d_plan_f32(size_t size1, size_t size2);
+
+/**
+ * @brief Create a 3D complex DFT plan (Single precision).
+ * @param size1 Size of the first dimension.
+ * @param size2 Size of the second dimension.
+ * @param size3 Size of the third dimension.
+ * @return Pointer to the created 3D DFT plan. Use `kfr_dft_delete_plan_f**` to free.
+ */
+KFR_API_SPEC KFR_DFT_PLAN_F32* kfr_dft_create_3d_plan_f32(size_t size1, size_t size2, size_t size3);
+
+/**
+ * @brief Create an N-dimensional complex DFT plan (Single precision).
+ * @param dims Number of dimensions.
+ * @param shape Array of sizes for each dimension.
+ * @return Pointer to the created N-dimensional DFT plan. Use `kfr_dft_delete_plan_f**` to free.
+ */
+KFR_API_SPEC KFR_DFT_PLAN_F32* kfr_dft_create_md_plan_f32(size_t dims, const unsigned* shape);
+
+/**
* @brief Create a complex DFT plan (Double precision).
* @param size Size of the DFT.
* @return Pointer to the created DFT plan. Use `kfr_dft_delete_plan_f**` to free.
@@ -192,6 +217,31 @@ KFR_API_SPEC KFR_DFT_PLAN_F32* kfr_dft_create_plan_f32(size_t size);
KFR_API_SPEC KFR_DFT_PLAN_F64* kfr_dft_create_plan_f64(size_t size);
/**
+ * @brief Create a 2D complex DFT plan (Double precision).
+ * @param size1 Size of the first dimension.
+ * @param size2 Size of the second dimension.
+ * @return Pointer to the created 2D DFT plan. Use `kfr_dft_delete_plan_f**` to free.
+ */
+KFR_API_SPEC KFR_DFT_PLAN_F64* kfr_dft_create_2d_plan_f64(size_t size1, size_t size2);
+
+/**
+ * @brief Create a 3D complex DFT plan (Double precision).
+ * @param size1 Size of the first dimension.
+ * @param size2 Size of the second dimension.
+ * @param size3 Size of the third dimension.
+ * @return Pointer to the created 3D DFT plan. Use `kfr_dft_delete_plan_f**` to free.
+ */
+KFR_API_SPEC KFR_DFT_PLAN_F64* kfr_dft_create_3d_plan_f64(size_t size1, size_t size2, size_t size3);
+
+/**
+ * @brief Create an N-dimensional complex DFT plan (Double precision).
+ * @param dims Number of dimensions.
+ * @param shape Array of sizes for each dimension.
+ * @return Pointer to the created N-dimensional DFT plan. Use `kfr_dft_delete_plan_f**` to free.
+ */
+KFR_API_SPEC KFR_DFT_PLAN_F64* kfr_dft_create_md_plan_f64(size_t dims, const unsigned* shape);
+
+/**
* @brief Dump details of the DFT plan to stdout for inspection.
* @param plan Pointer to the DFT plan.
*/
@@ -308,6 +358,13 @@ KFR_API_SPEC void kfr_dft_delete_plan_f64(KFR_DFT_PLAN_F64* plan);
KFR_API_SPEC KFR_DFT_REAL_PLAN_F32* kfr_dft_real_create_plan_f32(size_t size,
KFR_DFT_PACK_FORMAT pack_format);
+KFR_API_SPEC KFR_DFT_REAL_PLAN_F32* kfr_dft_real_create_2d_plan_f32(size_t size1, size_t size2,
+ bool real_out_is_enough);
+KFR_API_SPEC KFR_DFT_REAL_PLAN_F32* kfr_dft_real_create_3d_plan_f32(size_t size1, size_t size2, size_t size3,
+ bool real_out_is_enough);
+KFR_API_SPEC KFR_DFT_REAL_PLAN_F32* kfr_dft_real_create_md_plan_f32(size_t dims, const unsigned* shape,
+ bool real_out_is_enough);
+
/**
* @brief Create a real DFT plan (Double precision).
* @param size Size of the real DFT. Must be even.
@@ -317,6 +374,13 @@ KFR_API_SPEC KFR_DFT_REAL_PLAN_F32* kfr_dft_real_create_plan_f32(size_t size,
KFR_API_SPEC KFR_DFT_REAL_PLAN_F64* kfr_dft_real_create_plan_f64(size_t size,
KFR_DFT_PACK_FORMAT pack_format);
+KFR_API_SPEC KFR_DFT_REAL_PLAN_F64* kfr_dft_real_create_2d_plan_f64(size_t size1, size_t size2,
+ bool real_out_is_enough);
+KFR_API_SPEC KFR_DFT_REAL_PLAN_F64* kfr_dft_real_create_3d_plan_f64(size_t size1, size_t size2, size_t size3,
+ bool real_out_is_enough);
+KFR_API_SPEC KFR_DFT_REAL_PLAN_F64* kfr_dft_real_create_md_plan_f64(size_t dims, const unsigned* shape,
+ bool real_out_is_enough);
+
/**
* @brief Dump details of the real DFT plan to stdout for inspection.
* @param plan Pointer to the DFT plan.
diff --git a/include/kfr/dft/fft.hpp b/include/kfr/dft/fft.hpp
@@ -354,7 +354,7 @@ struct dft_plan_real : dft_plan<T>
};
/// @brief Multidimensional DFT
-template <typename T, index_t Dims>
+template <typename T, index_t Dims = dynamic_shape>
struct dft_plan_md
{
shape<Dims> size;
@@ -490,7 +490,7 @@ private:
};
/// @brief Multidimensional DFT
-template <typename T, index_t Dims>
+template <typename T, index_t Dims = dynamic_shape>
struct dft_plan_md_real
{
shape<Dims> size;
diff --git a/src/capi/capi.cpp b/src/capi/capi.cpp
@@ -23,6 +23,7 @@
disclosing the source code of your own applications.
See https://www.kfrlib.com for details.
*/
+#include <variant>
#define KFR_NO_C_COMPLEX_TYPES 1
#include <kfr/capi.h>
@@ -81,399 +82,562 @@ static void try_fn(Fn&& fn)
}
}
-extern "C"
+template <typename T>
+class var_dft_plan
+{
+public:
+ virtual ~var_dft_plan() {}
+ virtual void dump() = 0;
+ virtual size_t size() = 0;
+ virtual size_t temp_size() = 0;
+ virtual void execute(T*, const T*, u8*) = 0;
+ virtual void execute_inverse(T*, const T*, u8*) = 0;
+};
+
+static shape<dynamic_shape> init_shape(size_t dims, const unsigned* shape_)
{
- KFR_API_SPEC const char* kfr_version_string()
+ shape<dynamic_shape> sh(dims);
+ for (size_t i = 0; i < dims; ++i)
{
- return "KFR " KFR_VERSION_STRING KFR_DEBUG_STR " " KFR_ENABLED_ARCHS_LIST " " CMT_ARCH_BITNESS_NAME
- " (" CMT_COMPILER_FULL_NAME "/" CMT_OS_NAME ")" KFR_BUILD_DETAILS_1 KFR_BUILD_DETAILS_2;
+ sh[i] = shape_[i];
}
- KFR_API_SPEC uint32_t kfr_version() { return KFR_VERSION; }
- KFR_API_SPEC const char* kfr_enabled_archs() { return KFR_ENABLED_ARCHS_LIST; }
- KFR_API_SPEC int kfr_current_arch() { return static_cast<int>(get_cpu()); }
+ return sh;
+}
- KFR_API_SPEC const char* kfr_last_error() { return error.data(); }
+template <typename T, size_t Dims>
+struct var_dft_plan_select
+{
+ using complex = dft_plan_md<T, dynamic_shape>;
+ using real = dft_plan_md_real<T, dynamic_shape>;
+
+ static size_t size(const complex& plan) { return plan.size.product(); }
+ static size_t size(const real& plan) { return plan.size.product(); }
+};
+template <typename T>
+struct var_dft_plan_select<T, 1>
+{
+ using complex = dft_plan<T>;
+ using real = dft_plan_real<T>;
- KFR_API_SPEC void* kfr_allocate(size_t size)
- {
- return details::aligned_malloc(size, KFR_DEFAULT_ALIGNMENT);
- }
- KFR_API_SPEC void* kfr_allocate_aligned(size_t size, size_t alignment)
- {
- return details::aligned_malloc(size, alignment);
- }
- KFR_API_SPEC void kfr_deallocate(void* ptr) { return details::aligned_free(ptr); }
- KFR_API_SPEC size_t kfr_allocated_size(void* ptr) { return details::aligned_size(ptr); }
+ static size_t size(const complex& plan) { return plan.size; }
+ static size_t size(const real& plan) { return plan.size; }
+};
- KFR_API_SPEC void* kfr_add_ref(void* ptr)
+template <typename T, size_t Dims>
+class var_dft_plan_impl final : public var_dft_plan<T>
+{
+public:
+ template <typename... Args>
+ CMT_ALWAYS_INLINE var_dft_plan_impl(Args&&... args) : plan(std::forward<Args>(args)...)
{
- details::aligned_add_ref(ptr);
- return ptr;
}
- KFR_API_SPEC void kfr_release(void* ptr) { details::aligned_release(ptr); }
-
- KFR_API_SPEC void* kfr_reallocate(void* ptr, size_t new_size)
+ typename var_dft_plan_select<T, Dims>::complex plan;
+ void dump() { plan.dump(); }
+ size_t size() { return var_dft_plan_select<T, Dims>::size(plan); }
+ size_t temp_size() { return plan.temp_size; }
+ void execute(T* out, const T* in, u8* temp)
{
- return details::aligned_reallocate(ptr, new_size, KFR_DEFAULT_ALIGNMENT);
+ plan.execute(reinterpret_cast<complex<T>*>(out), reinterpret_cast<const complex<T>*>(in), temp,
+ cfalse);
}
- KFR_API_SPEC void* kfr_reallocate_aligned(void* ptr, size_t new_size, size_t alignment)
+ void execute_inverse(T* out, const T* in, u8* temp)
{
- return details::aligned_reallocate(ptr, new_size, alignment);
+ plan.execute(reinterpret_cast<complex<T>*>(out), reinterpret_cast<const complex<T>*>(in), temp,
+ ctrue);
}
+};
- KFR_API_SPEC KFR_DFT_PLAN_F32* kfr_dft_create_plan_f32(size_t size)
- {
- return try_fn([&]() { return reinterpret_cast<KFR_DFT_PLAN_F32*>(new kfr::dft_plan<float>(size)); },
- nullptr);
- }
- KFR_API_SPEC KFR_DFT_PLAN_F64* kfr_dft_create_plan_f64(size_t size)
+template <typename T, size_t Dims>
+class var_dft_plan_real_impl final : public var_dft_plan<T>
+{
+public:
+ template <typename... Args>
+ CMT_ALWAYS_INLINE var_dft_plan_real_impl(Args&&... args) : plan(std::forward<Args>(args)...)
{
- return try_fn([&]() { return reinterpret_cast<KFR_DFT_PLAN_F64*>(new kfr::dft_plan<double>(size)); },
- nullptr);
}
-
- KFR_API_SPEC void kfr_dft_dump_f32(KFR_DFT_PLAN_F32* plan)
+ typename var_dft_plan_select<T, Dims>::real plan;
+ void dump() { plan.dump(); }
+ size_t size() { return var_dft_plan_select<T, Dims>::size(plan); }
+ size_t temp_size() { return plan.temp_size; }
+ void execute(T* out, const T* in, u8* temp)
{
- try_fn([&] { reinterpret_cast<kfr::dft_plan<float>*>(plan)->dump(); });
+ plan.execute(reinterpret_cast<complex<T>*>(out), reinterpret_cast<const T*>(in), temp, cfalse);
}
- KFR_API_SPEC void kfr_dft_dump_f64(KFR_DFT_PLAN_F64* plan)
+ void execute_inverse(T* out, const T* in, u8* temp)
{
- try_fn([&] { reinterpret_cast<kfr::dft_plan<double>*>(plan)->dump(); });
+ plan.execute(reinterpret_cast<T*>(out), reinterpret_cast<const complex<T>*>(in), temp, ctrue);
}
+};
- KFR_API_SPEC size_t kfr_dft_get_size_f32(KFR_DFT_PLAN_F32* plan)
- {
- return try_fn([&]() { return reinterpret_cast<kfr::dft_plan<float>*>(plan)->size; }, 0);
- }
- KFR_API_SPEC size_t kfr_dft_get_size_f64(KFR_DFT_PLAN_F64* plan)
- {
- return try_fn([&]() { return reinterpret_cast<kfr::dft_plan<double>*>(plan)->size; }, 0);
- }
+extern "C"
+{
+KFR_API_SPEC const char* kfr_version_string()
+{
+ return "KFR " KFR_VERSION_STRING KFR_DEBUG_STR " " KFR_ENABLED_ARCHS_LIST " " CMT_ARCH_BITNESS_NAME
+ " (" CMT_COMPILER_FULL_NAME "/" CMT_OS_NAME ")" KFR_BUILD_DETAILS_1 KFR_BUILD_DETAILS_2;
+}
+KFR_API_SPEC uint32_t kfr_version() { return KFR_VERSION; }
+KFR_API_SPEC const char* kfr_enabled_archs() { return KFR_ENABLED_ARCHS_LIST; }
+KFR_API_SPEC int kfr_current_arch() { return static_cast<int>(get_cpu()); }
- KFR_API_SPEC size_t kfr_dft_get_temp_size_f32(KFR_DFT_PLAN_F32* plan)
- {
- return try_fn([&]() { return reinterpret_cast<kfr::dft_plan<float>*>(plan)->temp_size; }, 0);
- }
- KFR_API_SPEC size_t kfr_dft_get_temp_size_f64(KFR_DFT_PLAN_F64* plan)
- {
- return try_fn([&]() { return reinterpret_cast<kfr::dft_plan<double>*>(plan)->temp_size; }, 0);
- }
+KFR_API_SPEC const char* kfr_last_error() { return error.data(); }
- KFR_API_SPEC void kfr_dft_execute_f32(KFR_DFT_PLAN_F32* plan, kfr_c32* out, const kfr_c32* in,
- uint8_t* temp)
- {
- try_fn(
- [&]()
- {
- reinterpret_cast<kfr::dft_plan<float>*>(plan)->execute(
- reinterpret_cast<kfr::complex<float>*>(out),
- reinterpret_cast<const kfr::complex<float>*>(in), temp, kfr::cfalse);
- });
- }
- KFR_API_SPEC void kfr_dft_execute_f64(KFR_DFT_PLAN_F64* plan, kfr_c64* out, const kfr_c64* in,
- uint8_t* temp)
- {
- try_fn(
- [&]()
- {
- reinterpret_cast<kfr::dft_plan<double>*>(plan)->execute(
- reinterpret_cast<kfr::complex<double>*>(out),
- reinterpret_cast<const kfr::complex<double>*>(in), temp, kfr::cfalse);
- });
- }
- KFR_API_SPEC void kfr_dft_execute_inverse_f32(KFR_DFT_PLAN_F32* plan, kfr_c32* out, const kfr_c32* in,
- uint8_t* temp)
- {
- try_fn(
- [&]()
- {
- reinterpret_cast<kfr::dft_plan<float>*>(plan)->execute(
- reinterpret_cast<kfr::complex<float>*>(out),
- reinterpret_cast<const kfr::complex<float>*>(in), temp, kfr::ctrue);
- });
- }
- KFR_API_SPEC void kfr_dft_execute_inverse_f64(KFR_DFT_PLAN_F64* plan, kfr_c64* out, const kfr_c64* in,
- uint8_t* temp)
- {
- try_fn(
- [&]()
- {
- reinterpret_cast<kfr::dft_plan<double>*>(plan)->execute(
- reinterpret_cast<kfr::complex<double>*>(out),
- reinterpret_cast<const kfr::complex<double>*>(in), temp, kfr::ctrue);
- });
- }
+KFR_API_SPEC void* kfr_allocate(size_t size) { return details::aligned_malloc(size, KFR_DEFAULT_ALIGNMENT); }
+KFR_API_SPEC void* kfr_allocate_aligned(size_t size, size_t alignment)
+{
+ return details::aligned_malloc(size, alignment);
+}
+KFR_API_SPEC void kfr_deallocate(void* ptr) { return details::aligned_free(ptr); }
+KFR_API_SPEC size_t kfr_allocated_size(void* ptr) { return details::aligned_size(ptr); }
- KFR_API_SPEC void kfr_dft_delete_plan_f32(KFR_DFT_PLAN_F32* plan)
- {
- try_fn([&]() { delete reinterpret_cast<kfr::dft_plan<float>*>(plan); });
- }
- KFR_API_SPEC void kfr_dft_delete_plan_f64(KFR_DFT_PLAN_F64* plan)
- {
- try_fn([&]() { delete reinterpret_cast<kfr::dft_plan<double>*>(plan); });
- }
+KFR_API_SPEC void* kfr_add_ref(void* ptr)
+{
+ details::aligned_add_ref(ptr);
+ return ptr;
+}
+KFR_API_SPEC void kfr_release(void* ptr) { details::aligned_release(ptr); }
- // Real DFT plans
+KFR_API_SPEC void* kfr_reallocate(void* ptr, size_t new_size)
+{
+ return details::aligned_reallocate(ptr, new_size, KFR_DEFAULT_ALIGNMENT);
+}
+KFR_API_SPEC void* kfr_reallocate_aligned(void* ptr, size_t new_size, size_t alignment)
+{
+ return details::aligned_reallocate(ptr, new_size, alignment);
+}
- KFR_API_SPEC KFR_DFT_REAL_PLAN_F32* kfr_dft_real_create_plan_f32(size_t size,
- KFR_DFT_PACK_FORMAT pack_format)
- {
- return try_fn(
- [&]()
- {
- return reinterpret_cast<KFR_DFT_REAL_PLAN_F32*>(
- new kfr::dft_plan_real<float>(size, static_cast<dft_pack_format>(pack_format)));
- },
- nullptr);
- }
- KFR_API_SPEC KFR_DFT_REAL_PLAN_F64* kfr_dft_real_create_plan_f64(size_t size,
- KFR_DFT_PACK_FORMAT pack_format)
- {
- return try_fn(
- [&]()
- {
- return reinterpret_cast<KFR_DFT_REAL_PLAN_F64*>(
- new kfr::dft_plan_real<double>(size, static_cast<dft_pack_format>(pack_format)));
- },
- nullptr);
- }
+KFR_API_SPEC KFR_DFT_PLAN_F32* kfr_dft_create_plan_f32(size_t size)
+{
+ return try_fn([&]()
+ { return reinterpret_cast<KFR_DFT_PLAN_F32*>(new var_dft_plan_impl<float, 1>(size)); },
+ nullptr);
+}
- KFR_API_SPEC void kfr_dft_real_dump_f32(KFR_DFT_REAL_PLAN_F32* plan)
- {
- try_fn([&]() { reinterpret_cast<kfr::dft_plan_real<float>*>(plan)->dump(); });
- }
- KFR_API_SPEC void kfr_dft_real_dump_f64(KFR_DFT_REAL_PLAN_F64* plan)
- {
- try_fn([&]() { reinterpret_cast<kfr::dft_plan_real<double>*>(plan)->dump(); });
- }
+KFR_API_SPEC KFR_DFT_PLAN_F32* kfr_dft_create_2d_plan_f32(size_t size1, size_t size2)
+{
+ return try_fn(
+ [&]() {
+ return reinterpret_cast<KFR_DFT_PLAN_F32*>(
+ new var_dft_plan_impl<float, 2>(shape{ size1, size2 }));
+ },
+ nullptr);
+}
+KFR_API_SPEC KFR_DFT_PLAN_F32* kfr_dft_create_3d_plan_f32(size_t size1, size_t size2, size_t size3)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_DFT_PLAN_F32*>(
+ new var_dft_plan_impl<float, 3>(shape{ size1, size2, size3 }));
+ },
+ nullptr);
+}
+KFR_API_SPEC KFR_DFT_PLAN_F32* kfr_dft_create_md_plan_f32(size_t dims, const unsigned* shape)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_DFT_PLAN_F32*>(
+ new var_dft_plan_impl<float, dynamic_shape>(init_shape(dims, shape)));
+ },
+ nullptr);
+}
+KFR_API_SPEC KFR_DFT_PLAN_F64* kfr_dft_create_plan_f64(size_t size)
+{
+ return try_fn([&]()
+ { return reinterpret_cast<KFR_DFT_PLAN_F64*>(new var_dft_plan_impl<double, 1>(size)); },
+ nullptr);
+}
+KFR_API_SPEC KFR_DFT_PLAN_F64* kfr_dft_create_2d_plan_f64(size_t size1, size_t size2)
+{
+ return try_fn(
+ [&]() {
+ return reinterpret_cast<KFR_DFT_PLAN_F64*>(
+ new var_dft_plan_impl<double, 2>(shape{ size1, size2 }));
+ },
+ nullptr);
+}
+KFR_API_SPEC KFR_DFT_PLAN_F64* kfr_dft_create_3d_plan_f64(size_t size1, size_t size2, size_t size3)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_DFT_PLAN_F64*>(
+ new var_dft_plan_impl<double, 3>(shape{ size1, size2, size3 }));
+ },
+ nullptr);
+}
+KFR_API_SPEC KFR_DFT_PLAN_F64* kfr_dft_create_md_plan_f64(size_t dims, const unsigned* shape)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_DFT_PLAN_F64*>(
+ new var_dft_plan_impl<double, dynamic_shape>(init_shape(dims, shape)));
+ },
+ nullptr);
+}
- KFR_API_SPEC size_t kfr_dft_real_get_size_f32(KFR_DFT_REAL_PLAN_F32* plan)
- {
- return try_fn([&]() { return reinterpret_cast<kfr::dft_plan<float>*>(plan)->size; }, 0);
- }
- KFR_API_SPEC size_t kfr_dft_real_get_size_f64(KFR_DFT_REAL_PLAN_F64* plan)
- {
- return try_fn([&]() { return reinterpret_cast<kfr::dft_plan<double>*>(plan)->size; }, 0);
- }
+KFR_API_SPEC void kfr_dft_dump_f32(KFR_DFT_PLAN_F32* plan)
+{
+ try_fn([&] { reinterpret_cast<var_dft_plan<float>*>(plan)->dump(); });
+}
+KFR_API_SPEC void kfr_dft_dump_f64(KFR_DFT_PLAN_F64* plan)
+{
+ try_fn([&] { reinterpret_cast<var_dft_plan<double>*>(plan)->dump(); });
+}
- KFR_API_SPEC size_t kfr_dft_real_get_temp_size_f32(KFR_DFT_REAL_PLAN_F32* plan)
- {
- return try_fn([&]() { return reinterpret_cast<kfr::dft_plan<float>*>(plan)->temp_size; }, 0);
- }
- KFR_API_SPEC size_t kfr_dft_real_get_temp_size_f64(KFR_DFT_REAL_PLAN_F64* plan)
- {
- return try_fn([&]() { return reinterpret_cast<kfr::dft_plan<double>*>(plan)->temp_size; }, 0);
- }
+KFR_API_SPEC size_t kfr_dft_get_size_f32(KFR_DFT_PLAN_F32* plan)
+{
+ return try_fn([&]() { return reinterpret_cast<var_dft_plan<float>*>(plan)->size(); }, 0);
+}
+KFR_API_SPEC size_t kfr_dft_get_size_f64(KFR_DFT_PLAN_F64* plan)
+{
+ return try_fn([&]() { return reinterpret_cast<var_dft_plan<double>*>(plan)->size(); }, 0);
+}
- KFR_API_SPEC void kfr_dft_real_execute_f32(KFR_DFT_REAL_PLAN_F32* plan, kfr_c32* out, const float* in,
- uint8_t* temp)
- {
- try_fn(
- [&]()
- {
- reinterpret_cast<kfr::dft_plan_real<float>*>(plan)->execute(
- reinterpret_cast<kfr::complex<float>*>(out), in, temp);
- });
- }
- KFR_API_SPEC void kfr_dft_real_execute_f64(KFR_DFT_REAL_PLAN_F64* plan, kfr_c64* out, const double* in,
- uint8_t* temp)
- {
- try_fn(
- [&]()
- {
- reinterpret_cast<kfr::dft_plan_real<double>*>(plan)->execute(
- reinterpret_cast<kfr::complex<double>*>(out), in, temp);
- });
- }
- KFR_API_SPEC void kfr_dft_real_execute_inverse_f32(KFR_DFT_REAL_PLAN_F32* plan, float* out,
- const kfr_c32* in, uint8_t* temp)
- {
- try_fn(
- [&]()
- {
- reinterpret_cast<kfr::dft_plan_real<float>*>(plan)->execute(
- out, reinterpret_cast<const kfr::complex<float>*>(in), temp);
- });
- }
- KFR_API_SPEC void kfr_dft_real_execute_inverse_f64(KFR_DFT_REAL_PLAN_F64* plan, double* out,
- const kfr_c64* in, uint8_t* temp)
- {
- try_fn(
- [&]()
- {
- reinterpret_cast<kfr::dft_plan_real<double>*>(plan)->execute(
- out, reinterpret_cast<const kfr::complex<double>*>(in), temp);
- });
- }
+KFR_API_SPEC size_t kfr_dft_get_temp_size_f32(KFR_DFT_PLAN_F32* plan)
+{
+ return try_fn([&]() { return reinterpret_cast<var_dft_plan<float>*>(plan)->temp_size(); }, 0);
+}
+KFR_API_SPEC size_t kfr_dft_get_temp_size_f64(KFR_DFT_PLAN_F64* plan)
+{
+ return try_fn([&]() { return reinterpret_cast<var_dft_plan<double>*>(plan)->temp_size(); }, 0);
+}
- KFR_API_SPEC void kfr_dft_real_delete_plan_f32(KFR_DFT_REAL_PLAN_F32* plan)
- {
- try_fn([&]() { delete reinterpret_cast<kfr::dft_plan_real<float>*>(plan); });
- }
- KFR_API_SPEC void kfr_dft_real_delete_plan_f64(KFR_DFT_REAL_PLAN_F64* plan)
- {
- try_fn([&]() { delete reinterpret_cast<kfr::dft_plan_real<double>*>(plan); });
- }
+KFR_API_SPEC void kfr_dft_execute_f32(KFR_DFT_PLAN_F32* plan, kfr_c32* out, const kfr_c32* in, uint8_t* temp)
+{
+ try_fn(
+ [&]()
+ {
+ reinterpret_cast<var_dft_plan<float>*>(plan)->execute(reinterpret_cast<float*>(out),
+ reinterpret_cast<const float*>(in), temp);
+ });
+}
+KFR_API_SPEC void kfr_dft_execute_f64(KFR_DFT_PLAN_F64* plan, kfr_c64* out, const kfr_c64* in, uint8_t* temp)
+{
+ try_fn(
+ [&]()
+ {
+ reinterpret_cast<var_dft_plan<double>*>(plan)->execute(reinterpret_cast<double*>(out),
+ reinterpret_cast<const double*>(in), temp);
+ });
+}
+KFR_API_SPEC void kfr_dft_execute_inverse_f32(KFR_DFT_PLAN_F32* plan, kfr_c32* out, const kfr_c32* in,
+ uint8_t* temp)
+{
+ try_fn(
+ [&]()
+ {
+ reinterpret_cast<var_dft_plan<float>*>(plan)->execute_inverse(
+ reinterpret_cast<float*>(out), reinterpret_cast<const float*>(in), temp);
+ });
+}
+KFR_API_SPEC void kfr_dft_execute_inverse_f64(KFR_DFT_PLAN_F64* plan, kfr_c64* out, const kfr_c64* in,
+ uint8_t* temp)
+{
+ try_fn(
+ [&]()
+ {
+ reinterpret_cast<var_dft_plan<double>*>(plan)->execute_inverse(
+ reinterpret_cast<double*>(out), reinterpret_cast<const double*>(in), temp);
+ });
+}
- // Discrete Cosine Transform
+KFR_API_SPEC void kfr_dft_delete_plan_f32(KFR_DFT_PLAN_F32* plan)
+{
+ try_fn([&]() { delete reinterpret_cast<var_dft_plan<float>*>(plan); });
+}
+KFR_API_SPEC void kfr_dft_delete_plan_f64(KFR_DFT_PLAN_F64* plan)
+{
+ try_fn([&]() { delete reinterpret_cast<var_dft_plan<double>*>(plan); });
+}
- KFR_API_SPEC KFR_DCT_PLAN_F32* kfr_dct_create_plan_f32(size_t size)
- {
- return try_fn([&]() { return reinterpret_cast<KFR_DCT_PLAN_F32*>(new kfr::dct_plan<float>(size)); },
- nullptr);
- }
- KFR_API_SPEC KFR_DCT_PLAN_F64* kfr_dct_create_plan_f64(size_t size)
- {
- return try_fn([&]() { return reinterpret_cast<KFR_DCT_PLAN_F64*>(new kfr::dct_plan<double>(size)); },
- nullptr);
- }
+// Real DFT plans
- KFR_API_SPEC void kfr_dct_dump_f32(KFR_DCT_PLAN_F32* plan)
- {
- try_fn([&]() { reinterpret_cast<kfr::dct_plan<float>*>(plan)->dump(); });
- }
- KFR_API_SPEC void kfr_dct_dump_f64(KFR_DCT_PLAN_F64* plan)
- {
- try_fn([&]() { reinterpret_cast<kfr::dct_plan<double>*>(plan)->dump(); });
- }
+KFR_API_SPEC KFR_DFT_REAL_PLAN_F32* kfr_dft_real_create_plan_f32(size_t size, KFR_DFT_PACK_FORMAT pack_format)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_DFT_REAL_PLAN_F32*>(
+ new var_dft_plan_real_impl<float, 1>(size, static_cast<dft_pack_format>(pack_format)));
+ },
+ nullptr);
+}
- KFR_API_SPEC size_t kfr_dct_get_size_f32(KFR_DCT_PLAN_F32* plan)
- {
- return try_fn([&]() { return reinterpret_cast<kfr::dft_plan<float>*>(plan)->size; }, 0);
- }
- KFR_API_SPEC size_t kfr_dct_get_size_f64(KFR_DCT_PLAN_F64* plan)
- {
- return try_fn([&]() { return reinterpret_cast<kfr::dft_plan<double>*>(plan)->size; }, 0);
- }
+KFR_API_SPEC KFR_DFT_REAL_PLAN_F32* kfr_dft_real_create_2d_plan_f32(size_t size1, size_t size2,
+ bool real_out_is_enough)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_DFT_REAL_PLAN_F32*>(
+ new var_dft_plan_real_impl<float, 2>(shape{ size1, size2 }, real_out_is_enough));
+ },
+ nullptr);
+}
+KFR_API_SPEC KFR_DFT_REAL_PLAN_F32* kfr_dft_real_create_3d_plan_f32(size_t size1, size_t size2, size_t size3,
+ bool real_out_is_enough)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_DFT_REAL_PLAN_F32*>(
+ new var_dft_plan_real_impl<float, 3>(shape{ size1, size2, size3 }, real_out_is_enough));
+ },
+ nullptr);
+}
+KFR_API_SPEC KFR_DFT_REAL_PLAN_F32* kfr_dft_real_create_md_plan_f32(size_t dims, const unsigned* shape,
+ bool real_out_is_enough)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_DFT_REAL_PLAN_F32*>(
+ new var_dft_plan_real_impl<float, dynamic_shape>(init_shape(dims, shape), real_out_is_enough));
+ },
+ nullptr);
+}
- KFR_API_SPEC size_t kfr_dct_get_temp_size_f32(KFR_DCT_PLAN_F32* plan)
- {
- return try_fn([&]() { return reinterpret_cast<kfr::dft_plan<float>*>(plan)->temp_size; }, 0);
- }
- KFR_API_SPEC size_t kfr_dct_get_temp_size_f64(KFR_DCT_PLAN_F64* plan)
- {
- return try_fn([&]() { return reinterpret_cast<kfr::dft_plan<double>*>(plan)->temp_size; }, 0);
- }
+KFR_API_SPEC KFR_DFT_REAL_PLAN_F64* kfr_dft_real_create_plan_f64(size_t size, KFR_DFT_PACK_FORMAT pack_format)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_DFT_REAL_PLAN_F64*>(
+ new var_dft_plan_real_impl<double, 1>(size, static_cast<dft_pack_format>(pack_format)));
+ },
+ nullptr);
+}
- KFR_API_SPEC void kfr_dct_execute_f32(KFR_DCT_PLAN_F32* plan, float* out, const float* in, uint8_t* temp)
- {
- try_fn([&]() { reinterpret_cast<kfr::dct_plan<float>*>(plan)->execute(out, in, temp, kfr::cfalse); });
- }
- KFR_API_SPEC void kfr_dct_execute_f64(KFR_DCT_PLAN_F64* plan, double* out, const double* in,
- uint8_t* temp)
- {
- try_fn([&]()
- { reinterpret_cast<kfr::dct_plan<double>*>(plan)->execute(out, in, temp, kfr::cfalse); });
- }
- KFR_API_SPEC void kfr_dct_execute_inverse_f32(KFR_DCT_PLAN_F32* plan, float* out, const float* in,
- uint8_t* temp)
- {
- try_fn([&]() { reinterpret_cast<kfr::dct_plan<float>*>(plan)->execute(out, in, temp, kfr::ctrue); });
- }
- KFR_API_SPEC void kfr_dct_execute_inverse_f64(KFR_DCT_PLAN_F64* plan, double* out, const double* in,
- uint8_t* temp)
- {
- try_fn([&]() { reinterpret_cast<kfr::dct_plan<double>*>(plan)->execute(out, in, temp, kfr::ctrue); });
- }
+KFR_API_SPEC void kfr_dft_real_dump_f32(KFR_DFT_REAL_PLAN_F32* plan)
+{
+ try_fn([&]() { reinterpret_cast<var_dft_plan<float>*>(plan)->dump(); });
+}
+KFR_API_SPEC void kfr_dft_real_dump_f64(KFR_DFT_REAL_PLAN_F64* plan)
+{
+ try_fn([&]() { reinterpret_cast<var_dft_plan<double>*>(plan)->dump(); });
+}
- KFR_API_SPEC void kfr_dct_delete_plan_f32(KFR_DCT_PLAN_F32* plan)
- {
- try_fn([&]() { delete reinterpret_cast<kfr::dct_plan<float>*>(plan); });
- }
- KFR_API_SPEC void kfr_dct_delete_plan_f64(KFR_DCT_PLAN_F64* plan)
- {
- try_fn([&]() { delete reinterpret_cast<kfr::dct_plan<double>*>(plan); });
- }
+KFR_API_SPEC size_t kfr_dft_real_get_size_f32(KFR_DFT_REAL_PLAN_F32* plan)
+{
+ return try_fn([&]() { return reinterpret_cast<var_dft_plan<float>*>(plan)->size(); }, 0);
+}
+KFR_API_SPEC size_t kfr_dft_real_get_size_f64(KFR_DFT_REAL_PLAN_F64* plan)
+{
+ return try_fn([&]() { return reinterpret_cast<var_dft_plan<double>*>(plan)->size(); }, 0);
+}
- // Filters
+KFR_API_SPEC size_t kfr_dft_real_get_temp_size_f32(KFR_DFT_REAL_PLAN_F32* plan)
+{
+ return try_fn([&]() { return reinterpret_cast<var_dft_plan<float>*>(plan)->temp_size(); }, 0);
+}
+KFR_API_SPEC size_t kfr_dft_real_get_temp_size_f64(KFR_DFT_REAL_PLAN_F64* plan)
+{
+ return try_fn([&]() { return reinterpret_cast<var_dft_plan<double>*>(plan)->temp_size(); }, 0);
+}
- KFR_API_SPEC KFR_FILTER_F32* kfr_filter_create_fir_plan_f32(const kfr_f32* taps, size_t size)
- {
- return try_fn(
- [&]()
- { return reinterpret_cast<KFR_FILTER_F32*>(new fir_filter<float>(make_univector(taps, size))); },
- nullptr);
- }
- KFR_API_SPEC KFR_FILTER_F64* kfr_filter_create_fir_plan_f64(const kfr_f64* taps, size_t size)
- {
- return try_fn(
- [&]()
- { return reinterpret_cast<KFR_FILTER_F64*>(new fir_filter<double>(make_univector(taps, size))); },
- nullptr);
- }
+KFR_API_SPEC void kfr_dft_real_execute_f32(KFR_DFT_REAL_PLAN_F32* plan, kfr_c32* out, const float* in,
+ uint8_t* temp)
+{
+ try_fn(
+ [&]()
+ { reinterpret_cast<var_dft_plan<float>*>(plan)->execute(reinterpret_cast<float*>(out), in, temp); });
+}
+KFR_API_SPEC void kfr_dft_real_execute_f64(KFR_DFT_REAL_PLAN_F64* plan, kfr_c64* out, const double* in,
+ uint8_t* temp)
+{
+ try_fn(
+ [&]() {
+ reinterpret_cast<var_dft_plan<double>*>(plan)->execute(reinterpret_cast<double*>(out), in, temp);
+ });
+}
+KFR_API_SPEC void kfr_dft_real_execute_inverse_f32(KFR_DFT_REAL_PLAN_F32* plan, float* out, const kfr_c32* in,
+ uint8_t* temp)
+{
+ try_fn(
+ [&]()
+ {
+ reinterpret_cast<var_dft_plan<float>*>(plan)->execute_inverse(
+ out, reinterpret_cast<const float*>(in), temp);
+ });
+}
+KFR_API_SPEC void kfr_dft_real_execute_inverse_f64(KFR_DFT_REAL_PLAN_F64* plan, double* out,
+ const kfr_c64* in, uint8_t* temp)
+{
+ try_fn(
+ [&]()
+ {
+ reinterpret_cast<var_dft_plan<double>*>(plan)->execute_inverse(
+ out, reinterpret_cast<const double*>(in), temp);
+ });
+}
- KFR_API_SPEC KFR_FILTER_F32* kfr_filter_create_convolution_plan_f32(const kfr_f32* taps, size_t size,
- size_t block_size)
- {
- return try_fn(
- [&]()
- {
- return reinterpret_cast<KFR_FILTER_F32*>(
- new convolve_filter<float>(make_univector(taps, size), block_size ? block_size : 1024));
- },
- nullptr);
- }
- KFR_API_SPEC KFR_FILTER_F64* kfr_filter_create_convolution_plan_f64(const kfr_f64* taps, size_t size,
- size_t block_size)
- {
- return try_fn(
- [&]()
- {
- return reinterpret_cast<KFR_FILTER_F64*>(
- new convolve_filter<double>(make_univector(taps, size), block_size ? block_size : 1024));
- },
- nullptr);
- }
+KFR_API_SPEC void kfr_dft_real_delete_plan_f32(KFR_DFT_REAL_PLAN_F32* plan)
+{
+ try_fn([&]() { delete reinterpret_cast<var_dft_plan<float>*>(plan); });
+}
+KFR_API_SPEC void kfr_dft_real_delete_plan_f64(KFR_DFT_REAL_PLAN_F64* plan)
+{
+ try_fn([&]() { delete reinterpret_cast<var_dft_plan<double>*>(plan); });
+}
- KFR_API_SPEC KFR_FILTER_F32* kfr_filter_create_iir_plan_f32(const kfr_f32* sos, size_t sos_count)
- {
- return try_fn(
- [&]()
- {
- return reinterpret_cast<KFR_FILTER_F32*>(new iir_filter<float>(
- iir_params{ reinterpret_cast<const biquad_section<float>*>(sos), sos_count }));
- },
- nullptr);
- }
- KFR_API_SPEC KFR_FILTER_F64* kfr_filter_create_iir_plan_f64(const kfr_f64* sos, size_t sos_count)
- {
- return try_fn(
- [&]()
- {
- return reinterpret_cast<KFR_FILTER_F64*>(new iir_filter<double>(
- iir_params{ reinterpret_cast<const biquad_section<double>*>(sos), sos_count }));
- },
- nullptr);
- }
+// Discrete Cosine Transform
- KFR_API_SPEC void kfr_filter_process_f32(KFR_FILTER_F32* plan, kfr_f32* output, const kfr_f32* input,
- size_t size)
- {
- try_fn([&]() { reinterpret_cast<filter<float>*>(plan)->apply(output, input, size); });
- }
- KFR_API_SPEC void kfr_filter_process_f64(KFR_FILTER_F64* plan, kfr_f64* output, const kfr_f64* input,
- size_t size)
- {
- try_fn([&]() { reinterpret_cast<filter<double>*>(plan)->apply(output, input, size); });
- }
+KFR_API_SPEC KFR_DCT_PLAN_F32* kfr_dct_create_plan_f32(size_t size)
+{
+ return try_fn([&]() { return reinterpret_cast<KFR_DCT_PLAN_F32*>(new dct_plan<float>(size)); }, nullptr);
+}
+KFR_API_SPEC KFR_DCT_PLAN_F64* kfr_dct_create_plan_f64(size_t size)
+{
+ return try_fn([&]() { return reinterpret_cast<KFR_DCT_PLAN_F64*>(new dct_plan<double>(size)); }, nullptr);
+}
- KFR_API_SPEC void kfr_filter_reset_f32(KFR_FILTER_F32* plan)
- {
- try_fn([&]() { reinterpret_cast<filter<float>*>(plan)->reset(); });
- }
- KFR_API_SPEC void kfr_filter_reset_f64(KFR_FILTER_F64* plan)
- {
- try_fn([&]() { reinterpret_cast<filter<double>*>(plan)->reset(); });
- }
+KFR_API_SPEC void kfr_dct_dump_f32(KFR_DCT_PLAN_F32* plan)
+{
+ try_fn([&]() { reinterpret_cast<dct_plan<float>*>(plan)->dump(); });
+}
+KFR_API_SPEC void kfr_dct_dump_f64(KFR_DCT_PLAN_F64* plan)
+{
+ try_fn([&]() { reinterpret_cast<dct_plan<double>*>(plan)->dump(); });
+}
- KFR_API_SPEC void kfr_filter_delete_plan_f32(KFR_FILTER_F32* plan)
- {
- try_fn([&]() { delete reinterpret_cast<filter<f32>*>(plan); });
- }
- KFR_API_SPEC void kfr_filter_delete_plan_f64(KFR_FILTER_F64* plan)
- {
- try_fn([&]() { delete reinterpret_cast<filter<f64>*>(plan); });
- }
+KFR_API_SPEC size_t kfr_dct_get_size_f32(KFR_DCT_PLAN_F32* plan)
+{
+ return try_fn([&]() { return reinterpret_cast<dft_plan<float>*>(plan)->size; }, 0);
+}
+KFR_API_SPEC size_t kfr_dct_get_size_f64(KFR_DCT_PLAN_F64* plan)
+{
+ return try_fn([&]() { return reinterpret_cast<dft_plan<double>*>(plan)->size; }, 0);
+}
+
+KFR_API_SPEC size_t kfr_dct_get_temp_size_f32(KFR_DCT_PLAN_F32* plan)
+{
+ return try_fn([&]() { return reinterpret_cast<dft_plan<float>*>(plan)->temp_size; }, 0);
+}
+KFR_API_SPEC size_t kfr_dct_get_temp_size_f64(KFR_DCT_PLAN_F64* plan)
+{
+ return try_fn([&]() { return reinterpret_cast<dft_plan<double>*>(plan)->temp_size; }, 0);
+}
+
+KFR_API_SPEC void kfr_dct_execute_f32(KFR_DCT_PLAN_F32* plan, float* out, const float* in, uint8_t* temp)
+{
+ try_fn([&]() { reinterpret_cast<dct_plan<float>*>(plan)->execute(out, in, temp, cfalse); });
+}
+KFR_API_SPEC void kfr_dct_execute_f64(KFR_DCT_PLAN_F64* plan, double* out, const double* in, uint8_t* temp)
+{
+ try_fn([&]() { reinterpret_cast<dct_plan<double>*>(plan)->execute(out, in, temp, cfalse); });
+}
+KFR_API_SPEC void kfr_dct_execute_inverse_f32(KFR_DCT_PLAN_F32* plan, float* out, const float* in,
+ uint8_t* temp)
+{
+ try_fn([&]() { reinterpret_cast<dct_plan<float>*>(plan)->execute(out, in, temp, ctrue); });
+}
+KFR_API_SPEC void kfr_dct_execute_inverse_f64(KFR_DCT_PLAN_F64* plan, double* out, const double* in,
+ uint8_t* temp)
+{
+ try_fn([&]() { reinterpret_cast<dct_plan<double>*>(plan)->execute(out, in, temp, ctrue); });
+}
+
+KFR_API_SPEC void kfr_dct_delete_plan_f32(KFR_DCT_PLAN_F32* plan)
+{
+ try_fn([&]() { delete reinterpret_cast<dct_plan<float>*>(plan); });
+}
+KFR_API_SPEC void kfr_dct_delete_plan_f64(KFR_DCT_PLAN_F64* plan)
+{
+ try_fn([&]() { delete reinterpret_cast<dct_plan<double>*>(plan); });
+}
+
+// Filters
+
+KFR_API_SPEC KFR_FILTER_F32* kfr_filter_create_fir_plan_f32(const kfr_f32* taps, size_t size)
+{
+ return try_fn(
+ [&]()
+ { return reinterpret_cast<KFR_FILTER_F32*>(new fir_filter<float>(make_univector(taps, size))); },
+ nullptr);
+}
+KFR_API_SPEC KFR_FILTER_F64* kfr_filter_create_fir_plan_f64(const kfr_f64* taps, size_t size)
+{
+ return try_fn(
+ [&]()
+ { return reinterpret_cast<KFR_FILTER_F64*>(new fir_filter<double>(make_univector(taps, size))); },
+ nullptr);
+}
+
+KFR_API_SPEC KFR_FILTER_F32* kfr_filter_create_convolution_plan_f32(const kfr_f32* taps, size_t size,
+ size_t block_size)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_FILTER_F32*>(
+ new convolve_filter<float>(make_univector(taps, size), block_size ? block_size : 1024));
+ },
+ nullptr);
+}
+KFR_API_SPEC KFR_FILTER_F64* kfr_filter_create_convolution_plan_f64(const kfr_f64* taps, size_t size,
+ size_t block_size)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_FILTER_F64*>(
+ new convolve_filter<double>(make_univector(taps, size), block_size ? block_size : 1024));
+ },
+ nullptr);
+}
+
+KFR_API_SPEC KFR_FILTER_F32* kfr_filter_create_iir_plan_f32(const kfr_f32* sos, size_t sos_count)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_FILTER_F32*>(new iir_filter<float>(
+ iir_params{ reinterpret_cast<const biquad_section<float>*>(sos), sos_count }));
+ },
+ nullptr);
+}
+KFR_API_SPEC KFR_FILTER_F64* kfr_filter_create_iir_plan_f64(const kfr_f64* sos, size_t sos_count)
+{
+ return try_fn(
+ [&]()
+ {
+ return reinterpret_cast<KFR_FILTER_F64*>(new iir_filter<double>(
+ iir_params{ reinterpret_cast<const biquad_section<double>*>(sos), sos_count }));
+ },
+ nullptr);
+}
+
+KFR_API_SPEC void kfr_filter_process_f32(KFR_FILTER_F32* plan, kfr_f32* output, const kfr_f32* input,
+ size_t size)
+{
+ try_fn([&]() { reinterpret_cast<filter<float>*>(plan)->apply(output, input, size); });
+}
+KFR_API_SPEC void kfr_filter_process_f64(KFR_FILTER_F64* plan, kfr_f64* output, const kfr_f64* input,
+ size_t size)
+{
+ try_fn([&]() { reinterpret_cast<filter<double>*>(plan)->apply(output, input, size); });
+}
+
+KFR_API_SPEC void kfr_filter_reset_f32(KFR_FILTER_F32* plan)
+{
+ try_fn([&]() { reinterpret_cast<filter<float>*>(plan)->reset(); });
+}
+KFR_API_SPEC void kfr_filter_reset_f64(KFR_FILTER_F64* plan)
+{
+ try_fn([&]() { reinterpret_cast<filter<double>*>(plan)->reset(); });
+}
+
+KFR_API_SPEC void kfr_filter_delete_plan_f32(KFR_FILTER_F32* plan)
+{
+ try_fn([&]() { delete reinterpret_cast<filter<f32>*>(plan); });
+}
+KFR_API_SPEC void kfr_filter_delete_plan_f64(KFR_FILTER_F64* plan)
+{
+ try_fn([&]() { delete reinterpret_cast<filter<f64>*>(plan); });
+}
}
} // namespace kfr
diff --git a/tests/capi_test.cpp b/tests/capi_test.cpp
@@ -47,24 +47,26 @@ void test_memory()
#define DFT_SIZE 32
-void test_dft_f32()
+static void init_data_f32(kfr_f32* buf)
{
- printf("[TEST] DFT f32\n");
- const float eps = 0.00001f;
- KFR_DFT_PLAN_F32* plan = kfr_dft_create_plan_f32(DFT_SIZE);
- // kfr_dft_dump_f32(plan);
- kfr_f32 buf[DFT_SIZE * 2];
for (int i = 0; i < DFT_SIZE; i++)
{
buf[i * 2 + 0] = (float)(i) / DFT_SIZE;
buf[i * 2 + 1] = (float)(-i) / DFT_SIZE;
}
- uint8_t* tmp = (uint8_t*)kfr_allocate(kfr_dft_get_temp_size_f32(plan));
- kfr_dft_execute_f32(plan, buf, buf, tmp);
- kfr_dft_execute_inverse_f32(plan, buf, buf, tmp);
- kfr_dft_delete_plan_f32(plan);
- kfr_deallocate(tmp);
+}
+static void init_data_f64(kfr_f64* buf)
+{
+ for (int i = 0; i < DFT_SIZE; i++)
+ {
+ buf[i * 2 + 0] = (double)(i) / DFT_SIZE;
+ buf[i * 2 + 1] = (double)(-i) / DFT_SIZE;
+ }
+}
+static void test_data_f32(kfr_f32* buf)
+{
+ const float eps = 0.00001f;
for (int i = 0; i < DFT_SIZE; i++)
{
CHECK(fabsf(buf[i * 2 + 0] - (float)(i)) < eps, "DFT: wrong result at %d: re = %f", i,
@@ -73,32 +75,108 @@ void test_dft_f32()
buf[i * 2 + 1]);
}
}
+static void test_data_f64(kfr_f64* buf)
+{
+ const double eps = 0.00001;
+ for (int i = 0; i < DFT_SIZE; i++)
+ {
+ CHECK(fabs(buf[i * 2 + 0] - (double)(i)) < eps, "DFT: wrong result at %d: re = %f", i,
+ buf[i * 2 + 0]);
+ CHECK(fabs(buf[i * 2 + 1] - (double)(-i)) < eps, "DFT: wrong result at %d: im = %f", i,
+ buf[i * 2 + 1]);
+ }
+}
+
+void test_dft_f32()
+{
+ printf("[TEST] DFT f32\n");
+ // kfr_dft_dump_f32(plan);
+ kfr_f32 buf[DFT_SIZE * 2];
+ KFR_DFT_PLAN_F32* plan;
+ uint8_t* tmp;
+
+ init_data_f32(buf);
+ plan = kfr_dft_create_plan_f32(DFT_SIZE);
+ tmp = (uint8_t*)kfr_allocate(kfr_dft_get_temp_size_f32(plan));
+ kfr_dft_execute_f32(plan, buf, buf, tmp);
+ kfr_dft_execute_inverse_f32(plan, buf, buf, tmp);
+ kfr_deallocate(tmp);
+ kfr_dft_delete_plan_f32(plan);
+ test_data_f32(buf);
+
+ init_data_f32(buf);
+ plan = kfr_dft_create_2d_plan_f32(8, 4);
+ tmp = (uint8_t*)kfr_allocate(kfr_dft_get_temp_size_f32(plan));
+ kfr_dft_execute_f32(plan, buf, buf, tmp);
+ kfr_dft_execute_inverse_f32(plan, buf, buf, tmp);
+ kfr_deallocate(tmp);
+ kfr_dft_delete_plan_f32(plan);
+ test_data_f32(buf);
+
+ init_data_f32(buf);
+ plan = kfr_dft_create_3d_plan_f32(4, 4, 2);
+ tmp = (uint8_t*)kfr_allocate(kfr_dft_get_temp_size_f32(plan));
+ kfr_dft_execute_f32(plan, buf, buf, tmp);
+ kfr_dft_execute_inverse_f32(plan, buf, buf, tmp);
+ kfr_deallocate(tmp);
+ kfr_dft_delete_plan_f32(plan);
+ test_data_f32(buf);
+
+ unsigned sizes[4] = { 2, 2, 2, 4 };
+ init_data_f32(buf);
+ plan = kfr_dft_create_md_plan_f32(4, sizes);
+ tmp = (uint8_t*)kfr_allocate(kfr_dft_get_temp_size_f32(plan));
+ kfr_dft_execute_f32(plan, buf, buf, tmp);
+ kfr_dft_execute_inverse_f32(plan, buf, buf, tmp);
+ kfr_deallocate(tmp);
+ kfr_dft_delete_plan_f32(plan);
+ test_data_f32(buf);
+}
void test_dft_f64()
{
printf("[TEST] DFT f64\n");
- const float eps = 0.00000001;
- KFR_DFT_PLAN_F64* plan = kfr_dft_create_plan_f64(DFT_SIZE);
// kfr_dft_dump_f64(plan);
kfr_f64 buf[DFT_SIZE * 2];
- for (int i = 0; i < DFT_SIZE; i++)
- {
- buf[i * 2 + 0] = (double)(i) / DFT_SIZE;
- buf[i * 2 + 1] = (double)(-i) / DFT_SIZE;
- }
- uint8_t* tmp = (uint8_t*)kfr_allocate(kfr_dft_get_temp_size_f64(plan));
+ KFR_DFT_PLAN_F64* plan;
+ uint8_t* tmp;
+
+ init_data_f64(buf);
+ plan = kfr_dft_create_plan_f64(DFT_SIZE);
+ tmp = (uint8_t*)kfr_allocate(kfr_dft_get_temp_size_f64(plan));
+ kfr_dft_execute_f64(plan, buf, buf, tmp);
+ kfr_dft_execute_inverse_f64(plan, buf, buf, tmp);
+ kfr_deallocate(tmp);
+ kfr_dft_delete_plan_f64(plan);
+ test_data_f64(buf);
+
+ init_data_f64(buf);
+ plan = kfr_dft_create_2d_plan_f64(8, 4);
+ tmp = (uint8_t*)kfr_allocate(kfr_dft_get_temp_size_f64(plan));
kfr_dft_execute_f64(plan, buf, buf, tmp);
kfr_dft_execute_inverse_f64(plan, buf, buf, tmp);
+ kfr_deallocate(tmp);
kfr_dft_delete_plan_f64(plan);
+ test_data_f64(buf);
+
+ init_data_f64(buf);
+ plan = kfr_dft_create_3d_plan_f64(4, 4, 2);
+ tmp = (uint8_t*)kfr_allocate(kfr_dft_get_temp_size_f64(plan));
+ kfr_dft_execute_f64(plan, buf, buf, tmp);
+ kfr_dft_execute_inverse_f64(plan, buf, buf, tmp);
kfr_deallocate(tmp);
+ kfr_dft_delete_plan_f64(plan);
+ test_data_f64(buf);
- for (int i = 0; i < DFT_SIZE; i++)
- {
- CHECK(fabs(buf[i * 2 + 0] - (double)(i)) < eps, "DFT: wrong result at %d: re = %f", i,
- buf[i * 2 + 0]);
- CHECK(fabs(buf[i * 2 + 1] - (double)(-i)) < eps, "DFT: wrong result at %d: im = %f", i,
- buf[i * 2 + 1]);
- }
+ unsigned sizes[4] = { 2, 2, 2, 4 };
+ init_data_f64(buf);
+ plan = kfr_dft_create_md_plan_f64(4, sizes);
+ tmp = (uint8_t*)kfr_allocate(kfr_dft_get_temp_size_f64(plan));
+ kfr_dft_execute_f64(plan, buf, buf, tmp);
+ kfr_dft_execute_inverse_f64(plan, buf, buf, tmp);
+ kfr_deallocate(tmp);
+ kfr_dft_delete_plan_f64(plan);
+ test_data_f64(buf);
}
#define FILTER_SIZE 256