Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_FabArrayUtility.H
Go to the documentation of this file.
1#ifndef AMREX_FABARRAY_UTILITY_H_
2#define AMREX_FABARRAY_UTILITY_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_Concepts.H>
6#include <AMReX_FabArray.H>
7#include <AMReX_LayoutData.H>
8#include <AMReX_Print.H>
9#include <AMReX_ParReduce.H>
10#include <limits>
11
12namespace amrex {
13
14template <BaseFabType FAB, class F>
15typename FAB::value_type
16ReduceSum (FabArray<FAB> const& fa, int nghost, F&& f)
17{
18 return ReduceSum(fa, IntVect(nghost), std::forward<F>(f));
19}
20
22namespace fudetail {
23template <BaseFabType FAB, class F>
24typename FAB::value_type
25ReduceSum_host (FabArray<FAB> const& fa, IntVect const& nghost, F const& f)
26{
27 using value_type = typename FAB::value_type;
28 value_type sm = 0;
29
30#ifdef AMREX_USE_OMP
31#pragma omp parallel if (!system::regtest_reduction) reduction(+:sm)
32#endif
33 for (MFIter mfi(fa,true); mfi.isValid(); ++mfi)
34 {
35 const Box& bx = mfi.growntilebox(nghost);
36 auto const& arr = fa.const_array(mfi);
37 sm += f(bx, arr);
38 }
39
40 return sm;
41}
42}
44
45#ifdef AMREX_USE_GPU
47namespace fudetail {
48template <class OP, BaseFabType FAB, class F>
49std::conditional_t<std::is_same_v<OP,ReduceOpLogicalAnd> ||
50 std::is_same_v<OP,ReduceOpLogicalOr>,
51 int, typename FAB::value_type>
52ReduceMF (FabArray<FAB> const& fa, IntVect const& nghost, F const& f)
53{
54 using T = std::conditional_t<std::is_same_v<OP,ReduceOpLogicalAnd> ||
55 std::is_same_v<OP,ReduceOpLogicalOr>,
56 int, typename FAB::value_type>;
57 auto typ = fa.ixType();
58 auto const& ma = fa.const_arrays();
59 return ParReduce(TypeList<OP>{}, TypeList<T>{}, fa, nghost,
60 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept -> GpuTuple<T>
61 {
62 return { static_cast<T>(f(amrex::makeSingleCellBox(i,j,k,typ), ma[box_no])) };
63 });
64}
65
66template <class OP, BaseFabType FAB1, BaseFabType FAB2, class F>
67std::conditional_t<std::is_same_v<OP,ReduceOpLogicalAnd> ||
68 std::is_same_v<OP,ReduceOpLogicalOr>,
69 int, typename FAB1::value_type>
70ReduceMF (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2, IntVect const& nghost, F const& f)
71{
72 using T = std::conditional_t<std::is_same_v<OP,ReduceOpLogicalAnd> ||
73 std::is_same_v<OP,ReduceOpLogicalOr>,
74 int, typename FAB1::value_type>;
75 auto typ = fa1.ixType();
76 auto const& ma1 = fa1.const_arrays();
77 auto const& ma2 = fa2.const_arrays();
78 return ParReduce(TypeList<OP>{}, TypeList<T>{}, fa1, nghost,
79 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept -> GpuTuple<T>
80 {
81 return { static_cast<T>(f(amrex::makeSingleCellBox(i,j,k,typ),
82 ma1[box_no], ma2[box_no])) };
83 });
84}
85
86template <class OP, BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
87std::conditional_t<std::is_same_v<OP,ReduceOpLogicalAnd> ||
88 std::is_same_v<OP,ReduceOpLogicalOr>,
89 int, typename FAB1::value_type>
90ReduceMF (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
91 FabArray<FAB3> const& fa3, IntVect const& nghost, F const& f)
92{
93 using T = std::conditional_t<std::is_same_v<OP,ReduceOpLogicalAnd> ||
94 std::is_same_v<OP,ReduceOpLogicalOr>,
95 int, typename FAB1::value_type>;
96 auto typ = fa1.ixType();
97 auto const& ma1 = fa1.const_arrays();
98 auto const& ma2 = fa2.const_arrays();
99 auto const& ma3 = fa3.const_arrays();
100 return ParReduce(TypeList<OP>{}, TypeList<T>{}, fa1, nghost,
101 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept -> GpuTuple<T>
102 {
103 return { static_cast<T>(f(amrex::makeSingleCellBox(i,j,k,typ),
104 ma1[box_no], ma2[box_no], ma3[box_no])) };
105 });
106}
107
108template <BaseFabType FAB, class F>
110typename FAB::value_type ReduceSum_host_wrapper (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
111{
112 return ReduceSum_host(fa,nghost,std::forward<F>(f));
113}
114
115template <BaseFabType FAB, class F>
117typename FAB::value_type ReduceSum_host_wrapper (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
118{
119 amrex::ignore_unused(fa,nghost,f);
120 amrex::Abort("ReduceSum: Launch Region is off. Device lambda cannot be called by host.");
121 return 0;
122}
123}
125
126template <BaseFabType FAB, class F>
127typename FAB::value_type
128ReduceSum (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
129{
130 if (Gpu::inLaunchRegion()) {
131 return fudetail::ReduceMF<ReduceOpSum>(fa, nghost, std::forward<F>(f));
132 } else {
133 return fudetail::ReduceSum_host_wrapper(fa, nghost, std::forward<F>(f));
134 }
135}
136#else
137template <BaseFabType FAB, class F>
138typename FAB::value_type
139ReduceSum (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
140{
141 return fudetail::ReduceSum_host(fa, nghost, std::forward<F>(f));
142}
143#endif
144
145template <BaseFabType FAB1, BaseFabType FAB2, class F>
146typename FAB1::value_type
148 int nghost, F&& f)
149{
150 return ReduceSum(fa1, fa2, IntVect(nghost), std::forward<F>(f));
151}
152
154namespace fudetail {
155template <BaseFabType FAB1, BaseFabType FAB2, class F>
156typename FAB1::value_type
157ReduceSum_host (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
158 IntVect const& nghost, F const& f)
159{
160 using value_type = typename FAB1::value_type;
161 value_type sm = 0;
162
163#ifdef AMREX_USE_OMP
164#pragma omp parallel if (!system::regtest_reduction) reduction(+:sm)
165#endif
166 for (MFIter mfi(fa1,true); mfi.isValid(); ++mfi)
167 {
168 const Box& bx = mfi.growntilebox(nghost);
169 const auto& arr1 = fa1.const_array(mfi);
170 const auto& arr2 = fa2.const_array(mfi);
171 sm += f(bx, arr1, arr2);
172 }
173
174 return sm;
175}
176}
178
179#ifdef AMREX_USE_GPU
181namespace fudetail {
182template <BaseFabType FAB1, BaseFabType FAB2, class F>
184typename FAB1::value_type ReduceSum_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
185 IntVect const& nghost, F&& f)
186{
187 return ReduceSum_host(fa1,fa2,nghost,std::forward<F>(f));
188}
189
190template <BaseFabType FAB1, BaseFabType FAB2, class F>
192typename FAB1::value_type ReduceSum_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
193 IntVect const& nghost, F&& f)
194{
195 amrex::ignore_unused(fa1,fa2,nghost,f);
196 amrex::Abort("ReduceSum: Launch Region is off. Device lambda cannot be called by host.");
197 return 0;
198}
199}
201
202template <BaseFabType FAB1, BaseFabType FAB2, class F>
203typename FAB1::value_type
205 IntVect const& nghost, F&& f)
206{
207 if (Gpu::inLaunchRegion()) {
208 return fudetail::ReduceMF<ReduceOpSum>(fa1,fa2,nghost,std::forward<F>(f));
209 } else {
210 return fudetail::ReduceSum_host_wrapper(fa1,fa2,nghost, std::forward<F>(f));
211 }
212}
213#else
214template <BaseFabType FAB1, BaseFabType FAB2, class F>
215typename FAB1::value_type
216ReduceSum (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
217 IntVect const& nghost, F&& f)
218{
219 return fudetail::ReduceSum_host(fa1,fa2,nghost,std::forward<F>(f));
220}
221#endif
222
223template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
224typename FAB1::value_type
225ReduceSum (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2, FabArray<FAB3> const& fa3,
226 int nghost, F&& f)
227{
228 return ReduceSum(fa1, fa2, fa3, IntVect(nghost), std::forward<F>(f));
229}
230
232namespace fudetail {
233template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
234typename FAB1::value_type
235ReduceSum_host (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
236 FabArray<FAB3> const& fa3, IntVect const& nghost, F const& f)
237{
238 using value_type = typename FAB1::value_type;
239 value_type sm = 0;
240
241#ifdef AMREX_USE_OMP
242#pragma omp parallel if (!system::regtest_reduction) reduction(+:sm)
243#endif
244 for (MFIter mfi(fa1,true); mfi.isValid(); ++mfi)
245 {
246 const Box& bx = mfi.growntilebox(nghost);
247 const auto& arr1 = fa1.const_array(mfi);
248 const auto& arr2 = fa2.const_array(mfi);
249 const auto& arr3 = fa3.const_array(mfi);
250 sm += f(bx, arr1, arr2, arr3);
251 }
252
253 return sm;
254}
255}
257
258#ifdef AMREX_USE_GPU
260namespace fudetail {
261template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
263typename FAB1::value_type ReduceSum_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
264 FabArray<FAB3> const& fa3, IntVect const& nghost, F&& f)
265{
266 return fudetail::ReduceSum_host(fa1,fa2,fa3,nghost,std::forward<F>(f));
267}
268
269template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
271typename FAB1::value_type ReduceSum_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
272 FabArray<FAB3> const& fa3, IntVect const& nghost, F&& f)
273{
274 amrex::ignore_unused(fa1,fa2,fa3,nghost,f);
275 amrex::Abort("ReduceSum: Launch Region is off. Device lambda cannot be called by host.");
276 return 0;
277}
278}
280
281template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
282typename FAB1::value_type
284 FabArray<FAB3> const& fa3, IntVect const& nghost, F&& f)
285{
286 if (Gpu::inLaunchRegion()) {
287 return fudetail::ReduceMF<ReduceOpSum>(fa1,fa2,fa3,nghost,std::forward<F>(f));
288 } else {
289 return fudetail::ReduceSum_host_wrapper(fa1,fa2,fa3,nghost,std::forward<F>(f));
290 }
291}
292#else
293template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
294typename FAB1::value_type
295ReduceSum (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
296 FabArray<FAB3> const& fa3, IntVect const& nghost, F&& f)
297{
298 return fudetail::ReduceSum_host(fa1,fa2,fa3,nghost,std::forward<F>(f));
299}
300#endif
301
302template <BaseFabType FAB, class F>
303typename FAB::value_type
304ReduceMin (FabArray<FAB> const& fa, int nghost, F&& f)
305{
306 return ReduceMin(fa, IntVect(nghost), std::forward<F>(f));
307}
308
310namespace fudetail {
311template <BaseFabType FAB, class F>
312typename FAB::value_type
313ReduceMin_host (FabArray<FAB> const& fa, IntVect const& nghost, F const& f)
314{
315 using value_type = typename FAB::value_type;
316 value_type r = std::numeric_limits<value_type>::max();
317
318#ifdef AMREX_USE_OMP
319#pragma omp parallel reduction(min:r)
320#endif
321 for (MFIter mfi(fa,true); mfi.isValid(); ++mfi)
322 {
323 const Box& bx = mfi.growntilebox(nghost);
324 const auto& arr = fa.const_array(mfi);
325 r = std::min(r, f(bx, arr));
326 }
327 return r;
328}
329}
331
332#ifdef AMREX_USE_GPU
334namespace fudetail {
335template <BaseFabType FAB, class F>
337typename FAB::value_type ReduceMin_host_wrapper (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
338{
339 return ReduceMin_host(fa,nghost,std::forward<F>(f));
340}
341
342template <BaseFabType FAB, class F>
344typename FAB::value_type ReduceMin_host_wrapper (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
345{
346 amrex::ignore_unused(fa,nghost,f);
347 amrex::Abort("ReduceMin: Launch Region is off. Device lambda cannot be called by host.");
348 return 0;
349}
350}
352
353template <BaseFabType FAB, class F>
354typename FAB::value_type
355ReduceMin (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
356{
357 if (Gpu::inLaunchRegion()) {
358 return fudetail::ReduceMF<ReduceOpMin>(fa, nghost, std::forward<F>(f));
359 } else {
360 return fudetail::ReduceMin_host_wrapper(fa, nghost, std::forward<F>(f));
361 }
362}
363#else
364template <BaseFabType FAB, class F>
365typename FAB::value_type
366ReduceMin (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
367{
368 return fudetail::ReduceMin_host(fa, nghost, std::forward<F>(f));
369}
370#endif
371
372template <BaseFabType FAB1, BaseFabType FAB2, class F>
373typename FAB1::value_type
374ReduceMin (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2, int nghost, F&& f)
375{
376 return ReduceMin(fa1, fa2, IntVect(nghost), std::forward<F>(f));
377}
378
380namespace fudetail {
381template <BaseFabType FAB1, BaseFabType FAB2, class F>
382typename FAB1::value_type
383ReduceMin_host (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
384 IntVect const& nghost, F const& f)
385{
386 using value_type = typename FAB1::value_type;
387 value_type r = std::numeric_limits<value_type>::max();
388
389#ifdef AMREX_USE_OMP
390#pragma omp parallel reduction(min:r)
391#endif
392 for (MFIter mfi(fa1,true); mfi.isValid(); ++mfi)
393 {
394 const Box& bx = mfi.growntilebox(nghost);
395 const auto& arr1 = fa1.const_array(mfi);
396 const auto& arr2 = fa2.const_array(mfi);
397 r = std::min(r, f(bx, arr1, arr2));
398 }
399
400 return r;
401}
402}
404
405#ifdef AMREX_USE_GPU
407namespace fudetail {
408template <BaseFabType FAB1, BaseFabType FAB2, class F>
410typename FAB1::value_type ReduceMin_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
411 IntVect const& nghost, F&& f)
412{
413 return fudetail::ReduceMin_host(fa1,fa2,nghost,std::forward<F>(f));
414}
415
416template <BaseFabType FAB1, BaseFabType FAB2, class F>
418typename FAB1::value_type ReduceMin_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
419 IntVect const& nghost, F&& f)
420{
421 amrex::ignore_unused(fa1,fa2,nghost,f);
422 amrex::Abort("ReduceMin: Launch Region is off. Device lambda cannot be called by host.");
423 return 0;
424}
425}
427
428template <BaseFabType FAB1, BaseFabType FAB2, class F>
429typename FAB1::value_type
431 IntVect const& nghost, F&& f)
432{
433 if (Gpu::inLaunchRegion()) {
434 return fudetail::ReduceMF<ReduceOpMin>(fa1,fa2,nghost,std::forward<F>(f));
435 } else {
436 return fudetail::ReduceMin_host_wrapper(fa1,fa2,nghost,std::forward<F>(f));
437 }
438}
439#else
440template <BaseFabType FAB1, BaseFabType FAB2, class F>
441typename FAB1::value_type
442ReduceMin (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
443 IntVect const& nghost, F&& f)
444{
445 return fudetail::ReduceMin_host(fa1,fa2,nghost,std::forward<F>(f));
446}
447#endif
448
449template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
450typename FAB1::value_type
451ReduceMin (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2, FabArray<FAB3> const& fa3,
452 int nghost, F&& f)
453{
454 return ReduceMin(fa1, fa2, fa3, IntVect(nghost), std::forward<F>(f));
455}
456
458namespace fudetail {
459template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
460typename FAB1::value_type
461ReduceMin_host (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
462 FabArray<FAB3> const& fa3, IntVect const& nghost, F const& f)
463{
464 using value_type = typename FAB1::value_type;
465 value_type r = std::numeric_limits<value_type>::max();
466
467#ifdef AMREX_USE_OMP
468#pragma omp parallel reduction(min:r)
469#endif
470 for (MFIter mfi(fa1,true); mfi.isValid(); ++mfi)
471 {
472 const Box& bx = mfi.growntilebox(nghost);
473 const auto& arr1 = fa1.const_array(mfi);
474 const auto& arr2 = fa2.const_array(mfi);
475 const auto& arr3 = fa3.const_array(mfi);
476 r = std::min(r, f(bx, arr1, arr2, arr3));
477 }
478
479 return r;
480}
481}
483
484#ifdef AMREX_USE_GPU
486namespace fudetail {
487template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
489typename FAB1::value_type ReduceMin_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
490 FabArray<FAB3> const& fa3, IntVect const& nghost, F&& f)
491{
492 return fudetail::ReduceMin_host(fa1,fa2,fa3,nghost,std::forward<F>(f));
493}
494
495template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
497typename FAB1::value_type ReduceMin_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
498 FabArray<FAB3> const& fa3, IntVect const& nghost, F&& f)
499{
500 amrex::ignore_unused(fa1,fa2,fa3,nghost,f);
501 amrex::Abort("ReduceMin: Launch Region is off. Device lambda cannot be called by host.");
502 return 0;
503}
504}
506
507template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
508typename FAB1::value_type
510 FabArray<FAB3> const& fa3, IntVect const& nghost, F&& f)
511{
512 if (Gpu::inLaunchRegion()) {
513 return fudetail::ReduceMF<ReduceOpMin>(fa1,fa2,fa3,nghost,std::forward<F>(f));
514 } else {
515 return fudetail::ReduceMin_host_wrapper(fa1,fa2,fa3,nghost,std::forward<F>(f));
516 }
517}
518#else
519template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
520typename FAB1::value_type
521ReduceMin (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
522 FabArray<FAB3> const& fa3, IntVect const& nghost, F&& f)
523{
524 return fudetail::ReduceMin_host(fa1,fa2,fa3,nghost,std::forward<F>(f));
525}
526#endif
527
528template <BaseFabType FAB, class F>
529typename FAB::value_type
530ReduceMax (FabArray<FAB> const& fa, int nghost, F&& f)
531{
532 return ReduceMax(fa, IntVect(nghost), std::forward<F>(f));
533}
534
536namespace fudetail {
537template <BaseFabType FAB, class F>
538typename FAB::value_type
539ReduceMax_host (FabArray<FAB> const& fa, IntVect const& nghost, F const& f)
540{
541 using value_type = typename FAB::value_type;
542 value_type r = std::numeric_limits<value_type>::lowest();
543
544#ifdef AMREX_USE_OMP
545#pragma omp parallel reduction(max:r)
546#endif
547 for (MFIter mfi(fa,true); mfi.isValid(); ++mfi)
548 {
549 const Box& bx = mfi.growntilebox(nghost);
550 const auto& arr = fa.const_array(mfi);
551 r = std::max(r, f(bx, arr));
552 }
553
554 return r;
555}
556}
558
559#ifdef AMREX_USE_GPU
561namespace fudetail {
562template <BaseFabType FAB, class F>
564typename FAB::value_type ReduceMax_host_wrapper (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
565{
566 return ReduceMax_host(fa,nghost,std::forward<F>(f));
567}
568
569template <BaseFabType FAB, class F>
571typename FAB::value_type ReduceMax_host_wrapper (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
572{
573 amrex::ignore_unused(fa,nghost,f);
574 amrex::Abort("ReduceMax: Launch Region is off. Device lambda cannot be called by host.");
575 return 0;
576}
577}
579
580template <BaseFabType FAB, class F>
581typename FAB::value_type
582ReduceMax (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
583{
584 if (Gpu::inLaunchRegion()) {
585 return fudetail::ReduceMF<ReduceOpMax>(fa,nghost,std::forward<F>(f));
586 } else {
587 return fudetail::ReduceMax_host_wrapper(fa,nghost,std::forward<F>(f));
588 }
589}
590#else
591template <BaseFabType FAB, class F>
592typename FAB::value_type
593ReduceMax (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
594{
595 return fudetail::ReduceMax_host(fa,nghost,std::forward<F>(f));
596}
597#endif
598
599template <BaseFabType FAB1, BaseFabType FAB2, class F>
600typename FAB1::value_type
601ReduceMax (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2, int nghost, F&& f)
602{
603 return ReduceMax(fa1, fa2, IntVect(nghost), std::forward<F>(f));
604}
605
607namespace fudetail {
608template <BaseFabType FAB1, BaseFabType FAB2, class F>
609typename FAB1::value_type
610ReduceMax_host (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
611 IntVect const& nghost, F const& f)
612{
613 using value_type = typename FAB1::value_type;
614 value_type r = std::numeric_limits<value_type>::lowest();
615
616#ifdef AMREX_USE_OMP
617#pragma omp parallel reduction(max:r)
618#endif
619 for (MFIter mfi(fa1,true); mfi.isValid(); ++mfi)
620 {
621 const Box& bx = mfi.growntilebox(nghost);
622 const auto& arr1 = fa1.const_array(mfi);
623 const auto& arr2 = fa2.const_array(mfi);
624 r = std::max(r, f(bx, arr1, arr2));
625 }
626
627 return r;
628}
629}
631
632#ifdef AMREX_USE_GPU
634namespace fudetail {
635template <BaseFabType FAB1, BaseFabType FAB2, class F>
637typename FAB1::value_type ReduceMax_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
638 IntVect const& nghost, F&& f)
639{
640 return ReduceMax_host(fa1,fa2,nghost,std::forward<F>(f));
641}
642
643template <BaseFabType FAB1, BaseFabType FAB2, class F>
645typename FAB1::value_type ReduceMax_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
646 IntVect const& nghost, F&& f)
647{
648 amrex::ignore_unused(fa1,fa2,nghost,f);
649 amrex::Abort("ReduceMax: Launch Region is off. Device lambda cannot be called by host.");
650 return 0;
651}
652}
654
655template <BaseFabType FAB1, BaseFabType FAB2, class F>
656typename FAB1::value_type
658 IntVect const& nghost, F&& f)
659{
660 if (Gpu::inLaunchRegion()) {
661 return fudetail::ReduceMF<ReduceOpMax>(fa1,fa2,nghost,std::forward<F>(f));
662 } else {
663 return fudetail::ReduceMax_host_wrapper(fa1,fa2,nghost,std::forward<F>(f));
664 }
665}
666#else
667template <BaseFabType FAB1, BaseFabType FAB2, class F>
668typename FAB1::value_type
669ReduceMax (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
670 IntVect const& nghost, F&& f)
671{
672 return fudetail::ReduceMax_host(fa1,fa2,nghost,std::forward<F>(f));
673}
674#endif
675
676template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
677typename FAB1::value_type
678ReduceMax (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2, FabArray<FAB3> const& fa3,
679 int nghost, F&& f)
680{
681 return ReduceMax(fa1, fa2, fa3, IntVect(nghost), std::forward<F>(f));
682}
683
685namespace fudetail {
686template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
687typename FAB1::value_type
688ReduceMax_host (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
689 FabArray<FAB3> const& fa3, IntVect const& nghost, F const& f)
690{
691 using value_type = typename FAB1::value_type;
692 value_type r = std::numeric_limits<value_type>::lowest();
693
694#ifdef AMREX_USE_OMP
695#pragma omp parallel reduction(max:r)
696#endif
697 for (MFIter mfi(fa1,true); mfi.isValid(); ++mfi)
698 {
699 const Box& bx = mfi.growntilebox(nghost);
700 const auto& arr1 = fa1.const_array(mfi);
701 const auto& arr2 = fa2.const_array(mfi);
702 const auto& arr3 = fa3.const_array(mfi);
703 r = std::max(r, f(bx, arr1, arr2, arr3));
704 }
705
706 return r;
707}
708}
710
711#ifdef AMREX_USE_GPU
713namespace fudetail {
714template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
716typename FAB1::value_type ReduceMax_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
717 FabArray<FAB3> const& fa3, IntVect const& nghost, F&& f)
718{
719 return fudetail::ReduceMax_host(fa1,fa2,fa3,nghost,std::forward<F>(f));
720}
721
722template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
724typename FAB1::value_type ReduceMax_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
725 FabArray<FAB3> const& fa3, IntVect const& nghost, F&& f)
726{
727 amrex::ignore_unused(fa1,fa2,fa3,nghost,f);
728 amrex::Abort("ReduceMax: Launch Region is off. Device lambda cannot be called by host.");
729 return 0;
730}
731}
733
734template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
735typename FAB1::value_type
737 FabArray<FAB3> const& fa3, IntVect const& nghost, F&& f)
738{
739 if (Gpu::inLaunchRegion()) {
740 return fudetail::ReduceMF<ReduceOpMax>(fa1,fa2,fa3,nghost,std::forward<F>(f));
741 } else {
742 return fudetail::ReduceMax_host_wrapper(fa1,fa2,fa3,nghost,std::forward<F>(f));
743 }
744}
745#else
746template <BaseFabType FAB1, BaseFabType FAB2, BaseFabType FAB3, class F>
747typename FAB1::value_type
748ReduceMax (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
749 FabArray<FAB3> const& fa3, IntVect const& nghost, F&& f)
750{
751 return fudetail::ReduceMax_host(fa1,fa2,fa3,nghost,std::forward<F>(f));
752}
753#endif
754
755template <BaseFabType FAB, class F>
756bool
757ReduceLogicalAnd (FabArray<FAB> const& fa, int nghost, F&& f)
758{
759 return ReduceLogicalAnd(fa, IntVect(nghost), std::forward<F>(f));
760}
761
763namespace fudetail {
764template <BaseFabType FAB, class F>
765bool
766ReduceLogicalAnd_host (FabArray<FAB> const& fa, IntVect const& nghost, F const& f)
767{
768 int r = true;
769
770#ifdef AMREX_USE_OMP
771#pragma omp parallel reduction(&&:r)
772#endif
773 for (MFIter mfi(fa,true); mfi.isValid(); ++mfi)
774 {
775 const Box& bx = mfi.growntilebox(nghost);
776 const auto& arr = fa.const_array(mfi);
777 r = r && f(bx, arr);
778 }
779
780 return r;
781}
782}
784
785#ifdef AMREX_USE_GPU
787namespace fudetail {
788template <BaseFabType FAB, class F>
790bool ReduceLogicalAnd_host_wrapper (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
791{
792 return ReduceLogicalAnd_host(fa,nghost,std::forward<F>(f));
793}
794
795template <BaseFabType FAB, class F>
797bool ReduceLogicalAnd_host_wrapper (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
798{
799 amrex::ignore_unused(fa,nghost,f);
800 amrex::Abort("ReduceLogicalAnd: Launch Region is off. Device lambda cannot be called by host.");
801 return false;
802}
803}
805
806template <BaseFabType FAB, class F>
807bool
808ReduceLogicalAnd (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
809{
810 if (Gpu::inLaunchRegion()) {
811 return fudetail::ReduceMF<ReduceOpLogicalAnd>(fa,nghost,std::forward<F>(f));
812 } else {
813 return fudetail::ReduceLogicalAnd_host_wrapper(fa,nghost,std::forward<F>(f));
814 }
815}
816#else
817template <BaseFabType FAB, class F>
818bool
819ReduceLogicalAnd (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
820{
821 return fudetail::ReduceLogicalAnd_host(fa,nghost,std::forward<F>(f));
822}
823#endif
824
825template <BaseFabType FAB1, BaseFabType FAB2, class F>
826bool
828 int nghost, F&& f)
829{
830 return ReduceLogicalAnd(fa1, fa2, IntVect(nghost), std::forward<F>(f));
831}
832
834namespace fudetail {
835template <BaseFabType FAB1, BaseFabType FAB2, class F>
836bool
837ReduceLogicalAnd_host (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
838 IntVect const& nghost, F const& f)
839{
840 int r = true;
841
842#ifdef AMREX_USE_OMP
843#pragma omp parallel reduction(&&:r)
844#endif
845 for (MFIter mfi(fa1,true); mfi.isValid(); ++mfi)
846 {
847 const Box& bx = mfi.growntilebox(nghost);
848 const auto& arr1 = fa1.const_array(mfi);
849 const auto& arr2 = fa2.const_array(mfi);
850 r = r && f(bx, arr1, arr2);
851 }
852
853 return r;
854}
855}
857
858#ifdef AMREX_USE_GPU
860namespace fudetail {
861template <BaseFabType FAB1, BaseFabType FAB2, class F>
863bool ReduceLogicalAnd_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
864 IntVect const& nghost, F&& f)
865{
866 return ReduceLogicalAnd_host(fa1,fa2,nghost,std::forward<F>(f));
867}
868
869template <BaseFabType FAB1, BaseFabType FAB2, class F>
871bool ReduceLogicalAnd_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
872 IntVect const& nghost, F&& f)
873{
874 amrex::ignore_unused(fa1,fa2,nghost,f);
875 amrex::Abort("ReduceLogicalAnd: Luanch Region is off. Device lambda cannot be called by host.");
876 return false;
877}
878}
880
881template <BaseFabType FAB1, BaseFabType FAB2, class F>
882bool
884 IntVect const& nghost, F&& f)
885{
886 if (Gpu::inLaunchRegion()) {
887 return fudetail::ReduceMF<ReduceOpLogicalAnd>(fa1,fa2,nghost,std::forward<F>(f));
888 } else {
889 return fudetail::ReduceLogicalAnd_host_wrapper(fa1,fa2,nghost,std::forward<F>(f));
890 }
891}
892#else
893template <BaseFabType FAB1, BaseFabType FAB2, class F>
894bool
895ReduceLogicalAnd (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
896 IntVect const& nghost, F&& f)
897{
898 return fudetail::ReduceLogicalAnd_host(fa1,fa2,nghost,std::forward<F>(f));
899}
900#endif
901
902template <BaseFabType FAB, class F>
903bool
904ReduceLogicalOr (FabArray<FAB> const& fa, int nghost, F&& f)
905{
906 return ReduceLogicalOr(fa, IntVect(nghost), std::forward<F>(f));
907}
908
910namespace fudetail {
911template <BaseFabType FAB, class F>
912bool
913ReduceLogicalOr_host (FabArray<FAB> const& fa, IntVect const& nghost, F const& f)
914{
915 int r = false;
916
917#ifdef AMREX_USE_OMP
918#pragma omp parallel reduction(||:r)
919#endif
920 for (MFIter mfi(fa,true); mfi.isValid(); ++mfi)
921 {
922 const Box& bx = mfi.growntilebox(nghost);
923 const auto& arr = fa.const_array(mfi);
924 r = r || f(bx, arr);
925 }
926
927 return r;
928}
929}
931
932#ifdef AMREX_USE_GPU
934namespace fudetail {
935template <BaseFabType FAB, class F>
937bool ReduceLogicalOr_host_wrapper (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
938{
939 return ReduceLogicalOr_host(fa,nghost,std::forward<F>(f));
940}
941
942template <BaseFabType FAB, class F>
944bool ReduceLogicalOr_host (FabArray<FAB> const& fa, IntVect const& nghost, F&& /*f*/)
945{
946 amrex::ignore_unused(fa,nghost);
947 amrex::Abort("ReduceLogicalOr: Launch Region is off. Device lambda cannot be called by host.");
948 return 0;
949}
950}
952
953template <BaseFabType FAB, class F>
954bool
955ReduceLogicalOr (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
956{
957 if (Gpu::inLaunchRegion()) {
958 return fudetail::ReduceMF<ReduceOpLogicalOr>(fa,nghost,std::forward<F>(f));
959 } else {
960 return fudetail::ReduceLogicalOr_host_wrapper(fa,nghost,std::forward<F>(f));
961 }
962}
963#else
964template <BaseFabType FAB, class F>
965bool
966ReduceLogicalOr (FabArray<FAB> const& fa, IntVect const& nghost, F&& f)
967{
968 return fudetail::ReduceLogicalOr_host(fa,nghost,std::forward<F>(f));
969}
970#endif
971
972template <BaseFabType FAB1, BaseFabType FAB2, class F>
973bool
975 int nghost, F&& f)
976{
977 return ReduceLogicalOr(fa1, fa2, IntVect(nghost), std::forward<F>(f));
978}
979
981namespace fudetail {
982template <BaseFabType FAB1, BaseFabType FAB2, class F>
983bool
984ReduceLogicalOr_host (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
985 IntVect const& nghost, F const& f)
986{
987 int r = false;
988
989#ifdef AMREX_USE_OMP
990#pragma omp parallel reduction(||:r)
991#endif
992 for (MFIter mfi(fa1,true); mfi.isValid(); ++mfi)
993 {
994 const Box& bx = mfi.growntilebox(nghost);
995 const auto& arr1 = fa1.const_array(mfi);
996 const auto& arr2 = fa2.const_array(mfi);
997 r = r || f(bx, arr1, arr2);
998 }
999
1000 return r;
1001}
1002}
1004
1005#ifdef AMREX_USE_GPU
1007namespace fudetail {
1008template <BaseFabType FAB1, BaseFabType FAB2, class F>
1010bool ReduceLogicalOr_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
1011 IntVect const& nghost, F&& f)
1012{
1013 return fudetail::ReduceLogicalOr_host(fa1,fa2,nghost,std::forward<F>(f));
1014}
1015
1016template <BaseFabType FAB1, BaseFabType FAB2, class F>
1018bool ReduceLogicalOr_host_wrapper (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
1019 IntVect const& nghost, F&& f)
1020{
1021 amrex::ignore_unused(fa1,fa2,nghost,f);
1022 amrex::Abort("ReeuceLogicalOr: Launch Region is off. Device lambda cannot be called by host.");
1023 return false;
1024}
1025}
1027
1028template <BaseFabType FAB1, BaseFabType FAB2, class F>
1029bool
1031 IntVect const& nghost, F&& f)
1032{
1033 if (Gpu::inLaunchRegion()) {
1034 return fudetail::ReduceMF<ReduceOpLogicalOr>(fa1,fa2,nghost,std::forward<F>(f));
1035 } else {
1036 return fudetail::ReduceLogicalOr_host_wrapper(fa1,fa2,nghost,std::forward<F>(f));
1037 }
1038}
1039#else
1040template <BaseFabType FAB1, BaseFabType FAB2, class F>
1041bool
1042ReduceLogicalOr (FabArray<FAB1> const& fa1, FabArray<FAB2> const& fa2,
1043 IntVect const& nghost, F&& f)
1044{
1045 return fudetail::ReduceLogicalOr_host(fa1,fa2,nghost,std::forward<F>(f));
1046}
1047#endif
1048
1049template <BaseFabType FAB>
1050void
1051printCell (FabArray<FAB> const& mf, const IntVect& cell, int comp = -1,
1052 const IntVect& ng = IntVect::TheZeroVector())
1053{
1054 for (MFIter mfi(mf); mfi.isValid(); ++mfi)
1055 {
1056 const Box& bx = amrex::grow(mfi.validbox(), ng);
1057 if (bx.contains(cell)) {
1058 int n = (comp >= 0) ? 1 : mf.nComp();
1059 auto const& fab = mf.const_array(mfi);
1061 auto* dp = pv.data();
1062 auto f = [=] AMREX_GPU_HOST_DEVICE ()
1063 {
1064 if (comp >= 0) {
1065 *dp = fab(cell, comp);
1066 } else {
1067 for (int i = 0; i < n; ++i) {
1068 dp[i] = fab(cell,i);
1069 }
1070 }
1071 };
1072
1073#ifdef AMREX_USE_GPU
1074 if (mf.arena()->isManaged() || mf.arena()->isDevice()) {
1077 } else
1078#endif
1079 {
1080 f();
1081 }
1082
1083 if (comp >= 0) {
1084 amrex::AllPrint().SetPrecision(17) << " At cell " << cell << " in Box " << bx
1085 << ": " << *dp << '\n';
1086 } else {
1087 std::ostringstream ss;
1088 ss.precision(17);
1089 for (int i = 0; i < n-1; ++i)
1090 {
1091 ss << dp[i] << ", ";
1092 }
1093 ss << dp[n-1];
1094 amrex::AllPrint() << " At cell " << cell << " in Box " << bx
1095 << ": " << ss.view() << '\n';
1096 }
1097 }
1098 }
1099}
1100
1101template <BaseFabType FAB>
1102void
1103Swap (FabArray<FAB>& dst, FabArray<FAB>& src, int srccomp, int dstcomp, int numcomp, int nghost)
1104{
1105 Swap(dst,src,srccomp,dstcomp,numcomp,IntVect(nghost));
1106}
1107
1108template <BaseFabType FAB>
1109void
1110Swap (FabArray<FAB>& dst, FabArray<FAB>& src, int srccomp, int dstcomp, int numcomp, const IntVect& nghost)
1111{
1112 // We can take a shortcut and do a std::swap if we're swapping all of the data
1113 // and they are allocated in the same Arena.
1114
1115 bool explicit_swap = true;
1116
1117 if (srccomp == dstcomp && dstcomp == 0 && src.nComp() == dst.nComp() &&
1118 src.nGrowVect() == nghost && src.nGrowVect() == dst.nGrowVect() &&
1119 src.arena() == dst.arena() && src.hasEBFabFactory() == dst.hasEBFabFactory()) {
1120 explicit_swap = false;
1121 }
1122
1123 if (!explicit_swap) {
1124
1125 std::swap(dst, src);
1126
1127 } else {
1128#ifdef AMREX_USE_GPU
1129 if (Gpu::inLaunchRegion() && dst.isFusingCandidate()) {
1130 auto const& dstma = dst.arrays();
1131 auto const& srcma = src.arrays();
1132 ParallelFor(dst, nghost, numcomp,
1133 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1134 {
1135 const auto tmp = dstma[box_no](i,j,k,n+dstcomp);
1136 dstma[box_no](i,j,k,n+dstcomp) = srcma[box_no](i,j,k,n+srccomp);
1137 srcma[box_no](i,j,k,n+srccomp) = tmp;
1138 });
1139 if (!Gpu::inNoSyncRegion()) {
1141 }
1142 } else
1143#endif
1144 {
1145#ifdef AMREX_USE_OMP
1146#pragma omp parallel if (Gpu::notInLaunchRegion())
1147#endif
1148 for (MFIter mfi(dst,TilingIfNotGPU()); mfi.isValid(); ++mfi)
1149 {
1150 const Box& bx = mfi.growntilebox(nghost);
1151 if (bx.ok()) {
1152 auto sfab = src.array(mfi);
1153 auto dfab = dst.array(mfi);
1154 AMREX_HOST_DEVICE_PARALLEL_FOR_4D( bx, numcomp, i, j, k, n,
1155 {
1156 const auto tmp = dfab(i,j,k,n+dstcomp);
1157 dfab(i,j,k,n+dstcomp) = sfab(i,j,k,n+srccomp);
1158 sfab(i,j,k,n+srccomp) = tmp;
1159 });
1160 }
1161 }
1162 }
1163 }
1164}
1165
1166template <BaseFabType FAB>
1167void
1168Subtract (FabArray<FAB>& dst, FabArray<FAB> const& src, int srccomp, int dstcomp, int numcomp, int nghost)
1169{
1170 Subtract(dst,src,srccomp,dstcomp,numcomp,IntVect(nghost));
1171}
1172
1173template <BaseFabType FAB>
1174void
1175Subtract (FabArray<FAB>& dst, FabArray<FAB> const& src, int srccomp, int dstcomp, int numcomp, const IntVect& nghost)
1176{
1177#ifdef AMREX_USE_GPU
1178 if (Gpu::inLaunchRegion() && dst.isFusingCandidate()) {
1179 auto const& dstfa = dst.arrays();
1180 auto const& srcfa = src.const_arrays();
1181 ParallelFor(dst, nghost, numcomp,
1182 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1183 {
1184 dstfa[box_no](i,j,k,n+dstcomp) -= srcfa[box_no](i,j,k,n+srccomp);
1185 });
1186 if (!Gpu::inNoSyncRegion()) {
1188 }
1189 } else
1190#endif
1191 {
1192#ifdef AMREX_USE_OMP
1193#pragma omp parallel if (Gpu::notInLaunchRegion())
1194#endif
1195 for (MFIter mfi(dst,TilingIfNotGPU()); mfi.isValid(); ++mfi)
1196 {
1197 const Box& bx = mfi.growntilebox(nghost);
1198 if (bx.ok())
1199 {
1200 auto const srcFab = src.array(mfi);
1201 auto dstFab = dst.array(mfi);
1202 AMREX_HOST_DEVICE_PARALLEL_FOR_4D( bx, numcomp, i, j, k, n,
1203 {
1204 dstFab(i,j,k,n+dstcomp) -= srcFab(i,j,k,n+srccomp);
1205 });
1206 }
1207 }
1208 }
1209}
1210
1211
1212template <BaseFabType FAB>
1213void
1214Multiply (FabArray<FAB>& dst, FabArray<FAB> const& src, int srccomp, int dstcomp, int numcomp, int nghost)
1215{
1216 Multiply(dst,src,srccomp,dstcomp,numcomp,IntVect(nghost));
1217}
1218
1219template <BaseFabType FAB>
1220void
1221Multiply (FabArray<FAB>& dst, FabArray<FAB> const& src, int srccomp, int dstcomp, int numcomp, const IntVect& nghost)
1222{
1223#ifdef AMREX_USE_GPU
1224 if (Gpu::inLaunchRegion() && dst.isFusingCandidate()) {
1225 auto const& dstfa = dst.arrays();
1226 auto const& srcfa = src.const_arrays();
1227 ParallelFor(dst, nghost, numcomp,
1228 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1229 {
1230 dstfa[box_no](i,j,k,n+dstcomp) *= srcfa[box_no](i,j,k,n+srccomp);
1231 });
1232 if (!Gpu::inNoSyncRegion()) {
1234 }
1235 } else
1236#endif
1237 {
1238#ifdef AMREX_USE_OMP
1239#pragma omp parallel if (Gpu::notInLaunchRegion())
1240#endif
1241 for (MFIter mfi(dst,TilingIfNotGPU()); mfi.isValid(); ++mfi)
1242 {
1243 const Box& bx = mfi.growntilebox(nghost);
1244 if (bx.ok())
1245 {
1246 auto const srcFab = src.array(mfi);
1247 auto dstFab = dst.array(mfi);
1248 AMREX_HOST_DEVICE_PARALLEL_FOR_4D( bx, numcomp, i, j, k, n,
1249 {
1250 dstFab(i,j,k,n+dstcomp) *= srcFab(i,j,k,n+srccomp);
1251 });
1252 }
1253 }
1254 }
1255}
1256
1257
1258template <BaseFabType FAB>
1259void
1260Divide (FabArray<FAB>& dst, FabArray<FAB> const& src, int srccomp, int dstcomp, int numcomp, int nghost)
1261{
1262 Divide(dst,src,srccomp,dstcomp,numcomp,IntVect(nghost));
1263}
1264
1265template <BaseFabType FAB>
1266void
1267Divide (FabArray<FAB>& dst, FabArray<FAB> const& src, int srccomp, int dstcomp, int numcomp, const IntVect& nghost)
1268{
1269#ifdef AMREX_USE_GPU
1270 if (Gpu::inLaunchRegion() && dst.isFusingCandidate()) {
1271 auto const& dstfa = dst.arrays();
1272 auto const& srcfa = src.const_arrays();
1273 ParallelFor(dst, nghost, numcomp,
1274 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1275 {
1276 dstfa[box_no](i,j,k,n+dstcomp) /= srcfa[box_no](i,j,k,n+srccomp);
1277 });
1278 if (!Gpu::inNoSyncRegion()) {
1280 }
1281 } else
1282#endif
1283 {
1284#ifdef AMREX_USE_OMP
1285#pragma omp parallel if (Gpu::notInLaunchRegion())
1286#endif
1287 for (MFIter mfi(dst,TilingIfNotGPU()); mfi.isValid(); ++mfi)
1288 {
1289 const Box& bx = mfi.growntilebox(nghost);
1290 if (bx.ok())
1291 {
1292 auto const srcFab = src.array(mfi);
1293 auto dstFab = dst.array(mfi);
1294 AMREX_HOST_DEVICE_PARALLEL_FOR_4D( bx, numcomp, i, j, k, n,
1295 {
1296 dstFab(i,j,k,n+dstcomp) /= srcFab(i,j,k,n+srccomp);
1297 });
1298 }
1299 }
1300 }
1301}
1302
1303template <BaseFabType FAB>
1304void
1305Abs (FabArray<FAB>& fa, int icomp, int numcomp, int nghost)
1306{
1307 Abs(fa,icomp,numcomp,IntVect(nghost));
1308}
1309
1310template <BaseFabType FAB>
1311void
1312Abs (FabArray<FAB>& fa, int icomp, int numcomp, const IntVect& nghost)
1313{
1314#ifdef AMREX_USE_GPU
1315 if (Gpu::inLaunchRegion() && fa.isFusingCandidate()) {
1316 auto const& fabarr = fa.arrays();
1317 ParallelFor(fa, nghost, numcomp,
1318 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1319 {
1320 fabarr[box_no](i,j,k,n+icomp) = std::abs(fabarr[box_no](i,j,k,n+icomp));
1321 });
1322 if (!Gpu::inNoSyncRegion()) {
1324 }
1325 } else
1326#endif
1327 {
1328#ifdef AMREX_USE_OMP
1329#pragma omp parallel if (Gpu::notInLaunchRegion())
1330#endif
1331 for (MFIter mfi(fa,TilingIfNotGPU()); mfi.isValid(); ++mfi)
1332 {
1333 const Box& bx = mfi.growntilebox(nghost);
1334 if (bx.ok())
1335 {
1336 auto const& fab = fa.array(mfi);
1337 AMREX_HOST_DEVICE_PARALLEL_FOR_4D( bx, numcomp, i, j, k, n,
1338 {
1339 fab(i,j,k,n+icomp) = std::abs(fab(i,j,k,n+icomp));
1340 });
1341 }
1342 }
1343 }
1344}
1345
1346template <BaseFabType FAB>
1347void
1348prefetchToHost (FabArray<FAB> const& fa, const bool synchronous = true)
1349{
1350#ifdef AMREX_USE_GPU
1351 if (fa.arena()->isManaged()) {
1352 for (MFIter mfi(fa, MFItInfo().SetDeviceSync(synchronous)); mfi.isValid(); ++mfi) {
1353 fa.prefetchToHost(mfi);
1354 }
1355 }
1356#else
1357 amrex::ignore_unused(fa,synchronous);
1358#endif
1359}
1360
1361template <BaseFabType FAB>
1362void
1363prefetchToDevice (FabArray<FAB> const& fa, const bool synchronous = true)
1364{
1365#ifdef AMREX_USE_GPU
1366 if (fa.arena()->isManaged()) {
1367 for (MFIter mfi(fa, MFItInfo().SetDeviceSync(synchronous)); mfi.isValid(); ++mfi) {
1368 fa.prefetchToDevice(mfi);
1369 }
1370 }
1371#else
1372 amrex::ignore_unused(fa,synchronous);
1373#endif
1374}
1375
1376
1377template <BaseFabType FAB, BaseFabType IFAB>
1378void
1379OverrideSync (FabArray<FAB> & fa, FabArray<IFAB> const& msk, const Periodicity& period)
1380{
1381 BL_PROFILE("OverrideSync()");
1382
1383 OverrideSync_nowait(fa, msk, period);
1385}
1386
1387
1388template <BaseFabType FAB, BaseFabType IFAB>
1389void
1391{
1392 BL_PROFILE("OverrideSync_nowait()");
1393 AMREX_ASSERT_WITH_MESSAGE(!fa.os_temp, "OverrideSync_nowait() called when already in progress.");
1394
1395 if (fa.ixType().cellCentered()) { return; }
1396
1397 const int ncomp = fa.nComp();
1398
1399#ifdef AMREX_USE_GPU
1400 if (Gpu::inLaunchRegion() && fa.isFusingCandidate()) {
1401 auto const& fabarr = fa.arrays();
1402 auto const& ifabarr = msk.const_arrays();
1403 ParallelFor(fa, IntVect(0), ncomp,
1404 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k, int n) noexcept
1405 {
1406 if (!ifabarr[box_no](i,j,k)) { fabarr[box_no](i,j,k,n) = 0; }
1407 });
1408 if (!Gpu::inNoSyncRegion()) {
1410 }
1411 } else
1412#endif
1413 {
1414#ifdef AMREX_USE_OMP
1415#pragma omp parallel if (Gpu::notInLaunchRegion())
1416#endif
1417 for (MFIter mfi(fa,TilingIfNotGPU()); mfi.isValid(); ++mfi)
1418 {
1419 const Box& bx = mfi.tilebox();
1420 auto fab = fa.array(mfi);
1421 auto const ifab = msk.array(mfi);
1422 AMREX_HOST_DEVICE_PARALLEL_FOR_4D( bx, ncomp, i, j, k, n,
1423 {
1424 if (!ifab(i,j,k)) { fab(i,j,k,n) = 0; }
1425 });
1426 }
1427 }
1428
1429 fa.os_temp = std::make_unique< FabArray<FAB> > ( fa.boxArray(), fa.DistributionMap(),
1430 ncomp, 0, MFInfo(), fa.Factory() );
1431 fa.os_temp->setVal(0);
1432 fa.os_temp->ParallelCopy_nowait(fa, period, FabArrayBase::ADD);
1433}
1434
1435template <BaseFabType FAB>
1436void
1438{
1439 BL_PROFILE("OverrideSync_finish()");
1440
1441 if (fa.ixType().cellCentered()) { return; }
1442
1443 fa.os_temp->ParallelCopy_finish();
1444 amrex::Copy(fa, *(fa.os_temp), 0, 0, fa.nComp(), 0);
1445
1446 fa.os_temp.reset();
1447}
1448
1449template <BaseFabType FAB>
1450void
1452 int scomp, int dcomp, int ncomp)
1453{
1454 AMREX_ASSERT(isMFIterSafe(dst, src));
1455 AMREX_ASSERT(dst.nGrowVect() == src.nGrowVect());
1456#ifdef AMREX_USE_GPU
1457 for (MFIter mfi(dst); mfi.isValid(); ++mfi) {
1458 void* pdst = dst[mfi].dataPtr(dcomp);
1459 void const* psrc = src[mfi].dataPtr(scomp);
1460 Gpu::dtoh_memcpy_async(pdst, psrc, dst[mfi].nBytes(mfi.fabbox(), ncomp));
1461 }
1462#else
1463 Copy(dst, src, scomp, dcomp, ncomp, dst.nGrowVect());
1464#endif
1465}
1466
1467template <BaseFabType FAB>
1468void
1470{
1471 dtoh_memcpy(dst, src, 0, 0, dst.nComp());
1472}
1473
1474template <BaseFabType FAB>
1475void
1477 int scomp, int dcomp, int ncomp)
1478{
1479 AMREX_ASSERT(isMFIterSafe(dst, src));
1480 AMREX_ASSERT(dst.nGrowVect() == src.nGrowVect());
1481#ifdef AMREX_USE_GPU
1482 for (MFIter mfi(dst); mfi.isValid(); ++mfi) {
1483 void* pdst = dst[mfi].dataPtr(dcomp);
1484 void const* psrc = src[mfi].dataPtr(scomp);
1485 Gpu::htod_memcpy_async(pdst, psrc, dst[mfi].nBytes(mfi.fabbox(), ncomp));
1486 }
1487#else
1488 Copy(dst, src, scomp, dcomp, ncomp, dst.nGrowVect());
1489#endif
1490}
1491
1492template <BaseFabType FAB>
1493void
1495{
1496 htod_memcpy(dst, src, 0, 0, dst.nComp());
1497}
1498
1499template <BaseFabType FAB>
1500IntVect
1501indexFromValue (FabArray<FAB> const& mf, int comp, IntVect const& nghost,
1502 typename FAB::value_type value)
1503{
1504 IntVect loc;
1505
1506#ifdef AMREX_USE_GPU
1507 if (Gpu::inLaunchRegion())
1508 {
1510 int* p = aa.data();
1511 // This is a device ptr to 1+AMREX_SPACEDIM int zeros.
1512 // The first is used as an atomic bool and the others for intvect.
1513 if (mf.isFusingCandidate()) {
1514 auto const& ma = mf.const_arrays();
1515 ParallelFor(mf, nghost, [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept
1516 {
1517 int* flag = p;
1518 if (*flag == 0) {
1519 if (ma[box_no](i,j,k,comp) == value) {
1520 if (Gpu::Atomic::Exch(flag,1) == 0) {
1521 AMREX_D_TERM(p[1] = i;,
1522 p[2] = j;,
1523 p[3] = k;);
1524 }
1525 }
1526 }
1527 });
1528 } else {
1529 for (MFIter mfi(mf,MFItInfo().SetDeviceSync(false)); mfi.isValid(); ++mfi) {
1530 const Box& bx = amrex::grow(mfi.validbox(), nghost);
1531 auto const& arr = mf.const_array(mfi);
1532 amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
1533 {
1534 int* flag = p;
1535 if (*flag == 0) {
1536 if (arr(i,j,k,comp) == value) {
1537 if (Gpu::Atomic::Exch(flag,1) == 0) {
1538 AMREX_D_TERM(p[1] = i;,
1539 p[2] = j;,
1540 p[3] = k;);
1541 }
1542 }
1543 }
1544 });
1545 }
1546 }
1547 int const* tmp = aa.copyToHost();
1548 AMREX_D_TERM(loc[0] = tmp[1];,
1549 loc[1] = tmp[2];,
1550 loc[2] = tmp[3];);
1551 }
1552 else
1553#endif
1554 {
1555 bool f = false;
1556#ifdef AMREX_USE_OMP
1557#pragma omp parallel
1558#endif
1559 {
1560 IntVect priv_loc = IntVect::TheMinVector();
1561 for (MFIter mfi(mf,true); mfi.isValid(); ++mfi)
1562 {
1563 const Box& bx = mfi.growntilebox(nghost);
1564 auto const& fab = mf.const_array(mfi);
1565 AMREX_LOOP_3D(bx, i, j, k,
1566 {
1567 if (fab(i,j,k,comp) == value) {
1568 priv_loc = IntVect(AMREX_D_DECL(i,j,k));
1569 }
1570 });
1571 }
1572
1573 if (priv_loc.allGT(IntVect::TheMinVector())) {
1574 bool old;
1575// we should be able to test on _OPENMP < 201107 for capture (version 3.1)
1576// but we must work around a bug in gcc < 4.9
1577// And, with NVHPC 21.9 to <23.1, we saw an ICE with the atomic capture (NV bug: #3390723)
1578#if defined(AMREX_USE_OMP) && defined(_OPENMP) && (_OPENMP < 201307 || (defined(__NVCOMPILER) && __NVCOMPILER_MAJOR__ < 23)) // OpenMP 4.0
1579#pragma omp critical (amrex_indexfromvalue)
1580#elif defined(AMREX_USE_OMP)
1581#pragma omp atomic capture
1582#endif
1583 {
1584 old = f;
1585 f = true;
1586 }
1587
1588 if (old == false) { loc = priv_loc; }
1589 }
1590 }
1591 }
1592
1593 return loc;
1594}
1595
1607template <BaseFabType FAB>
1608typename FAB::value_type
1609Dot (FabArray<FAB> const& x, int xcomp, FabArray<FAB> const& y, int ycomp, int ncomp,
1610 IntVect const& nghost, bool local = false)
1611{
1612 BL_ASSERT(x.boxArray() == y.boxArray());
1613 BL_ASSERT(x.DistributionMap() == y.DistributionMap());
1614 BL_ASSERT(x.nGrowVect().allGE(nghost) && y.nGrowVect().allGE(nghost));
1615
1616 BL_PROFILE("amrex::Dot()");
1617
1618 using T = typename FAB::value_type;
1619 auto sm = T(0.0);
1620#ifdef AMREX_USE_GPU
1621 if (Gpu::inLaunchRegion()) {
1622 auto const& xma = x.const_arrays();
1623 auto const& yma = y.const_arrays();
1624 sm = ParReduce(TypeList<ReduceOpSum>{}, TypeList<T>{}, x, nghost,
1625 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept -> GpuTuple<T>
1626 {
1627 auto t = T(0.0);
1628 auto const& xfab = xma[box_no];
1629 auto const& yfab = yma[box_no];
1630 for (int n = 0; n < ncomp; ++n) {
1631 t += xfab(i,j,k,xcomp+n) * yfab(i,j,k,ycomp+n);
1632 }
1633 return t;
1634 });
1635 } else
1636#endif
1637 {
1638#ifdef AMREX_USE_OMP
1639#pragma omp parallel if (!system::regtest_reduction) reduction(+:sm)
1640#endif
1641 for (MFIter mfi(x,true); mfi.isValid(); ++mfi)
1642 {
1643 Box const& bx = mfi.growntilebox(nghost);
1644 auto const& xfab = x.const_array(mfi);
1645 auto const& yfab = y.const_array(mfi);
1646 AMREX_LOOP_4D(bx, ncomp, i, j, k, n,
1647 {
1648 sm += xfab(i,j,k,xcomp+n) * yfab(i,j,k,ycomp+n);
1649 });
1650 }
1651 }
1652
1653 if (!local) {
1655 }
1656
1657 return sm;
1658}
1659
1669template <BaseFabType FAB>
1670typename FAB::value_type
1671Dot (FabArray<FAB> const& x, int xcomp, int ncomp, IntVect const& nghost, bool local = false)
1672{
1673 BL_ASSERT(x.nGrowVect().allGE(nghost));
1674
1675 BL_PROFILE("amrex::Dot()");
1676
1677 using T = typename FAB::value_type;
1678 auto sm = T(0.0);
1679#ifdef AMREX_USE_GPU
1680 if (Gpu::inLaunchRegion()) {
1681 auto const& xma = x.const_arrays();
1682 sm = ParReduce(TypeList<ReduceOpSum>{}, TypeList<T>{}, x, nghost,
1683 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept -> GpuTuple<T>
1684 {
1685 auto t = T(0.0);
1686 auto const& xfab = xma[box_no];
1687 for (int n = 0; n < ncomp; ++n) {
1688 auto v = xfab(i,j,k,xcomp+n);
1689 t += v*v;
1690 }
1691 return t;
1692 });
1693 } else
1694#endif
1695 {
1696#ifdef AMREX_USE_OMP
1697#pragma omp parallel if (!system::regtest_reduction) reduction(+:sm)
1698#endif
1699 for (MFIter mfi(x,true); mfi.isValid(); ++mfi)
1700 {
1701 Box const& bx = mfi.growntilebox(nghost);
1702 auto const& xfab = x.const_array(mfi);
1703 AMREX_LOOP_4D(bx, ncomp, i, j, k, n,
1704 {
1705 auto v = xfab(i,j,k,xcomp+n);
1706 sm += v*v;
1707 });
1708 }
1709 }
1710
1711 if (!local) {
1713 }
1714
1715 return sm;
1716}
1717
1730template <BaseFabType IFAB, BaseFabType FAB>
1731typename FAB::value_type
1732Dot (FabArray<IFAB> const& mask, FabArray<FAB> const& x, int xcomp,
1733 FabArray<FAB> const& y, int ycomp, int ncomp, IntVect const& nghost,
1734 bool local = false)
1735{
1736 BL_ASSERT(x.boxArray() == y.boxArray() && x.boxArray() == mask.boxArray());
1737 BL_ASSERT(x.DistributionMap() == y.DistributionMap() && x.DistributionMap() == mask.DistributionMap());
1738 BL_ASSERT(x.nGrowVect().allGE(nghost) && y.nGrowVect().allGE(nghost) &&
1739 mask.nGrowVect().allGE(nghost));
1740
1741 BL_PROFILE("amrex::Dot()");
1742
1743 using T = typename FAB::value_type;
1744 auto sm = T(0.0);
1745#ifdef AMREX_USE_GPU
1746 if (Gpu::inLaunchRegion()) {
1747 auto const& mma = mask.const_arrays();
1748 auto const& xma = x.const_arrays();
1749 auto const& yma = y.const_arrays();
1750 sm = ParReduce(TypeList<ReduceOpSum>{}, TypeList<T>{}, x, nghost,
1751 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept -> GpuTuple<T>
1752 {
1753 auto t = T(0.0);
1754 auto m = T(mma[box_no](i,j,k));
1755 if (m != 0) {
1756 auto const& xfab = xma[box_no];
1757 auto const& yfab = yma[box_no];
1758 for (int n = 0; n < ncomp; ++n) {
1759 t += xfab(i,j,k,xcomp+n) * yfab(i,j,k,ycomp+n);
1760 }
1761 }
1762 return t*m;
1763 });
1764 } else
1765#endif
1766 {
1767#ifdef AMREX_USE_OMP
1768#pragma omp parallel if (!system::regtest_reduction) reduction(+:sm)
1769#endif
1770 for (MFIter mfi(x,true); mfi.isValid(); ++mfi)
1771 {
1772 Box const& bx = mfi.growntilebox(nghost);
1773 auto const& mfab = mask.const_array(mfi);
1774 auto const& xfab = x.const_array(mfi);
1775 auto const& yfab = y.const_array(mfi);
1776 AMREX_LOOP_4D(bx, ncomp, i, j, k, n,
1777 {
1778 auto m = T(mfab(i,j,k));
1779 sm += m * xfab(i,j,k,xcomp+n) * yfab(i,j,k,ycomp+n);
1780 });
1781 }
1782 }
1783
1784 if (!local) {
1786 }
1787
1788 return sm;
1789}
1790
1801template <BaseFabType IFAB, BaseFabType FAB>
1802typename FAB::value_type
1803Dot (FabArray<IFAB> const& mask, FabArray<FAB> const& x, int xcomp, int ncomp,
1804 IntVect const& nghost, bool local = false)
1805{
1806 BL_ASSERT(x.boxArray() == mask.boxArray());
1807 BL_ASSERT(x.DistributionMap() == mask.DistributionMap());
1808 BL_ASSERT(x.nGrowVect().allGE(nghost) && mask.nGrowVect().allGE(nghost));
1809
1810 BL_PROFILE("amrex::Dot()");
1811
1812 using T = typename FAB::value_type;
1813 auto sm = T(0.0);
1814#ifdef AMREX_USE_GPU
1815 if (Gpu::inLaunchRegion()) {
1816 auto const& mma = mask.const_arrays();
1817 auto const& xma = x.const_arrays();
1818 sm = ParReduce(TypeList<ReduceOpSum>{}, TypeList<T>{}, x, nghost,
1819 [=] AMREX_GPU_DEVICE (int box_no, int i, int j, int k) noexcept -> GpuTuple<T>
1820 {
1821 auto t = T(0.0);
1822 auto m = T(mma[box_no](i,j,k));
1823 if (m != 0) {
1824 auto const& xfab = xma[box_no];
1825 for (int n = 0; n < ncomp; ++n) {
1826 auto v = xfab(i,j,k,xcomp+n);
1827 t += v*v;
1828 }
1829 }
1830 return t*m;
1831 });
1832 } else
1833#endif
1834 {
1835#ifdef AMREX_USE_OMP
1836#pragma omp parallel if (!system::regtest_reduction) reduction(+:sm)
1837#endif
1838 for (MFIter mfi(x,true); mfi.isValid(); ++mfi)
1839 {
1840 Box const& bx = mfi.growntilebox(nghost);
1841 auto const& mfab = mask.const_array(mfi);
1842 auto const& xfab = x.const_array(mfi);
1843 AMREX_LOOP_4D(bx, ncomp, i, j, k, n,
1844 {
1845 auto m = T(mfab(i,j,k));
1846 auto v = xfab(i,j,k,xcomp+n);
1847 sm += m*v*v;
1848 });
1849 }
1850 }
1851
1852 if (!local) {
1854 }
1855
1856 return sm;
1857}
1858
1860template <MultiFabLike MF>
1861void setVal (MF& dst, typename MF::value_type val)
1862{
1863 dst.setVal(val);
1864}
1865
1867template <MultiFabLike MF>
1868void setBndry (MF& dst, typename MF::value_type val, int scomp, int ncomp)
1869{
1870 dst.setBndry(val, scomp, ncomp);
1871}
1872
1874template <MultiFabLike MF>
1875void Scale (MF& dst, typename MF::value_type val, int scomp, int ncomp, int nghost)
1876{
1877 dst.mult(val, scomp, ncomp, nghost);
1878}
1879
1881template <MultiFabLike DMF, MultiFabLike SMF>
1882void LocalCopy (DMF& dst, SMF const& src, int scomp, int dcomp,
1883 int ncomp, IntVect const& nghost)
1884{
1885 amrex::Copy(dst, src, scomp, dcomp, ncomp, nghost);
1886}
1887
1889template <MultiFabLike MF>
1890void LocalAdd (MF& dst, MF const& src, int scomp, int dcomp,
1891 int ncomp, IntVect const& nghost)
1892{
1893 amrex::Add(dst, src, scomp, dcomp, ncomp, nghost);
1894}
1895
1897template <MultiFabLike MF>
1898void Saxpy (MF& dst, typename MF::value_type a, MF const& src, int scomp, int dcomp,
1899 int ncomp, IntVect const& nghost)
1900{
1901 MF::Saxpy(dst, a, src, scomp, dcomp, ncomp, nghost);
1902}
1903
1905template <MultiFabLike MF>
1906void Xpay (MF& dst, typename MF::value_type a, MF const& src, int scomp, int dcomp,
1907 int ncomp, IntVect const& nghost)
1908{
1909 MF::Xpay(dst, a, src, scomp, dcomp, ncomp, nghost);
1910}
1911
1913template <MultiFabLike MF>
1914void Saxpy_Xpay (MF& dst, typename MF::value_type a_saxpy, MF const& src_saxpy,
1915 typename MF::value_type a_xpay, MF const& src_xpay, int scomp, int dcomp,
1916 int ncomp, IntVect const& nghost)
1917{
1918 MF::Saxpy_Xpay(dst, a_saxpy, src_saxpy, a_xpay, src_xpay, scomp, dcomp, ncomp, nghost);
1919}
1920
1922template <MultiFabLike MF>
1923void Saxpy_Saxpy (MF& dst1, typename MF::value_type a1, MF const& src1,
1924 MF& dst2, typename MF::value_type a2, MF const& src2, int scomp, int dcomp,
1925 int ncomp, IntVect const& nghost)
1926{
1927 MF::Saxpy_Saxpy(dst1, a1, src1, dst2, a2, src2, scomp, dcomp, ncomp, nghost);
1928}
1929
1931template <MultiFabLike MF>
1932void Saypy_Saxpy (MF& dst1, typename MF::value_type a1,
1933 MF& dst2, typename MF::value_type a2, MF const& src, int scomp, int dcomp,
1934 int ncomp, IntVect const& nghost)
1935{
1936 MF::Saypy_Saxpy(dst1, a1, dst2, a2, src, scomp, dcomp, ncomp, nghost);
1937}
1938
1940template <MultiFabLike MF>
1941void LinComb (MF& dst,
1942 typename MF::value_type a, MF const& src_a, int acomp,
1943 typename MF::value_type b, MF const& src_b, int bcomp,
1944 int dcomp, int ncomp, IntVect const& nghost)
1945{
1946 MF::LinComb(dst, a, src_a, acomp, b, src_b, bcomp, dcomp, ncomp, nghost);
1947}
1948
1950template <MultiFabLike MF>
1951void ParallelCopy (MF& dst, MF const& src, int scomp, int dcomp, int ncomp,
1952 IntVect const& ng_src = IntVect(0),
1953 IntVect const& ng_dst = IntVect(0),
1954 Periodicity const& period = Periodicity::NonPeriodic())
1955{
1956 dst.ParallelCopy(src, scomp, dcomp, ncomp, ng_src, ng_dst, period);
1957}
1958
1959template <MultiFabLike MF>
1960[[nodiscard]] typename MF::value_type
1961norminf (MF const& mf, int scomp, int ncomp, IntVect const& nghost,
1962 bool local = false)
1963{
1964 return mf.norminf(scomp, ncomp, nghost, local);
1965}
1966
1968template <MultiFabLike MF, std::size_t N>
1969void setVal (Array<MF,N>& dst, typename MF::value_type val)
1970{
1971 for (auto& mf: dst) {
1972 mf.setVal(val);
1973 }
1974}
1975
1977template <MultiFabLike MF, std::size_t N>
1978void setBndry (Array<MF,N>& dst, typename MF::value_type val, int scomp, int ncomp)
1979{
1980 for (auto& mf : dst) {
1981 mf.setBndry(val, scomp, ncomp);
1982 }
1983}
1984
1986template <MultiFabLike MF, std::size_t N>
1987void Scale (Array<MF,N>& dst, typename MF::value_type val, int scomp, int ncomp,
1988 int nghost)
1989{
1990 for (auto& mf : dst) {
1991 mf.mult(val, scomp, ncomp, nghost);
1992 }
1993}
1994
1996template <MultiFabLike DMF, MultiFabLike SMF, std::size_t N>
1997void LocalCopy (Array<DMF,N>& dst, Array<SMF,N> const& src, int scomp, int dcomp,
1998 int ncomp, IntVect const& nghost)
1999{
2000 for (std::size_t i = 0; i < N; ++i) {
2001 amrex::Copy(dst[i], src[i], scomp, dcomp, ncomp, nghost);
2002 }
2003}
2004
2006template <MultiFabLike MF, std::size_t N>
2007void LocalAdd (Array<MF,N>& dst, Array<MF,N> const& src, int scomp, int dcomp,
2008 int ncomp, IntVect const& nghost)
2009{
2010 for (std::size_t i = 0; i < N; ++i) {
2011 amrex::Add(dst[i], src[i], scomp, dcomp, ncomp, nghost);
2012 }
2013}
2014
2016template <MultiFabLike MF, std::size_t N>
2017void Saxpy (Array<MF,N>& dst, typename MF::value_type a,
2018 Array<MF,N> const& src, int scomp, int dcomp, int ncomp,
2019 IntVect const& nghost)
2020{
2021 for (std::size_t i = 0; i < N; ++i) {
2022 MF::Saxpy(dst[i], a, src[i], scomp, dcomp, ncomp, nghost);
2023 }
2024}
2025
2027template <MultiFabLike MF, std::size_t N>
2028void Xpay (Array<MF,N>& dst, typename MF::value_type a,
2029 Array<MF,N> const& src, int scomp, int dcomp, int ncomp,
2030 IntVect const& nghost)
2031{
2032 for (std::size_t i = 0; i < N; ++i) {
2033 MF::Xpay(dst[i], a, src[i], scomp, dcomp, ncomp, nghost);
2034 }
2035}
2036
2038template <MultiFabLike MF, std::size_t N>
2040 typename MF::value_type a, Array<MF,N> const& src_a, int acomp,
2041 typename MF::value_type b, Array<MF,N> const& src_b, int bcomp,
2042 int dcomp, int ncomp, IntVect const& nghost)
2043{
2044 for (std::size_t i = 0; i < N; ++i) {
2045 MF::LinComb(dst[i], a, src_a[i], acomp, b, src_b[i], bcomp, dcomp, ncomp, nghost);
2046 }
2047}
2048
2050template <MultiFabLike MF, std::size_t N>
2051void ParallelCopy (Array<MF,N>& dst, Array<MF,N> const& src,
2052 int scomp, int dcomp, int ncomp,
2053 IntVect const& ng_src = IntVect(0),
2054 IntVect const& ng_dst = IntVect(0),
2055 Periodicity const& period = Periodicity::NonPeriodic())
2056{
2057 for (std::size_t i = 0; i < N; ++i) {
2058 dst[i].ParallelCopy(src[i], scomp, dcomp, ncomp, ng_src, ng_dst, period);
2059 }
2060}
2061
2062template <MultiFabLike MF, std::size_t N>
2063[[nodiscard]] typename MF::value_type
2064norminf (Array<MF,N> const& mf, int scomp, int ncomp, IntVect const& nghost,
2065 bool local = false)
2066{
2067 auto r = typename MF::value_type(0);
2068 for (std::size_t i = 0; i < N; ++i) {
2069 auto tmp = mf[i].norminf(scomp, ncomp, nghost, true);
2070 r = std::max(r,tmp);
2071 }
2072 if (!local) {
2074 }
2075 return r;
2076}
2077
2078template <MultiFabLike MF, std::size_t N>
2079requires (N > 0)
2080[[nodiscard]] int nComp (Array<MF,N> const& mf)
2081{
2082 return mf[0].nComp();
2083}
2084
2085template <MultiFabLike MF, std::size_t N>
2086requires (N > 0)
2087[[nodiscard]] IntVect nGrowVect (Array<MF,N> const& mf)
2088{
2089 return mf[0].nGrowVect();
2090}
2091
2092template <MultiFabLike MF, std::size_t N>
2093requires (N > 0)
2094[[nodiscard]] BoxArray const&
2096{
2097 return mf[0].boxArray();
2098}
2099
2100template <MultiFabLike MF, std::size_t N>
2101requires (N > 0)
2102[[nodiscard]] DistributionMapping const&
2104{
2105 return mf[0].DistributionMap();
2106}
2107
2108/*
2109 * \brief Return a mask indicating how many duplicates are in each point
2110 *
2111 * \param fa input FabArray
2112 * \param nghost number of ghost cells included in counting
2113 * \param period periodicity
2114 */
2115template <class FAB>
2116FabArray<BaseFab<int>>
2117OverlapMask (FabArray<FAB> const& fa, IntVect const& nghost, Periodicity const& period)
2118{
2119 BL_PROFILE("OverlapMask()");
2120
2121 const BoxArray& ba = fa.boxArray();
2122 const DistributionMapping& dm = fa.DistributionMap();
2123
2124 FabArray<BaseFab<int>> mask(ba, dm, 1, nghost);
2125 mask.setVal(1);
2126
2127 const std::vector<IntVect>& pshifts = period.shiftIntVect();
2128
2130
2131 bool run_on_gpu = Gpu::inLaunchRegion();
2132 amrex::ignore_unused(run_on_gpu, tags);
2133#ifdef AMREX_USE_OMP
2134#pragma omp parallel if (!run_on_gpu)
2135#endif
2136 {
2137 std::vector< std::pair<int,Box> > isects;
2138
2139 for (MFIter mfi(mask); mfi.isValid(); ++mfi)
2140 {
2141 const Box& bx = mask[mfi].box();
2142 auto const& arr = mask.array(mfi);
2143
2144 for (const auto& iv : pshifts)
2145 {
2146 ba.intersections(bx+iv, isects, false, nghost);
2147 for (const auto& is : isects)
2148 {
2149 Box const& b = is.second-iv;
2150 if (iv == 0 && b == bx) { continue; }
2151#ifdef AMREX_USE_GPU
2152 if (run_on_gpu) {
2153 tags.push_back(Array4BoxTag<int>{.dfab = arr, .dbox = b});
2154 } else
2155#endif
2156 {
2157 amrex::LoopConcurrentOnCpu(b, [=] (int i, int j, int k) noexcept
2158 {
2159 arr(i,j,k) += 1;
2160 });
2161 }
2162 }
2163 }
2164 }
2165 }
2166
2167#ifdef AMREX_USE_GPU
2168 amrex::ParallelFor(tags, 1,
2169 [=] AMREX_GPU_DEVICE (int i, int j, int k, int n, Array4BoxTag<int> const& tag) noexcept
2170 {
2171 Gpu::Atomic::AddNoRet(tag.dfab.ptr(i,j,k,n), 1);
2172 });
2173#endif
2174
2175 return mask;
2176}
2177
2178}
2179
2180#endif
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define BL_ASSERT(EX)
Definition AMReX_BLassert.H:39
#define AMREX_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:37
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_HOST_DEVICE_PARALLEL_FOR_4D(...)
Definition AMReX_GpuLaunchMacrosC.nolint.H:111
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
Real * pdst
Definition AMReX_HypreMLABecLap.cpp:1140
Array4< int const > mask
Definition AMReX_InterpFaceRegister.cpp:93
#define AMREX_LOOP_3D(bx, i, j, k, block)
Definition AMReX_Loop.nolint.H:4
#define AMREX_LOOP_4D(bx, ncomp, i, j, k, n, block)
Definition AMReX_Loop.nolint.H:16
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Print on all processors of the default communicator.
Definition AMReX_Print.H:113
virtual bool isManaged() const
Check whether it is managed GPU memory.
Definition AMReX_Arena.cpp:88
virtual bool isDevice() const
Check whether it is non-managed GPU device memory.
Definition AMReX_Arena.cpp:100
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:564
std::vector< std::pair< int, Box > > intersections(const Box &bx) const
Return intersections of Box and BoxArray.
Definition AMReX_BoxArray.cpp:1186
__host__ __device__ bool contains(const IntVectND< dim > &p) const noexcept
Return true if argument is contained within BoxND.
Definition AMReX_Box.H:216
__host__ __device__ bool ok() const noexcept
Checks if it is a proper BoxND (including a valid type).
Definition AMReX_Box.H:212
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
IntVect nGrowVect() const noexcept
Definition AMReX_FabArrayBase.H:80
bool isFusingCandidate() const noexcept
Is this a good candidate for kernel fusing?
Definition AMReX_FabArrayBase.cpp:2705
IndexType ixType() const noexcept
Return index type.
Definition AMReX_FabArrayBase.H:86
const DistributionMapping & DistributionMap() const noexcept
Return constant reference to associated DistributionMapping.
Definition AMReX_FabArrayBase.H:131
@ ADD
Definition AMReX_FabArrayBase.H:394
int nComp() const noexcept
Return number of variables (aka components) associated with each point.
Definition AMReX_FabArrayBase.H:83
const BoxArray & boxArray() const noexcept
Return a constant reference to the BoxArray that defines the valid region associated with this FabArr...
Definition AMReX_FabArrayBase.H:95
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:344
void mult(value_type val, int comp, int num_comp, int nghost=0)
Definition AMReX_FabArray.H:2968
void setBndry(value_type val)
Set all values in the boundary region to val.
Definition AMReX_FabArray.H:2487
std::unique_ptr< FabArray< FAB > > os_temp
Definition AMReX_FabArray.H:1661
const FabFactory< FAB > & Factory() const noexcept
Definition AMReX_FabArray.H:442
void prefetchToHost(const MFIter &mfi) const noexcept
Definition AMReX_FabArray.H:541
MultiArray4< typename FabArray< FAB >::value_type > arrays() noexcept
Definition AMReX_FabArray.H:633
void prefetchToDevice(const MFIter &mfi) const noexcept
Definition AMReX_FabArray.H:551
Arena * arena() const noexcept
Definition AMReX_FabArray.H:445
MultiArray4< typename FabArray< FAB >::value_type const > const_arrays() const noexcept
Definition AMReX_FabArray.H:647
void setVal(value_type val)
Set all components in the entire region of each FAB to val.
Definition AMReX_FabArray.H:2705
Array4< typename FabArray< FAB >::value_type const > array(const MFIter &mfi) const noexcept
Definition AMReX_FabArray.H:561
bool hasEBFabFactory() const noexcept
Definition AMReX_FabArray.H:449
Array4< typename FabArray< FAB >::value_type const > const_array(const MFIter &mfi) const noexcept
Definition AMReX_FabArray.H:585
GPU-compatible tuple.
Definition AMReX_Tuple.H:98
Definition AMReX_GpuBuffer.H:24
T const * data() const noexcept
Definition AMReX_GpuBuffer.H:51
__host__ __device__ bool cellCentered() const noexcept
True if the IndexTypeND is CELL based in all directions.
Definition AMReX_IndexType.H:102
__host__ __device__ constexpr bool allGT(const IntVectND< dim > &rhs) const noexcept
Returns true if this is greater than argument for all components. NOTE: This is NOT a strict weak ord...
Definition AMReX_IntVect.H:517
__host__ static __device__ constexpr IntVectND< dim > TheZeroVector() noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:771
__host__ static __device__ constexpr IntVectND< dim > TheMinVector() noexcept
Definition AMReX_IntVect.H:819
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
T * data() noexcept
Definition AMReX_PODVector.H:666
This provides length of period for periodic domains. 0 means it is not periodic in that direction....
Definition AMReX_Periodicity.H:17
static const Periodicity & NonPeriodic() noexcept
Definition AMReX_Periodicity.cpp:52
std::vector< IntVect > shiftIntVect(IntVect const &nghost=IntVect(0)) const
Definition AMReX_Periodicity.cpp:8
Print & SetPrecision(int p)
Definition AMReX_Print.H:86
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:29
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1289
std::array< T, N > Array
Definition AMReX_Array.H:26
void Sum(T &v, MPI_Comm comm)
Definition AMReX_ParallelReduce.H:221
void Max(KeyValuePair< K, V > &vi, MPI_Comm comm)
Definition AMReX_ParallelReduce.H:133
__host__ __device__ AMREX_FORCE_INLINE T Exch(T *address, T val) noexcept
Definition AMReX_GpuAtomic.H:487
__host__ __device__ AMREX_FORCE_INLINE void AddNoRet(T *sum, T value) noexcept
Definition AMReX_GpuAtomic.H:283
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
void dtoh_memcpy_async(void *p_h, const void *p_d, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:435
bool inLaunchRegion() noexcept
Definition AMReX_GpuControl.H:88
bool inNoSyncRegion() noexcept
Definition AMReX_GpuControl.H:148
void htod_memcpy_async(void *p_d, const void *p_h, const std::size_t sz) noexcept
Definition AMReX_GpuDevice.H:421
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
Definition AMReX_Amr.cpp:50
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
__host__ __device__ BoxND< dim > makeSingleCellBox(int i, int j, int k, IndexTypeND< dim > typ=IndexTypeND< dim >::TheCellType())
Definition AMReX_Box.H:2152
int nComp(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2852
void Copy(FabArray< DFAB > &dst, FabArray< SFAB > const &src, int srccomp, int dstcomp, int numcomp, int nghost)
Definition AMReX_FabArray.H:180
void htod_memcpy(FabArray< FAB > &dst, FabArray< FAB > const &src, int scomp, int dcomp, int ncomp)
Definition AMReX_FabArrayUtility.H:1476
__host__ __device__ void Swap(T &t1, T &t2) noexcept
Definition AMReX_Algorithm.H:94
void OverrideSync_nowait(FabArray< FAB > &fa, FabArray< IFAB > const &msk, const Periodicity &period)
Definition AMReX_FabArrayUtility.H:1390
void Add(FabArray< FAB > &dst, FabArray< FAB > const &src, int srccomp, int dstcomp, int numcomp, int nghost)
Definition AMReX_FabArray.H:239
ReduceData< Ts... >::Type ParReduce(TypeList< Ops... > operation_list, TypeList< Ts... > type_list, FabArray< FAB > const &fa, IntVect const &nghost, F &&f)
Parallel reduce for MultiFab/FabArray. The reduce result is local and it's the user's responsibility ...
Definition AMReX_ParReduce.H:48
DistributionMapping const & DistributionMap(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2867
IntVect nGrowVect(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2857
void Scale(MF &dst, typename MF::value_type val, int scomp, int ncomp, int nghost)
dst *= val
Definition AMReX_FabArrayUtility.H:1875
FAB::value_type Dot(FabArray< FAB > const &x, int xcomp, FabArray< FAB > const &y, int ycomp, int ncomp, IntVect const &nghost, bool local=false)
Compute dot products of two FabArrays.
Definition AMReX_FabArrayUtility.H:1609
void Saxpy_Xpay(MF &dst, typename MF::value_type a_saxpy, MF const &src_saxpy, typename MF::value_type a_xpay, MF const &src_xpay, int scomp, int dcomp, int ncomp, IntVect const &nghost)
dst += a_saxpy * src_saxpy followed by dst = src_xpay + a_xpay * dst
Definition AMReX_FabArrayUtility.H:1914
void prefetchToDevice(FabArray< FAB > const &fa, const bool synchronous=true)
Definition AMReX_FabArrayUtility.H:1363
FAB::value_type ReduceMin(FabArray< FAB > const &fa, int nghost, F &&f)
Definition AMReX_FabArrayUtility.H:304
void ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:202
FabArray< BaseFab< int > > OverlapMask(FabArray< FAB > const &fa, IntVect const &nghost, Periodicity const &period)
Definition AMReX_FabArrayUtility.H:2117
void prefetchToHost(FabArray< FAB > const &fa, const bool synchronous=true)
Definition AMReX_FabArrayUtility.H:1348
bool isMFIterSafe(const FabArrayBase &x, const FabArrayBase &y)
Definition AMReX_MFIter.H:252
void Saxpy_Saxpy(MF &dst1, typename MF::value_type a1, MF const &src1, MF &dst2, typename MF::value_type a2, MF const &src2, int scomp, int dcomp, int ncomp, IntVect const &nghost)
dst1 += a1 * src1 followed by dst2 += a2 * src2
Definition AMReX_FabArrayUtility.H:1923
IntVect indexFromValue(FabArray< FAB > const &mf, int comp, IntVect const &nghost, typename FAB::value_type value)
Definition AMReX_FabArrayUtility.H:1501
void LinComb(MF &dst, typename MF::value_type a, MF const &src_a, int acomp, typename MF::value_type b, MF const &src_b, int bcomp, int dcomp, int ncomp, IntVect const &nghost)
dst = a*src_a + b*src_b
Definition AMReX_FabArrayUtility.H:1941
void ParallelCopy(MF &dst, MF const &src, int scomp, int dcomp, int ncomp, IntVect const &ng_src=IntVect(0), IntVect const &ng_dst=IntVect(0), Periodicity const &period=Periodicity::NonPeriodic())
dst = src w/ MPI communication
Definition AMReX_FabArrayUtility.H:1951
void Saypy_Saxpy(MF &dst1, typename MF::value_type a1, MF &dst2, typename MF::value_type a2, MF const &src, int scomp, int dcomp, int ncomp, IntVect const &nghost)
dst1 += a1 * dst2 followed by dst2 += a2 * src
Definition AMReX_FabArrayUtility.H:1932
void Multiply(FabArray< FAB > &dst, FabArray< FAB > const &src, int srccomp, int dstcomp, int numcomp, int nghost)
Definition AMReX_FabArrayUtility.H:1214
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void Abs(FabArray< FAB > &fa, int icomp, int numcomp, int nghost)
Definition AMReX_FabArrayUtility.H:1305
void OverrideSync_finish(FabArray< FAB > &fa)
Definition AMReX_FabArrayUtility.H:1437
void single_task(L &&f) noexcept
Definition AMReX_GpuLaunchFunctsC.H:1239
void dtoh_memcpy(FabArray< FAB > &dst, FabArray< FAB > const &src, int scomp, int dcomp, int ncomp)
Definition AMReX_FabArrayUtility.H:1451
void LocalCopy(DMF &dst, SMF const &src, int scomp, int dcomp, int ncomp, IntVect const &nghost)
dst = src
Definition AMReX_FabArrayUtility.H:1882
bool TilingIfNotGPU() noexcept
Definition AMReX_MFIter.H:12
void setBndry(MF &dst, typename MF::value_type val, int scomp, int ncomp)
dst = val in ghost cells.
Definition AMReX_FabArrayUtility.H:1868
void printCell(FabArray< FAB > const &mf, const IntVect &cell, int comp=-1, const IntVect &ng=IntVect::TheZeroVector())
Definition AMReX_FabArrayUtility.H:1051
MF::value_type norminf(MF const &mf, int scomp, int ncomp, IntVect const &nghost, bool local=false)
Definition AMReX_FabArrayUtility.H:1961
bool ReduceLogicalOr(FabArray< FAB > const &fa, int nghost, F &&f)
Definition AMReX_FabArrayUtility.H:904
void LoopConcurrentOnCpu(Dim3 lo, Dim3 hi, F const &f) noexcept
Definition AMReX_Loop.H:388
void Xpay(MF &dst, typename MF::value_type a, MF const &src, int scomp, int dcomp, int ncomp, IntVect const &nghost)
dst = src + a * dst
Definition AMReX_FabArrayUtility.H:1906
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241
const int[]
Definition AMReX_BLProfiler.cpp:1664
void OverrideSync(FabArray< FAB > &fa, FabArray< IFAB > const &msk, const Periodicity &period)
Definition AMReX_FabArrayUtility.H:1379
void Divide(FabArray< FAB > &dst, FabArray< FAB > const &src, int srccomp, int dstcomp, int numcomp, int nghost)
Definition AMReX_FabArrayUtility.H:1260
bool ReduceLogicalAnd(FabArray< FAB > const &fa, int nghost, F &&f)
Definition AMReX_FabArrayUtility.H:757
FAB::value_type ReduceMax(FabArray< FAB > const &fa, int nghost, F &&f)
Definition AMReX_FabArrayUtility.H:530
void Subtract(FabArray< FAB > &dst, FabArray< FAB > const &src, int srccomp, int dstcomp, int numcomp, int nghost)
Definition AMReX_FabArrayUtility.H:1168
void Saxpy(MF &dst, typename MF::value_type a, MF const &src, int scomp, int dcomp, int ncomp, IntVect const &nghost)
dst += a * src
Definition AMReX_FabArrayUtility.H:1898
void LocalAdd(MF &dst, MF const &src, int scomp, int dcomp, int ncomp, IntVect const &nghost)
dst += src
Definition AMReX_FabArrayUtility.H:1890
BoxArray const & boxArray(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2862
FAB::value_type ReduceSum(FabArray< FAB > const &fa, int nghost, F &&f)
Definition AMReX_FabArrayUtility.H:16
void setVal(MF &dst, typename MF::value_type val)
dst = val
Definition AMReX_FabArrayUtility.H:1861
Definition AMReX_TagParallelFor.H:58
Array4< T > dfab
Definition AMReX_TagParallelFor.H:59
Definition AMReX_TypeTraits.H:94
FabArray memory allocation information.
Definition AMReX_FabArray.H:68
Definition AMReX_MFIter.H:20
Struct for holding types.
Definition AMReX_TypeList.H:13