aboutsummaryrefslogtreecommitdiff
path: root/Source/ablastr/utils/Communication.H
blob: 9be9fec1e58b014c1ccb29d19d1c3cddf96bc656 (plain) (blame)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
/* Copyright 2019-2022 Andrew Myers
 *
 * This file is part of WarpX.
 *
 * License: BSD-3-Clause-LBNL
 */
#ifndef ABLASTR_UTILS_COMMUNICATION_H_
#define ABLASTR_UTILS_COMMUNICATION_H_

#include <AMReX_FabArrayBase.H>
#include <AMReX_GpuDevice.H>
#include <AMReX_GpuQualifiers.H>
#include <AMReX_Periodicity.H>
#include <AMReX_Vector.H>

#include <AMReX_BaseFwd.H>

#include <optional>


namespace ablastr::utils::communication
{

using comm_float_type = float;

template <class FAB1, class FAB2>
void
mixedCopy (amrex::FabArray<FAB1>& dst, amrex::FabArray<FAB2> const& src, int srccomp, int dstcomp, int numcomp, const amrex::IntVect& nghost)
{
    auto const& srcma = src.const_arrays();
    auto const& dstma = dst.arrays();
    ParallelFor(dst, nghost, numcomp,
    [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
    {
        dstma[box_no](i,j,k,dstcomp+n) = (typename FAB1::value_type) srcma[box_no](i,j,k,srccomp+n);
    });
    amrex::Gpu::synchronize();
}

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::Periodicity::NonPeriodic(),
                  amrex::FabArrayBase::CpOp op = amrex::FabArrayBase::COPY);

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 = amrex::Periodicity::NonPeriodic());

void FillBoundary (amrex::MultiFab &mf,
                   bool do_single_precision_comms,
                   const amrex::Periodicity &period = amrex::Periodicity::NonPeriodic());

void FillBoundary (amrex::MultiFab &mf,
                   amrex::IntVect ng,
                   bool do_single_precision_comms,
                   const amrex::Periodicity &period = amrex::Periodicity::NonPeriodic(),
                   std::optional<bool> nodal_sync = std::nullopt);

void FillBoundary (amrex::iMultiFab &mf,
                   const amrex::Periodicity &period = amrex::Periodicity::NonPeriodic());

void FillBoundary (amrex::iMultiFab&         mf,
                   amrex::IntVect            ng,
                   const amrex::Periodicity& period = amrex::Periodicity::NonPeriodic());

void
FillBoundary(amrex::Vector<amrex::MultiFab *> const &mf, bool do_single_precision_comms,
             const amrex::Periodicity &period);

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 = amrex::Periodicity::NonPeriodic());

void OverrideSync (amrex::MultiFab &mf,
                   bool do_single_precision_comms,
                   const amrex::Periodicity &period = amrex::Periodicity::NonPeriodic());
}

#endif // ABLASTR_UTILS_COMMUNICATION_H_