1 #ifndef TOROIDAL_AVERAGE_HPP
2 #define TOROIDAL_AVERAGE_HPP
12 if(pol_decomp.mpi.n_intpl_ranks==1)
return;
15 View<double*,CLayout, typename T::device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged>> f_1D((
double*)(f.data()), f.size());
25 MPI_Allreduce(MPI_IN_PLACE, f_mpi.data(), f_mpi.size(), MPI_DOUBLE, MPI_SUM, pol_decomp.mpi.intpl_comm);
37 if(pol_decomp.mpi.n_intpl_ranks==1)
return;
40 View<double*,CLayout, typename T::device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged>> f_1D((
double*)(f.data()), f.size());
50 MPI_Allreduce(MPI_IN_PLACE, f_mpi.data(), f_mpi.size(), MPI_DOUBLE, MPI_SUM, pol_decomp.mpi.intpl_comm);
57 double inv_n_planes = 1.0/pol_decomp.mpi.n_intpl_ranks;
60 Kokkos::parallel_for(
"toroidal_average", Kokkos::RangePolicy<typename T::execution_space>(0, f.size()), KOKKOS_LAMBDA(
const int i){
61 f_1D(i) *= inv_n_planes;
70 View<double*,CLayout, typename T::device_type, Kokkos::MemoryTraits<Kokkos::Unmanaged>> f_1D((
double*)(f.data()), f.size());
73 View<double*,CLayout, typename T::device_type> toroidal_avg(
NoInit(
"toroidal_avg"), f_1D.layout());
74 Kokkos::deep_copy(toroidal_avg, f_1D);
80 Kokkos::parallel_for(
"toroidal_average", Kokkos::RangePolicy<typename T::execution_space>(0, f.size()), KOKKOS_LAMBDA(
const int i){
81 f_1D(i) -= toroidal_avg(i);
void mirror_copy(T1 &view_dest, const T2 &view_src)
Definition: my_mirror_view.hpp:122
void toroidal_average_in_place(const DomainDecomposition< DeviceType > &pol_decomp, T &f)
Definition: toroidal_average.hpp:34
Kokkos::Device< HostExSpace, HostMemSpace > HostType
Definition: space_settings.hpp:57
subroutine remove_toroidal_average(n, x)
Removes the toroidal average of a quantity of type real(8) that is global on a poloidal plane but loc...
Definition: tiny_functions.F90:72
void toroidal_sum_in_place(const DomainDecomposition< DeviceType > &pol_decomp, T &f)
Definition: toroidal_average.hpp:9
View< T *, CLayout, Device > my_mirror_scratch_view(const View< T *, CLayout, Device, Kokkos::MemoryTraits< Kokkos::Unmanaged >> &view, Device nd)
Definition: my_mirror_view.hpp:97
View< T *, CLayout, Device > my_mirror_view(const View< T *, CLayout, Device > &view, Device nd)
Definition: my_mirror_view.hpp:14
void parallel_for(const std::string name, int n_ptl, Function func, Option option, HostAoSoA aosoa_h, DeviceAoSoA aosoa_d)
Definition: streamed_parallel_for.hpp:252
View< double *, CLayout, typename T::device_type > split_toroidal_average(const DomainDecomposition< DeviceType > &pol_decomp, T &f)
Definition: toroidal_average.hpp:68
Kokkos::ViewAllocateWithoutInitializing NoInit
Definition: space_settings.hpp:69