ARTEMIS
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_FabArray.H>
11 #include <AMReX_Gpu.H>
12 #include <AMReX_GpuDevice.H>
13 #include <AMReX_iMultiFab.H>
14 #include <AMReX_IntVect.H>
15 #include <AMReX_MultiFab.H>
16 #include <AMReX_Periodicity.H>
17 #include <AMReX_TypeTraits.H>
18 
19 #include "WarpX.H"
20 
21 #include <optional>
22 
23 
25 {
26 
28 
29 template <class FAB1, class FAB2>
30 void
31 mixedCopy (amrex::FabArray<FAB1>& dst, amrex::FabArray<FAB2> const& src, int srccomp, int dstcomp, int numcomp, const amrex::IntVect& nghost)
32 {
33  auto const& srcma = src.const_arrays();
34  auto const& dstma = dst.arrays();
35  ParallelFor(dst, nghost, numcomp,
36  [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
37  {
38  dstma[box_no](i,j,k,dstcomp+n) = (typename FAB1::value_type) srcma[box_no](i,j,k,srccomp+n);
39  });
41 }
42 
44  const amrex::MultiFab &src,
45  int src_comp,
46  int dst_comp,
47  int num_comp,
48  const amrex::IntVect &src_nghost,
49  const amrex::IntVect &dst_nghost,
50  bool do_single_precision_comms,
53 
54 void ParallelAdd (amrex::MultiFab &dst,
55  const amrex::MultiFab &src,
56  int src_comp,
57  int dst_comp,
58  int num_comp,
59  const amrex::IntVect &src_nghost,
60  const amrex::IntVect &dst_nghost,
61  bool do_single_precision_comms,
63 
65  bool do_single_precision_comms,
67 
69  amrex::IntVect ng,
70  bool do_single_precision_comms,
72  std::optional<bool> nodal_sync = std::nullopt);
73 
76 
78  amrex::IntVect ng,
80 
81 void
82 FillBoundary(amrex::Vector<amrex::MultiFab *> const &mf, bool do_single_precision_comms,
83  const amrex::Periodicity &period);
84 
86  bool do_single_precision_comms,
88 
90  int start_comp,
91  int num_comps, amrex::IntVect ng,
92  bool do_single_precision_comms,
94 
95 void
97  int start_comp,
98  int num_comps,
99  amrex::IntVect src_ng,
100  amrex::IntVect dst_ng,
101  bool do_single_precision_comms,
103 
104 void OverrideSync (amrex::MultiFab &mf,
105  bool do_single_precision_comms,
107 }
108 
109 #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:20
void OverrideSync(amrex::MultiFab &mf, bool do_single_precision_comms, const amrex::Periodicity &period)
Definition: Communication.cpp:221
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:65
void mixedCopy(amrex::FabArray< FAB1 > &dst, amrex::FabArray< FAB2 > const &src, int srccomp, int dstcomp, int numcomp, const amrex::IntVect &nghost)
Definition: Communication.H:31
float comm_float_type
Definition: Communication.H:27
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:57
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:22
void SumBoundary(amrex::MultiFab &mf, bool do_single_precision_comms, const amrex::Periodicity &period)
Definition: Communication.cpp:142
void synchronize() noexcept
i
Definition: check_interp_points_and_weights.py:174
float
Definition: plot_parallel.py:59
int n
Definition: run_libensemble_on_warpx.py:67