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
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
|
/* Copyright 2019 Andrew Myers, Maxence Thevenet, Weiqun Zhang
*
*
* This file is part of WarpX.
*
* License: BSD-3-Clause-LBNL
*/
#include <WarpX.H>
#include <Filter.H>
#ifdef _OPENMP
#include <omp.h>
#endif
using namespace amrex;
#ifdef AMREX_USE_CUDA
/* \brief Apply stencil on MultiFab (GPU version, 2D/3D).
* \param dstmf Destination MultiFab
* \param srcmf source MultiFab
* \param scomp first component of srcmf on which the filter is applied
* \param dcomp first component of dstmf on which the filter is applied
* \param ncomp Number of components on which the filter is applied.
*/
void
Filter::ApplyStencil (MultiFab& dstmf, const MultiFab& srcmf, int scomp, int dcomp, int ncomp)
{
BL_PROFILE("BilinearFilter::ApplyStencil(MultiFab)");
ncomp = std::min(ncomp, srcmf.nComp());
for (MFIter mfi(dstmf); mfi.isValid(); ++mfi)
{
const auto& src = srcmf.array(mfi);
const auto& dst = dstmf.array(mfi);
const Box& tbx = mfi.growntilebox();
const Box& gbx = amrex::grow(tbx,stencil_length_each_dir-1);
// tmpfab has enough ghost cells for the stencil
FArrayBox tmp_fab(gbx,ncomp);
Elixir tmp_eli = tmp_fab.elixir(); // Prevent the tmp data from being deleted too early
auto const& tmp = tmp_fab.array();
// Copy values in srcfab into tmpfab
const Box& ibx = gbx & srcmf[mfi].box();
AMREX_PARALLEL_FOR_4D ( gbx, ncomp, i, j, k, n,
{
if (ibx.contains(IntVect(AMREX_D_DECL(i,j,k)))) {
tmp(i,j,k,n) = src(i,j,k,n+scomp);
} else {
tmp(i,j,k,n) = 0.0;
}
});
// Apply filter
DoFilter(tbx, tmp, dst, 0, dcomp, ncomp);
}
}
/* \brief Apply stencil on FArrayBox (GPU version, 2D/3D).
* \param dstfab Destination FArrayBox
* \param srcmf source FArrayBox
* \param tbx Grown box on which srcfab is defined.
* \param scomp first component of srcfab on which the filter is applied
* \param dcomp first component of dstfab on which the filter is applied
* \param ncomp Number of components on which the filter is applied.
*/
void
Filter::ApplyStencil (FArrayBox& dstfab, const FArrayBox& srcfab,
const Box& tbx, int scomp, int dcomp, int ncomp)
{
BL_PROFILE("BilinearFilter::ApplyStencil(FArrayBox)");
ncomp = std::min(ncomp, srcfab.nComp());
const auto& src = srcfab.array();
const auto& dst = dstfab.array();
const Box& gbx = amrex::grow(tbx,stencil_length_each_dir-1);
// tmpfab has enough ghost cells for the stencil
FArrayBox tmp_fab(gbx,ncomp);
Elixir tmp_eli = tmp_fab.elixir(); // Prevent the tmp data from being deleted too early
auto const& tmp = tmp_fab.array();
// Copy values in srcfab into tmpfab
const Box& ibx = gbx & srcfab.box();
AMREX_PARALLEL_FOR_4D ( gbx, ncomp, i, j, k, n,
{
if (ibx.contains(IntVect(AMREX_D_DECL(i,j,k)))) {
tmp(i,j,k,n) = src(i,j,k,n+scomp);
} else {
tmp(i,j,k,n) = 0.0;
}
});
// Apply filter
DoFilter(tbx, tmp, dst, 0, dcomp, ncomp);
}
/* \brief Apply stencil (2D/3D, CPU/GPU)
*/
void Filter::DoFilter (const Box& tbx,
Array4<Real const> const& tmp,
Array4<Real > const& dst,
int scomp, int dcomp, int ncomp)
{
amrex::Real const* AMREX_RESTRICT sx = stencil_x.data();
amrex::Real const* AMREX_RESTRICT sy = stencil_y.data();
amrex::Real const* AMREX_RESTRICT sz = stencil_z.data();
Dim3 slen_local = slen;
AMREX_PARALLEL_FOR_4D ( tbx, ncomp, i, j, k, n,
{
Real d = 0.0;
for (int iz=0; iz < slen_local.z; ++iz){
for (int iy=0; iy < slen_local.y; ++iy){
for (int ix=0; ix < slen_local.x; ++ix){
#if (AMREX_SPACEDIM == 3)
Real sss = sx[ix]*sy[iy]*sz[iz];
#else
Real sss = sx[ix]*sz[iy];
#endif
#if (AMREX_SPACEDIM == 3)
d += sss*( tmp(i-ix,j-iy,k-iz,scomp+n)
+tmp(i+ix,j-iy,k-iz,scomp+n)
+tmp(i-ix,j+iy,k-iz,scomp+n)
+tmp(i+ix,j+iy,k-iz,scomp+n)
+tmp(i-ix,j-iy,k+iz,scomp+n)
+tmp(i+ix,j-iy,k+iz,scomp+n)
+tmp(i-ix,j+iy,k+iz,scomp+n)
+tmp(i+ix,j+iy,k+iz,scomp+n));
#else
d += sss*( tmp(i-ix,j-iy,k,scomp+n)
+tmp(i+ix,j-iy,k,scomp+n)
+tmp(i-ix,j+iy,k,scomp+n)
+tmp(i+ix,j+iy,k,scomp+n));
#endif
}
}
}
dst(i,j,k,dcomp+n) = d;
});
}
#else
/* \brief Apply stencil on MultiFab (CPU version, 2D/3D).
* \param dstmf Destination MultiFab
* \param srcmf source MultiFab
* \param scomp first component of srcmf on which the filter is applied
* \param dcomp first component of dstmf on which the filter is applied
* \param ncomp Number of components on which the filter is applied.
*/
void
Filter::ApplyStencil (amrex::MultiFab& dstmf, const amrex::MultiFab& srcmf, int scomp, int dcomp, int ncomp)
{
BL_PROFILE("BilinearFilter::ApplyStencil()");
ncomp = std::min(ncomp, srcmf.nComp());
#ifdef _OPENMP
#pragma omp parallel
#endif
{
FArrayBox tmpfab;
for (MFIter mfi(dstmf,true); mfi.isValid(); ++mfi){
const auto& srcfab = srcmf[mfi];
auto& dstfab = dstmf[mfi];
const Box& tbx = mfi.growntilebox();
const Box& gbx = amrex::grow(tbx,stencil_length_each_dir-1);
// tmpfab has enough ghost cells for the stencil
tmpfab.resize(gbx,ncomp);
tmpfab.setVal(0.0, gbx, 0, ncomp);
// Copy values in srcfab into tmpfab
const Box& ibx = gbx & srcfab.box();
tmpfab.copy(srcfab, ibx, scomp, ibx, 0, ncomp);
// Apply filter
DoFilter(tbx, tmpfab.array(), dstfab.array(), 0, dcomp, ncomp);
}
}
}
/* \brief Apply stencil on FArrayBox (CPU version, 2D/3D).
* \param dstfab Destination FArrayBox
* \param srcmf source FArrayBox
* \param tbx Grown box on which srcfab is defined.
* \param scomp first component of srcfab on which the filter is applied
* \param dcomp first component of dstfab on which the filter is applied
* \param ncomp Number of components on which the filter is applied.
*/
void
Filter::ApplyStencil (amrex::FArrayBox& dstfab, const amrex::FArrayBox& srcfab,
const amrex::Box& tbx, int scomp, int dcomp, int ncomp)
{
BL_PROFILE("BilinearFilter::ApplyStencil(FArrayBox)");
ncomp = std::min(ncomp, srcfab.nComp());
FArrayBox tmpfab;
const Box& gbx = amrex::grow(tbx,stencil_length_each_dir-1);
// tmpfab has enough ghost cells for the stencil
tmpfab.resize(gbx,ncomp);
tmpfab.setVal(0.0, gbx, 0, ncomp);
// Copy values in srcfab into tmpfab
const Box& ibx = gbx & srcfab.box();
tmpfab.copy(srcfab, ibx, scomp, ibx, 0, ncomp);
// Apply filter
DoFilter(tbx, tmpfab.array(), dstfab.array(), 0, dcomp, ncomp);
}
void Filter::DoFilter (const Box& tbx,
Array4<Real const> const& tmp,
Array4<Real > const& dst,
int scomp, int dcomp, int ncomp)
{
const auto lo = amrex::lbound(tbx);
const auto hi = amrex::ubound(tbx);
// tmp and dst are of type Array4 (Fortran ordering)
amrex::Real const* AMREX_RESTRICT sx = stencil_x.data();
amrex::Real const* AMREX_RESTRICT sy = stencil_y.data();
amrex::Real const* AMREX_RESTRICT sz = stencil_z.data();
for (int n = 0; n < ncomp; ++n) {
// Set dst value to 0.
for (int k = lo.z; k <= hi.z; ++k) {
for (int j = lo.y; j <= hi.y; ++j) {
for (int i = lo.x; i <= hi.x; ++i) {
dst(i,j,k,dcomp+n) = 0.0;
}
}
}
// 3 nested loop on 3D stencil
for (int iz=0; iz < slen.z; ++iz){
for (int iy=0; iy < slen.y; ++iy){
for (int ix=0; ix < slen.x; ++ix){
#if (AMREX_SPACEDIM == 3)
Real sss = sx[ix]*sy[iy]*sz[iz];
#else
Real sss = sx[ix]*sz[iy];
#endif
// 3 nested loop on 3D array
for (int k = lo.z; k <= hi.z; ++k) {
for (int j = lo.y; j <= hi.y; ++j) {
AMREX_PRAGMA_SIMD
for (int i = lo.x; i <= hi.x; ++i) {
#if (AMREX_SPACEDIM == 3)
dst(i,j,k,dcomp+n) += sss*(tmp(i-ix,j-iy,k-iz,scomp+n)
+tmp(i+ix,j-iy,k-iz,scomp+n)
+tmp(i-ix,j+iy,k-iz,scomp+n)
+tmp(i+ix,j+iy,k-iz,scomp+n)
+tmp(i-ix,j-iy,k+iz,scomp+n)
+tmp(i+ix,j-iy,k+iz,scomp+n)
+tmp(i-ix,j+iy,k+iz,scomp+n)
+tmp(i+ix,j+iy,k+iz,scomp+n));
#else
dst(i,j,k,dcomp+n) += sss*(tmp(i-ix,j-iy,k,scomp+n)
+tmp(i+ix,j-iy,k,scomp+n)
+tmp(i-ix,j+iy,k,scomp+n)
+tmp(i+ix,j+iy,k,scomp+n));
#endif
}
}
}
}
}
}
}
}
#endif // #ifdef AMREX_USE_CUDA
|