IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
CIC.h
Go to the documentation of this file.
1//
2// Class CIC
3// First order/cloud-in-cell grid interpolation. Currently implemented as
4// global functions, but in order to support higher or lower order interpolation,
5// these should be moved into structs.
6//
7#ifndef CIC_INTERPOLATION_H
8#define CIC_INTERPOLATION_H
9
10namespace ippl {
11 namespace detail {
21 template <unsigned long Point, unsigned long Index, typename Weights>
22 KOKKOS_INLINE_FUNCTION constexpr typename Weights::value_type interpolationWeight(
23 const Weights& wlo, const Weights& whi);
24
33 template <unsigned long Point, unsigned long Index, typename Indices>
34 KOKKOS_INLINE_FUNCTION constexpr typename Indices::value_type interpolationIndex(
35 const Indices& args);
36
50 template <unsigned long ScatterPoint, unsigned long... Index, typename View, typename T,
51 typename IndexType = size_t>
52 KOKKOS_INLINE_FUNCTION constexpr void scatterToPoint(
53 const std::index_sequence<Index...>&, const View& view,
54 const Vector<T, View::rank>& wlo, const Vector<T, View::rank>& whi,
55 const Vector<IndexType, View::rank>& args, const T& val);
56
77 template <unsigned long... ScatterPoint, typename View, typename T,
78 typename IndexType = size_t>
79 KOKKOS_INLINE_FUNCTION constexpr void scatterToField(
80 const std::index_sequence<ScatterPoint...>&, const View& view,
81 const Vector<T, View::rank>& wlo, const Vector<T, View::rank>& whi,
82 const Vector<IndexType, View::rank>& args, T val = 1);
83
97 template <unsigned long GatherPoint, unsigned long... Index, typename View, typename T,
98 typename IndexType = size_t>
99 KOKKOS_INLINE_FUNCTION constexpr typename View::value_type gatherFromPoint(
100 const std::index_sequence<Index...>&, const View& view,
101 const Vector<T, View::rank>& wlo, const Vector<T, View::rank>& whi,
103
116 template <unsigned long... GatherPoint, typename View, typename T,
117 typename IndexType = size_t>
118 KOKKOS_INLINE_FUNCTION constexpr typename View::value_type gatherFromField(
119 const std::index_sequence<GatherPoint...>&, const View& view,
120 const Vector<T, View::rank>& wlo, const Vector<T, View::rank>& whi,
122 } // namespace detail
123} // namespace ippl
124
125#include "Interpolation/CIC.hpp"
126
127#endif
ippl::Vector< T, Dim > Vector
Definition datatypes.h:26
Definition Archive.h:20
KOKKOS_INLINE_FUNCTION constexpr Weights::value_type interpolationWeight(const Weights &wlo, const Weights &whi)
Definition CIC.hpp:11
KOKKOS_INLINE_FUNCTION constexpr View::value_type gatherFromPoint(const std::index_sequence< Index... > &, const View &view, const Vector< T, View::rank > &wlo, const Vector< T, View::rank > &whi, const Vector< IndexType, View::rank > &args)
Definition CIC.hpp:59
KOKKOS_INLINE_FUNCTION constexpr View::value_type gatherFromField(const std::index_sequence< GatherPoint... > &, const View &view, const Vector< T, View::rank > &wlo, const Vector< T, View::rank > &whi, const Vector< IndexType, View::rank > &args)
Definition CIC.hpp:68
KOKKOS_INLINE_FUNCTION constexpr Indices::value_type interpolationIndex(const Indices &args)
Definition CIC.hpp:24
KOKKOS_INLINE_FUNCTION constexpr void scatterToPoint(const std::index_sequence< Index... > &, const View &view, const Vector< T, View::rank > &wlo, const Vector< T, View::rank > &whi, const Vector< IndexType, View::rank > &args, const T &val)
Definition CIC.hpp:38
KOKKOS_INLINE_FUNCTION constexpr void scatterToField(const std::index_sequence< ScatterPoint... > &, const View &view, const Vector< T, View::rank > &wlo, const Vector< T, View::rank > &whi, const Vector< IndexType, View::rank > &args, T val=1)
Definition CIC.hpp:47