WarpX
Communication.H
Go to the documentation of this file.
1 /* Copyright 2019-2022 Andrew Myers
2  *
3  * This file is part of WarpX.
4  *
5  * License: BSD-3-Clause-LBNL
6  */
7 #ifndef ABLASTR_UTILS_COMMUNICATION_H_
8 #define ABLASTR_UTILS_COMMUNICATION_H_
9 
10 #include <AMReX_FabArrayBase.H>
11 #include <AMReX_GpuDevice.H>
12 #include <AMReX_GpuQualifiers.H>
13 #include <AMReX_Periodicity.H>
14 #include <AMReX_Vector.H>
15 
16 #include <AMReX_BaseFwd.H>
17 
18 #include <optional>
19 
20 
22 {
23 
25 
26 template <class FAB1, class FAB2>
27 void
28 mixedCopy (amrex::FabArray<FAB1>& dst, amrex::FabArray<FAB2> const& src, int srccomp, int dstcomp, int numcomp, const amrex::IntVect& nghost)
29 {
30  auto const& srcma = src.const_arrays();
31  auto const& dstma = dst.arrays();
32  ParallelFor(dst, nghost, numcomp,
33  [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
34  {
35  dstma[box_no](i,j,k,dstcomp+n) = (typename FAB1::value_type) srcma[box_no](i,j,k,srccomp+n);
36  });
38 }
39 
41  const amrex::MultiFab &src,
42  int src_comp,
43  int dst_comp,
44  int num_comp,
45  const amrex::IntVect &src_nghost,
46  const amrex::IntVect &dst_nghost,
47  bool do_single_precision_comms,
50 
51 void ParallelAdd (amrex::MultiFab &dst,
52  const amrex::MultiFab &src,
53  int src_comp,
54  int dst_comp,
55  int num_comp,
56  const amrex::IntVect &src_nghost,
57  const amrex::IntVect &dst_nghost,
58  bool do_single_precision_comms,
60 
62  bool do_single_precision_comms,
64  std::optional<bool> nodal_sync=std::nullopt);
65 
67  amrex::IntVect ng,
68  bool do_single_precision_comms,
70  std::optional<bool> nodal_sync = std::nullopt);
71 
74 
76  amrex::IntVect ng,
78 
79 void
80 FillBoundary(amrex::Vector<amrex::MultiFab *> const &mf, bool do_single_precision_comms,
81  const amrex::Periodicity &period, std::optional<bool> nodal_sync=std::nullopt);
82 
83 void
85  int start_comp,
86  int num_comps,
87  amrex::IntVect src_ng,
88  amrex::IntVect dst_ng,
89  bool do_single_precision_comms,
91 
93  bool do_single_precision_comms,
95 }
96 
97 #endif // ABLASTR_UTILS_COMMUNICATION_H_
#define AMREX_GPU_DEVICE
MultiArray4< typename FabArray< FAB >::value_type > arrays() noexcept
MultiArray4< typename FabArray< FAB >::value_type const > const_arrays() const noexcept
static const Periodicity & NonPeriodic() noexcept
Definition: Communication.cpp:26
void OverrideSync(amrex::MultiFab &mf, bool do_single_precision_comms, const amrex::Periodicity &period)
Definition: Communication.cpp:177
void FillBoundary(amrex::MultiFab &mf, amrex::IntVect ng, bool do_single_precision_comms, const amrex::Periodicity &period, std::optional< bool > nodal_sync)
Definition: Communication.cpp:71
void mixedCopy(amrex::FabArray< FAB1 > &dst, amrex::FabArray< FAB2 > const &src, int srccomp, int dstcomp, int numcomp, const amrex::IntVect &nghost)
Definition: Communication.H:28
float comm_float_type
Definition: Communication.H:24
void ParallelAdd(amrex::MultiFab &dst, const amrex::MultiFab &src, int src_comp, int dst_comp, int num_comp, const amrex::IntVect &src_nghost, const amrex::IntVect &dst_nghost, bool do_single_precision_comms, const amrex::Periodicity &period)
Definition: Communication.cpp:63
void SumBoundary(amrex::MultiFab &mf, int start_comp, int num_comps, amrex::IntVect src_ng, amrex::IntVect dst_ng, bool do_single_precision_comms, const amrex::Periodicity &period)
Definition: Communication.cpp:149
void ParallelCopy(amrex::MultiFab &dst, const amrex::MultiFab &src, int src_comp, int dst_comp, int num_comp, const amrex::IntVect &src_nghost, const amrex::IntVect &dst_nghost, bool do_single_precision_comms, const amrex::Periodicity &period, amrex::FabArrayBase::CpOp op)
Definition: Communication.cpp:28
void synchronize() noexcept
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... >, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
i
Definition: check_interp_points_and_weights.py:174
float
Definition: plot_parallel.py:59