Line data Source code
1 : #ifndef AMREX_FillPatchUtil_H_
2 : #define AMREX_FillPatchUtil_H_
3 : #include <AMReX_Config.H>
4 :
5 : #include <AMReX_MultiFab.H>
6 : #include <AMReX_iMultiFab.H>
7 : #include <AMReX_Geometry.H>
8 : #include <AMReX_PhysBCFunct.H>
9 : #include <AMReX_Interpolater.H>
10 : #include <AMReX_MFInterpolater.H>
11 : #include <AMReX_Array.H>
12 : #include <AMReX_Utility.H>
13 :
14 : #ifdef AMREX_USE_EB
15 : #include <AMReX_EB2.H>
16 : #include <AMReX_EBFabFactory.H>
17 : #include <AMReX_EBInterpolater.H>
18 : #include <AMReX_EBMFInterpolater.H>
19 : #endif
20 :
21 : #ifdef AMREX_USE_OMP
22 : #include <omp.h>
23 : #endif
24 :
25 : #include <cmath>
26 : #include <limits>
27 :
28 : namespace amrex
29 : {
30 :
31 : template <typename MFFAB>
32 : struct NullInterpHook
33 : {
34 : template <class F=MFFAB, std::enable_if_t<IsBaseFab<F>::value,int> = 0>
35 432 : void operator() (MFFAB& /*fab*/, const Box& /*bx*/, int /*icomp*/, int /*ncomp*/) const {}
36 :
37 : template <class F=MFFAB, std::enable_if_t<IsBaseFab<F>::value,int> = 0>
38 : void operator() (Array<MFFAB*, AMREX_SPACEDIM> /*fab*/, const Box& /*bx*/, int /*icomp*/, int /*ncomp*/) const {}
39 :
40 : template <class F=MFFAB, std::enable_if_t<IsFabArray<F>::value,int> = 0>
41 : void operator() (MFFAB& /*mf*/, int /*icomp*/, int /*ncomp*/) const {}
42 : };
43 :
44 : /**
45 : * \brief Test if AMR grids are properly nested.
46 : *
47 : * If grids are not properly nested, FillPatch functions may fail.
48 : *
49 : * \tparam Interp Interpolater type
50 : *
51 : * \param ratio refinement ratio
52 : * \param blocking_factor blocking factor on the fine level
53 : * \param ngrow number of ghost cells of fine MultiFab
54 : * \param boxType index type
55 : * \param mapper an interpolater object
56 : */
57 : template <typename Interp>
58 : bool ProperlyNested (const IntVect& ratio, const IntVect& blocking_factor, int ngrow,
59 : const IndexType& boxType, Interp* mapper);
60 :
61 : /**
62 : * \brief FillPatch with data from the current level
63 : *
64 : * The destination MultiFab/FabArray is on the same AMR level as the
65 : * source MultiFab/FabArray. Usually this can only be used on AMR level
66 : * 0, because filling fine level MF usually requires coarse level
67 : * data. If needed, interpolation in time is performed.
68 : *
69 : * \tparam MF the MultiFab/FabArray type
70 : * \tparam BC functor for filling physical boundaries
71 : *
72 : * \param mf destination MF
73 : * \param nghost number of ghost cells of mf needed to be filled
74 : * \param time time associated with mf
75 : * \param smf source MFs
76 : * \param stime times associated smf
77 : * \param scomp starting component of the source MFs
78 : * \param dcomp starting component of the destination MF
79 : * \param ncomp number of components
80 : * \param geom Geometry for this level
81 : * \param physbcf functor for physical boundaries
82 : * \param bcfcomp starting component for physbcf
83 : */
84 : template <typename MF, typename BC>
85 : std::enable_if_t<IsFabArray<MF>::value>
86 : FillPatchSingleLevel (MF& mf, IntVect const& nghost, Real time,
87 : const Vector<MF*>& smf, const Vector<Real>& stime,
88 : int scomp, int dcomp, int ncomp,
89 : const Geometry& geom,
90 : BC& physbcf, int bcfcomp);
91 :
92 : /**
93 : * \brief FillPatch with data from the current level
94 : *
95 : * The destination MultiFab/FabArray is on the same AMR level as the
96 : * source MultiFab/FabArray. Usually this can only be used on AMR level
97 : * 0, because filling fine level MF usually requires coarse level
98 : * data. If needed, interpolation in time is performed. Ghost cells of
99 : * the destination MF are not filled.
100 : *
101 : * \tparam MF the MultiFab/FabArray type
102 : * \tparam BC functor for filling physical boundaries
103 : *
104 : * \param mf destination MF
105 : * \param time time associated with mf
106 : * \param smf source MFs
107 : * \param stime times associated smf
108 : * \param scomp starting component of the source MFs
109 : * \param dcomp starting component of the destination MF
110 : * \param ncomp number of components
111 : * \param geom Geometry for this level
112 : * \param physbcf functor for physical boundaries
113 : * \param bcfcomp starting component for physbcf
114 : */
115 : template <typename MF, typename BC>
116 : std::enable_if_t<IsFabArray<MF>::value>
117 : FillPatchSingleLevel (MF& mf, Real time,
118 : const Vector<MF*>& smf, const Vector<Real>& stime,
119 : int scomp, int dcomp, int ncomp,
120 : const Geometry& geom,
121 : BC& physbcf, int bcfcomp);
122 :
123 : /**
124 : * \brief FillPatch with data from the current level and the level below.
125 : *
126 : * First, we fill the destination MultiFab/FabArray with the current
127 : * level data as much as possible. This may include interpolation in
128 : * time. For the rest of the destination MF, we fill them with the
129 : * coarse level data using interpolation in space (and in time if
130 : * needed).
131 : *
132 : * \tparam MF the MultiFab/FabArray type
133 : * \tparam BC functor for filling physical boundaries
134 : * \tparam Interp spatial interpolater
135 : * \tparam PreInterpHook pre-interpolation hook
136 : * \tparam PostInterpHook post-interpolation hook
137 : *
138 : * \param mf destination MF on the fine level
139 : * \param nghost number of ghost cells of mf needed to be filled
140 : * \param time time associated with mf
141 : * \param cmf source MFs on the coarse level
142 : * \param ct times associated cmf
143 : * \param fmf source MFs on the fine level
144 : * \param ft times associated fmf
145 : * \param scomp starting component of the source MFs
146 : * \param dcomp starting component of the destination MF
147 : * \param ncomp number of components
148 : * \param cgeom Geometry for the coarse level
149 : * \param fgeom Geometry for the fine level
150 : * \param cbc functor for physical boundaries on the coarse level
151 : * \param cbccomp starting component for cbc
152 : * \param fbc functor for physical boundaries on the fine level
153 : * \param fbccomp starting component for fbc
154 : * \param ratio refinement ratio
155 : * \param mapper spatial interpolater
156 : * \param bcs boundary types for each component. We need this because some interpolaters need it.
157 : * \param bcscomp starting component for bcs
158 : * \param pre_interp pre-interpolation hook
159 : * \param post_interp post-interpolation hook
160 : */
161 : template <typename MF, typename BC, typename Interp,
162 : typename PreInterpHook=NullInterpHook<typename MF::FABType::value_type>,
163 : typename PostInterpHook=NullInterpHook<typename MF::FABType::value_type> >
164 : std::enable_if_t<IsFabArray<MF>::value>
165 : FillPatchTwoLevels (MF& mf, IntVect const& nghost, Real time,
166 : const Vector<MF*>& cmf, const Vector<Real>& ct,
167 : const Vector<MF*>& fmf, const Vector<Real>& ft,
168 : int scomp, int dcomp, int ncomp,
169 : const Geometry& cgeom, const Geometry& fgeom,
170 : BC& cbc, int cbccomp,
171 : BC& fbc, int fbccomp,
172 : const IntVect& ratio,
173 : Interp* mapper,
174 : const Vector<BCRec>& bcs, int bcscomp,
175 : const PreInterpHook& pre_interp = {},
176 : const PostInterpHook& post_interp = {});
177 :
178 : /**
179 : * \brief FillPatch with data from the current level and the level below.
180 : *
181 : * First, we fill the destination MultiFab/FabArray with the current
182 : * level data as much as possible. This may include interpolation in
183 : * time. For the rest of the destination MF, we fill them with the
184 : * coarse level data using interpolation in space (and in time if
185 : * needed). Ghost cells of the destination MF are not filled.
186 : *
187 : * \tparam MF the MultiFab/FabArray type
188 : * \tparam BC functor for filling physical boundaries
189 : * \tparam Interp spatial interpolater
190 : * \tparam PreInterpHook pre-interpolation hook
191 : * \tparam PostInterpHook post-interpolation hook
192 : *
193 : * \param mf destination MF on the fine level
194 : * \param time time associated with mf
195 : * \param cmf source MFs on the coarse level
196 : * \param ct times associated cmf
197 : * \param fmf source MFs on the fine level
198 : * \param ft times associated fmf
199 : * \param scomp starting component of the source MFs
200 : * \param dcomp starting component of the destination MF
201 : * \param ncomp number of components
202 : * \param cgeom Geometry for the coarse level
203 : * \param fgeom Geometry for the fine level
204 : * \param cbc functor for physical boundaries on the coarse level
205 : * \param cbccomp starting component for cbc
206 : * \param fbc functor for physical boundaries on the fine level
207 : * \param fbccomp starting component for fbc
208 : * \param ratio refinement ratio
209 : * \param mapper spatial interpolater
210 : * \param bcs boundary types for each component. We need this because some interpolaters need it.
211 : * \param bcscomp starting component for bcs
212 : * \param pre_interp pre-interpolation hook
213 : * \param post_interp post-interpolation hook
214 : */
215 : template <typename MF, typename BC, typename Interp,
216 : typename PreInterpHook=NullInterpHook<typename MF::FABType::value_type>,
217 : typename PostInterpHook=NullInterpHook<typename MF::FABType::value_type> >
218 : std::enable_if_t<IsFabArray<MF>::value>
219 : FillPatchTwoLevels (MF& mf, Real time,
220 : const Vector<MF*>& cmf, const Vector<Real>& ct,
221 : const Vector<MF*>& fmf, const Vector<Real>& ft,
222 : int scomp, int dcomp, int ncomp,
223 : const Geometry& cgeom, const Geometry& fgeom,
224 : BC& cbc, int cbccomp,
225 : BC& fbc, int fbccomp,
226 : const IntVect& ratio,
227 : Interp* mapper,
228 : const Vector<BCRec>& bcs, int bcscomp,
229 : const PreInterpHook& pre_interp = {},
230 : const PostInterpHook& post_interp = {});
231 :
232 : /**
233 : * \brief FillPatch for face variables with data from the current level
234 : * and the level below. Sometimes, we need to fillpatch all
235 : * AMREX_SPACEDIM face MultiFabs togother to satisfy certain constraint
236 : * such as divergence preserving.
237 : *
238 : * First, we fill the destination MultiFab/FabArray's with the current
239 : * level data as much as possible. This may include interpolation in
240 : * time. For the rest of the destination MFs, we fill them with the
241 : * coarse level data using interpolation in space (and in time if
242 : * needed).
243 : *
244 : * \tparam MF the MultiFab/FabArray type
245 : * \tparam BC functor for filling physical boundaries
246 : * \tparam Interp spatial interpolater
247 : * \tparam PreInterpHook pre-interpolation hook
248 : * \tparam PostInterpHook post-interpolation hook
249 : *
250 : * \param mf destination MFs on the fine level
251 : * \param nghost number of ghost cells of mf needed to be filled
252 : * \param time time associated with mf
253 : * \param cmf source MFs on the coarse level
254 : * \param ct times associated cmf
255 : * \param fmf source MFs on the fine level
256 : * \param ft times associated fmf
257 : * \param scomp starting component of the source MFs
258 : * \param dcomp starting component of the destination MFs
259 : * \param ncomp number of components
260 : * \param cgeom Geometry for the coarse level
261 : * \param fgeom Geometry for the fine level
262 : * \param cbc functor for physical boundaries on the coarse level
263 : * \param cbccomp starting component for cbc
264 : * \param fbc functor for physical boundaries on the fine level
265 : * \param fbccomp starting component for fbc
266 : * \param ratio refinement ratio
267 : * \param mapper spatial interpolater
268 : * \param bcs boundary types for each component. We need this because some interpolaters need it.
269 : * \param bcscomp starting component for bcs
270 : * \param pre_interp pre-interpolation hook
271 : * \param post_interp post-interpolation hook
272 : */
273 : template <typename MF, typename BC, typename Interp,
274 : typename PreInterpHook=NullInterpHook<typename MF::FABType::value_type>,
275 : typename PostInterpHook=NullInterpHook<typename MF::FABType::value_type> >
276 : std::enable_if_t<IsFabArray<MF>::value>
277 : FillPatchTwoLevels (Array<MF*, AMREX_SPACEDIM> const& mf, IntVect const& nghost, Real time,
278 : const Vector<Array<MF*, AMREX_SPACEDIM> >& cmf, const Vector<Real>& ct,
279 : const Vector<Array<MF*, AMREX_SPACEDIM> >& fmf, const Vector<Real>& ft,
280 : int scomp, int dcomp, int ncomp,
281 : const Geometry& cgeom, const Geometry& fgeom,
282 : Array<BC, AMREX_SPACEDIM>& cbc, const Array<int, AMREX_SPACEDIM>& cbccomp,
283 : Array<BC, AMREX_SPACEDIM>& fbc, const Array<int, AMREX_SPACEDIM>& fbccomp,
284 : const IntVect& ratio,
285 : Interp* mapper,
286 : const Array<Vector<BCRec>, AMREX_SPACEDIM>& bcs, const Array<int, AMREX_SPACEDIM>& bcscomp,
287 : const PreInterpHook& pre_interp = {},
288 : const PostInterpHook& post_interp = {});
289 :
290 : /**
291 : * \brief FillPatch for face variables with data from the current level
292 : * and the level below. Sometimes, we need to fillpatch all
293 : * AMREX_SPACEDIM face MultiFabs togother to satisfy certain constraint
294 : * such as divergence preserving.
295 : *
296 : * First, we fill the destination MultiFab/FabArray's with the current
297 : * level data as much as possible. This may include interpolation in
298 : * time. For the rest of the destination MFs, we fill them with the
299 : * coarse level data using interpolation in space (and in time if
300 : * needed).
301 : *
302 : * \tparam MF the MultiFab/FabArray type
303 : * \tparam BC functor for filling physical boundaries
304 : * \tparam Interp spatial interpolater
305 : * \tparam PreInterpHook pre-interpolation hook
306 : * \tparam PostInterpHook post-interpolation hook
307 : *
308 : * \param mf destination MFs on the fine level
309 : * \param nghost number of ghost cells of mf needed to be filled
310 : * \param time time associated with mf
311 : * \param cmf source MFs on the coarse level
312 : * \param ct times associated cmf
313 : * \param fmf source MFs on the fine level
314 : * \param ft times associated fmf
315 : * \param scomp starting component of the source MFs
316 : * \param dcomp starting component of the destination MFs
317 : * \param ncomp number of components
318 : * \param cgeom Geometry for the coarse level
319 : * \param fgeom Geometry for the fine level
320 : * \param cbc functor for physical boundaries on the coarse level
321 : * \param cbccomp starting component for cbc
322 : * \param fbc functor for physical boundaries on the fine level
323 : * \param fbccomp starting component for fbc
324 : * \param ratio refinement ratio
325 : * \param mapper spatial interpolater
326 : * \param bcs boundary types for each component. We need this because some interpolaters need it.
327 : * \param bcscomp starting component for bcs
328 : * \param pre_interp pre-interpolation hook
329 : * \param post_interp post-interpolation hook
330 : */
331 : template <typename MF, typename BC, typename Interp,
332 : typename PreInterpHook=NullInterpHook<typename MF::FABType::value_type>,
333 : typename PostInterpHook=NullInterpHook<typename MF::FABType::value_type> >
334 : std::enable_if_t<IsFabArray<MF>::value>
335 : FillPatchTwoLevels (Array<MF*, AMREX_SPACEDIM> const& mf, IntVect const& nghost, Real time,
336 : const Vector<Array<MF*, AMREX_SPACEDIM> >& cmf, const Vector<Real>& ct,
337 : const Vector<Array<MF*, AMREX_SPACEDIM> >& fmf, const Vector<Real>& ft,
338 : int scomp, int dcomp, int ncomp,
339 : const Geometry& cgeom, const Geometry& fgeom,
340 : Array<BC, AMREX_SPACEDIM>& cbc, int cbccomp,
341 : Array<BC, AMREX_SPACEDIM>& fbc, int fbccomp,
342 : const IntVect& ratio,
343 : Interp* mapper,
344 : const Array<Vector<BCRec>, AMREX_SPACEDIM>& bcs, int bcscomp,
345 : const PreInterpHook& pre_interp = {},
346 : const PostInterpHook& post_interp = {});
347 :
348 : /**
349 : * \brief FillPatch for face variables with data from the current level
350 : * and the level below. Sometimes, we need to fillpatch all
351 : * AMREX_SPACEDIM face MultiFabs togother to satisfy certain constraint
352 : * such as divergence preserving.
353 : *
354 : * First, we fill the destination MultiFab/FabArray's with the current
355 : * level data as much as possible. This may include interpolation in
356 : * time. For the rest of the destination MFs, we fill them with the
357 : * coarse level data using interpolation in space (and in time if
358 : * needed). Ghost cells of the destination MFs are not filled.
359 : *
360 : * \tparam MF the MultiFab/FabArray type
361 : * \tparam BC functor for filling physical boundaries
362 : * \tparam Interp spatial interpolater
363 : * \tparam PreInterpHook pre-interpolation hook
364 : * \tparam PostInterpHook post-interpolation hook
365 : *
366 : * \param mf destination MFs on the fine level
367 : * \param time time associated with mf
368 : * \param cmf source MFs on the coarse level
369 : * \param ct times associated cmf
370 : * \param fmf source MFs on the fine level
371 : * \param ft times associated fmf
372 : * \param scomp starting component of the source MFs
373 : * \param dcomp starting component of the destination MFs
374 : * \param ncomp number of components
375 : * \param cgeom Geometry for the coarse level
376 : * \param fgeom Geometry for the fine level
377 : * \param cbc functor for physical boundaries on the coarse level
378 : * \param cbccomp starting component for cbc
379 : * \param fbc functor for physical boundaries on the fine level
380 : * \param fbccomp starting component for fbc
381 : * \param ratio refinement ratio
382 : * \param mapper spatial interpolater
383 : * \param bcs boundary types for each component. We need this because some interpolaters need it.
384 : * \param bcscomp starting component for bcs
385 : * \param pre_interp pre-interpolation hook
386 : * \param post_interp post-interpolation hook
387 : */
388 : template <typename MF, typename BC, typename Interp,
389 : typename PreInterpHook=NullInterpHook<typename MF::FABType::value_type>,
390 : typename PostInterpHook=NullInterpHook<typename MF::FABType::value_type> >
391 : std::enable_if_t<IsFabArray<MF>::value>
392 : FillPatchTwoLevels (Array<MF*, AMREX_SPACEDIM> const& mf, Real time,
393 : const Vector<Array<MF*, AMREX_SPACEDIM> >& cmf, const Vector<Real>& ct,
394 : const Vector<Array<MF*, AMREX_SPACEDIM> >& fmf, const Vector<Real>& ft,
395 : int scomp, int dcomp, int ncomp,
396 : const Geometry& cgeom, const Geometry& fgeom,
397 : Array<BC, AMREX_SPACEDIM>& cbc, int cbccomp,
398 : Array<BC, AMREX_SPACEDIM>& fbc, int fbccomp,
399 : const IntVect& ratio,
400 : Interp* mapper,
401 : const Array<Vector<BCRec>, AMREX_SPACEDIM>& bcs, int bcscomp,
402 : const PreInterpHook& pre_interp = {},
403 : const PostInterpHook& post_interp = {});
404 :
405 : #ifdef AMREX_USE_EB
406 : /**
407 : * \brief FillPatch with data from the current level and the level below.
408 : *
409 : * First, we fill the destination MultiFab/FabArray with the current
410 : * level data as much as possible. This may include interpolation in
411 : * time. For the rest of the destination MF, we fill them with the
412 : * coarse level data using interpolation in space (and in time if
413 : * needed).
414 : *
415 : * \tparam MF the MultiFab/FabArray type
416 : * \tparam BC functor for filling physical boundaries
417 : * \tparam Interp spatial interpolater
418 : * \tparam PreInterpHook pre-interpolation hook
419 : * \tparam PostInterpHook post-interpolation hook
420 : *
421 : * \param mf destination MF on the fine level
422 : * \param nghost number of ghost cells of mf needed to be filled
423 : * \param time time associated with mf
424 : * \param index_space EB IndexSpace
425 : * \param cmf source MFs on the coarse level
426 : * \param ct times associated cmf
427 : * \param fmf source MFs on the fine level
428 : * \param ft times associated fmf
429 : * \param scomp starting component of the source MFs
430 : * \param dcomp starting component of the destination MF
431 : * \param ncomp number of components
432 : * \param cgeom Geometry for the coarse level
433 : * \param fgeom Geometry for the fine level
434 : * \param cbc functor for physical boundaries on the coarse level
435 : * \param cbccomp starting component for cbc
436 : * \param fbc functor for physical boundaries on the fine level
437 : * \param fbccomp starting component for fbc
438 : * \param ratio refinement ratio
439 : * \param mapper spatial interpolater
440 : * \param bcs boundary types for each component. We need this because some interpolaters need it.
441 : * \param bcscomp starting component for bcs
442 : * \param pre_interp pre-interpolation hook
443 : * \param post_interp post-interpolation hook
444 : */
445 : template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
446 : std::enable_if_t<IsFabArray<MF>::value>
447 : FillPatchTwoLevels (MF& mf, IntVect const& nghost, Real time,
448 : const EB2::IndexSpace& index_space,
449 : const Vector<MF*>& cmf, const Vector<Real>& ct,
450 : const Vector<MF*>& fmf, const Vector<Real>& ft,
451 : int scomp, int dcomp, int ncomp,
452 : const Geometry& cgeom, const Geometry& fgeom,
453 : BC& cbc, int cbccomp,
454 : BC& fbc, int fbccomp,
455 : const IntVect& ratio,
456 : Interp* mapper,
457 : const Vector<BCRec>& bcs, int bcscomp,
458 : const PreInterpHook& pre_interp,
459 : const PostInterpHook& post_interp);
460 :
461 : /**
462 : * \brief FillPatch with data from the current level and the level below.
463 : *
464 : * First, we fill the destination MultiFab/FabArray with the current
465 : * level data as much as possible. This may include interpolation in
466 : * time. For the rest of the destination MF, we fill them with the
467 : * coarse level data using interpolation in space (and in time if
468 : * needed). Ghost cells of the destination MF are not filled.
469 : *
470 : * \tparam MF the MultiFab/FabArray type
471 : * \tparam BC functor for filling physical boundaries
472 : * \tparam Interp spatial interpolater
473 : * \tparam PreInterpHook pre-interpolation hook
474 : * \tparam PostInterpHook post-interpolation hook
475 : *
476 : * \param mf destination MF on the fine level
477 : * \param time time associated with mf
478 : * \param index_space EB IndexSpace
479 : * \param cmf source MFs on the coarse level
480 : * \param ct times associated cmf
481 : * \param fmf source MFs on the fine level
482 : * \param ft times associated fmf
483 : * \param scomp starting component of the source MFs
484 : * \param dcomp starting component of the destination MF
485 : * \param ncomp number of components
486 : * \param cgeom Geometry for the coarse level
487 : * \param fgeom Geometry for the fine level
488 : * \param cbc functor for physical boundaries on the coarse level
489 : * \param cbccomp starting component for cbc
490 : * \param fbc functor for physical boundaries on the fine level
491 : * \param fbccomp starting component for fbc
492 : * \param ratio refinement ratio
493 : * \param mapper spatial interpolater
494 : * \param bcs boundary types for each component. We need this because some interpolaters need it.
495 : * \param bcscomp starting component for bcs
496 : * \param pre_interp pre-interpolation hook
497 : * \param post_interp post-interpolation hook
498 : */
499 : template <typename MF, typename BC, typename Interp, typename PreInterpHook, typename PostInterpHook>
500 : std::enable_if_t<IsFabArray<MF>::value>
501 : FillPatchTwoLevels (MF& mf, Real time,
502 : const EB2::IndexSpace& index_space,
503 : const Vector<MF*>& cmf, const Vector<Real>& ct,
504 : const Vector<MF*>& fmf, const Vector<Real>& ft,
505 : int scomp, int dcomp, int ncomp,
506 : const Geometry& cgeom, const Geometry& fgeom,
507 : BC& cbc, int cbccomp,
508 : BC& fbc, int fbccomp,
509 : const IntVect& ratio,
510 : Interp* mapper,
511 : const Vector<BCRec>& bcs, int bcscomp,
512 : const PreInterpHook& pre_interp,
513 : const PostInterpHook& post_interp);
514 : #endif
515 :
516 : /**
517 : * \brief Fill with interpolation of coarse level data
518 : *
519 : * No ghost cells of the destination MF are filled.
520 : *
521 : * \tparam MF the MultiFab/FabArray type
522 : * \tparam BC functor for filling physical boundaries
523 : * \tparam Interp spatial interpolater
524 : * \tparam PreInterpHook pre-interpolation hook
525 : * \tparam PostInterpHook post-interpolation hook
526 : *
527 : * \param mf destination MF on the fine level
528 : * \param time time associated with mf
529 : * \param cmf source MF on the coarse level
530 : * \param scomp starting component of the source MF
531 : * \param dcomp starting component of the destination MF
532 : * \param ncomp number of components
533 : * \param cgeom Geometry for the coarse level
534 : * \param fgeom Geometry for the fine level
535 : * \param cbc functor for physical boundaries on the coarse level
536 : * \param cbccomp starting component for cbc
537 : * \param fbc functor for physical boundaries on the fine level
538 : * \param fbccomp starting component for fbc
539 : * \param ratio refinement ratio
540 : * \param mapper spatial interpolater
541 : * \param bcs boundary types for each component. We need this because some interpolaters need it.
542 : * \param bcscomp starting component for bcs
543 : * \param pre_interp pre-interpolation hook
544 : * \param post_interp post-interpolation hook
545 : */
546 : template <typename MF, typename BC, typename Interp,
547 : typename PreInterpHook=NullInterpHook<typename MF::FABType::value_type>,
548 : typename PostInterpHook=NullInterpHook<typename MF::FABType::value_type> >
549 : std::enable_if_t<IsFabArray<MF>::value>
550 : InterpFromCoarseLevel (MF& mf, Real time,
551 : const MF& cmf, int scomp, int dcomp, int ncomp,
552 : const Geometry& cgeom, const Geometry& fgeom,
553 : BC& cbc, int cbccomp,
554 : BC& fbc, int fbccomp,
555 : const IntVect& ratio,
556 : Interp* mapper,
557 : const Vector<BCRec>& bcs, int bcscomp,
558 : const PreInterpHook& pre_interp = {},
559 : const PostInterpHook& post_interp = {});
560 :
561 : /**
562 : * \brief Fill with interpolation of coarse level data
563 : *
564 : * \tparam MF the MultiFab/FabArray type
565 : * \tparam BC functor for filling physical boundaries
566 : * \tparam Interp spatial interpolater
567 : * \tparam PreInterpHook pre-interpolation hook
568 : * \tparam PostInterpHook post-interpolation hook
569 : *
570 : * \param mf destination MF on the fine level
571 : * \param nghost number of ghost cells of mf needed to be filled
572 : * \param time time associated with mf
573 : * \param cmf source MF on the coarse level
574 : * \param scomp starting component of the source MF
575 : * \param dcomp starting component of the destination MF
576 : * \param ncomp number of components
577 : * \param cgeom Geometry for the coarse level
578 : * \param fgeom Geometry for the fine level
579 : * \param cbc functor for physical boundaries on the coarse level
580 : * \param cbccomp starting component for cbc
581 : * \param fbc functor for physical boundaries on the fine level
582 : * \param fbccomp starting component for fbc
583 : * \param ratio refinement ratio
584 : * \param mapper spatial interpolater
585 : * \param bcs boundary types for each component. We need this because some interpolaters need it.
586 : * \param bcscomp starting component for bcs
587 : * \param pre_interp pre-interpolation hook
588 : * \param post_interp post-interpolation hook
589 : */
590 : template <typename MF, typename BC, typename Interp,
591 : typename PreInterpHook=NullInterpHook<typename MF::FABType::value_type>,
592 : typename PostInterpHook=NullInterpHook<typename MF::FABType::value_type> >
593 : std::enable_if_t<IsFabArray<MF>::value>
594 : InterpFromCoarseLevel (MF& mf, IntVect const& nghost, Real time,
595 : const MF& cmf, int scomp, int dcomp, int ncomp,
596 : const Geometry& cgeom, const Geometry& fgeom,
597 : BC& cbc, int cbccomp,
598 : BC& fbc, int fbccomp,
599 : const IntVect& ratio,
600 : Interp* mapper,
601 : const Vector<BCRec>& bcs, int bcscomp,
602 : const PreInterpHook& pre_interp = {},
603 : const PostInterpHook& post_interp = {});
604 :
605 : /**
606 : * \brief Fill face variables with data from the coarse level.
607 : * Sometimes, we need to fillpatch all AMREX_SPACEDIM face MultiFabs
608 : * togother to satisfy certain constraint such as divergence preserving.
609 : *
610 : * \tparam MF the MultiFab/FabArray type
611 : * \tparam BC functor for filling physical boundaries
612 : * \tparam Interp spatial interpolater
613 : * \tparam PreInterpHook pre-interpolation hook
614 : * \tparam PostInterpHook post-interpolation hook
615 : *
616 : * \param mf destination MFs on the fine level
617 : * \param time time associated with mf
618 : * \param cmf source MFs on the coarse level
619 : * \param scomp starting component of the source MFs
620 : * \param dcomp starting component of the destination MFs
621 : * \param ncomp number of components
622 : * \param cgeom Geometry for the coarse level
623 : * \param fgeom Geometry for the fine level
624 : * \param cbc functor for physical boundaries on the coarse level
625 : * \param cbccomp starting component for cbc
626 : * \param fbc functor for physical boundaries on the fine level
627 : * \param fbccomp starting component for fbc
628 : * \param ratio refinement ratio
629 : * \param mapper spatial interpolater
630 : * \param bcs boundary types for each component. We need this because some interpolaters need it.
631 : * \param bcscomp starting component for bcs
632 : * \param pre_interp pre-interpolation hook
633 : * \param post_interp post-interpolation hook
634 : */
635 : template <typename MF, typename BC, typename Interp,
636 : typename PreInterpHook=NullInterpHook<typename MF::FABType::value_type>,
637 : typename PostInterpHook=NullInterpHook<typename MF::FABType::value_type> >
638 : std::enable_if_t<IsFabArray<MF>::value>
639 : InterpFromCoarseLevel (Array<MF*, AMREX_SPACEDIM> const& mf, Real time,
640 : const Array<MF*, AMREX_SPACEDIM>& cmf, int scomp, int dcomp, int ncomp,
641 : const Geometry& cgeom, const Geometry& fgeom,
642 : Array<BC, AMREX_SPACEDIM>& cbc, int cbccomp,
643 : Array<BC, AMREX_SPACEDIM>& fbc, int fbccomp,
644 : const IntVect& ratio,
645 : Interp* mapper,
646 : const Array<Vector<BCRec>, AMREX_SPACEDIM>& bcs, int bcscomp,
647 : const PreInterpHook& pre_interp = {},
648 : const PostInterpHook& post_interp = {});
649 :
650 : /**
651 : * \brief Fill face variables with data from the coarse level.
652 : * Sometimes, we need to fillpatch all AMREX_SPACEDIM face MultiFabs
653 : * togother to satisfy certain constraint such as divergence preserving.
654 : *
655 : * \tparam MF the MultiFab/FabArray type
656 : * \tparam BC functor for filling physical boundaries
657 : * \tparam Interp spatial interpolater
658 : * \tparam PreInterpHook pre-interpolation hook
659 : * \tparam PostInterpHook post-interpolation hook
660 : *
661 : * \param mf destination MFs on the fine level
662 : * \param nghost number of ghost cells of mf needed to be filled
663 : * \param time time associated with mf
664 : * \param cmf source MFs on the coarse level
665 : * \param scomp starting component of the source MFs
666 : * \param dcomp starting component of the destination MFs
667 : * \param ncomp number of components
668 : * \param cgeom Geometry for the coarse level
669 : * \param fgeom Geometry for the fine level
670 : * \param cbc functor for physical boundaries on the coarse level
671 : * \param cbccomp starting component for cbc
672 : * \param fbc functor for physical boundaries on the fine level
673 : * \param fbccomp starting component for fbc
674 : * \param ratio refinement ratio
675 : * \param mapper spatial interpolater
676 : * \param bcs boundary types for each component. We need this because some interpolaters need it.
677 : * \param bcscomp starting component for bcs
678 : * \param pre_interp pre-interpolation hook
679 : * \param post_interp post-interpolation hook
680 : */
681 : template <typename MF, typename BC, typename Interp,
682 : typename PreInterpHook=NullInterpHook<typename MF::FABType::value_type>,
683 : typename PostInterpHook=NullInterpHook<typename MF::FABType::value_type> >
684 : std::enable_if_t<IsFabArray<MF>::value>
685 : InterpFromCoarseLevel (Array<MF*, AMREX_SPACEDIM> const& mf, IntVect const& nghost, Real time,
686 : const Array<MF*, AMREX_SPACEDIM>& cmf, int scomp, int dcomp, int ncomp,
687 : const Geometry& cgeom, const Geometry& fgeom,
688 : Array<BC, AMREX_SPACEDIM>& cbc, int cbccomp,
689 : Array<BC, AMREX_SPACEDIM>& fbc, int fbccomp,
690 : const IntVect& ratio,
691 : Interp* mapper,
692 : const Array<Vector<BCRec>, AMREX_SPACEDIM>& bcs, int bcscomp,
693 : const PreInterpHook& pre_interp = {},
694 : const PostInterpHook& post_interp = {});
695 :
696 : #ifndef BL_NO_FORT
697 : enum InterpEM_t { InterpE, InterpB};
698 :
699 : void InterpCrseFineBndryEMfield (InterpEM_t interp_type,
700 : const Array<MultiFab,AMREX_SPACEDIM>& crse,
701 : Array<MultiFab,AMREX_SPACEDIM>& fine,
702 : const Geometry& cgeom, const Geometry& fgeom,
703 : int ref_ratio);
704 :
705 : void InterpCrseFineBndryEMfield (InterpEM_t interp_type,
706 : const Array<MultiFab const*,AMREX_SPACEDIM>& crse,
707 : const Array<MultiFab*,AMREX_SPACEDIM>& fine,
708 : const Geometry& cgeom, const Geometry& fgeom,
709 : int ref_ratio);
710 : #endif
711 :
712 : /**
713 : * \brief FillPatch with data from AMR levels.
714 : *
715 : * First, we try to fill the destination MultiFab/FabArray with this
716 : * level's data if it's available. For the unfilled region, we try to
717 : * fill with the coarse level below if it's available. Even coarser
718 : * levels will be used if necessary till all regions are filled. This
719 : * function is more expensive than FillPatchTwoLevels. So if one knows
720 : * FillPatchTwoLevels can do the job because grids are properly nested,
721 : * this function should be avoided.
722 : *
723 : * \tparam MF the MultiFab/FabArray type
724 : * \tparam BC functor for filling physical boundaries
725 : * \tparam Interp spatial interpolater
726 : *
727 : * \param mf destination MF
728 : * \param level AMR level associated with mf
729 : * \param nghost number of ghost cells of mf needed to be filled
730 : * \param time time associated with mf
731 : * \param smf source MFs. The outer Vector is for AMR levels, whereas the inner
732 : * Vector is for data at various times. It is not an error if the
733 : * level for the destination MF is finer than data in smf (i.e.,
734 : * `level >= smf.size()`).
735 : * \param st times associated smf
736 : * \param scomp starting component of the source MFs
737 : * \param dcomp starting component of the destination MF
738 : * \param ncomp number of components
739 : * \param geom Geometry objects for AMR levels. The size must be big enough
740 : * such that `level < geom.size()`.
741 : * \param bc functors for physical boundaries on AMR levels. The size must be
742 : * big enough such that `level < bc.size()`.
743 : * \param bccomp starting component for bc
744 : * \param ratio refinement ratio for AMR levels. The size must be big enough
745 : * such that `level < bc.size()-1`.
746 : * \param mapper spatial interpolater
747 : * \param bcr boundary types for each component. We need this because some interpolaters need it.
748 : * \param bcrcomp starting component for bcr
749 : */
750 : template <typename MF, typename BC, typename Interp>
751 : std::enable_if_t<IsFabArray<MF>::value>
752 : FillPatchNLevels (MF& mf, int level, const IntVect& nghost, Real time,
753 : const Vector<Vector<MF*>>& smf, const Vector<Vector<Real>>& st,
754 : int scomp, int dcomp, int ncomp,
755 : const Vector<Geometry>& geom,
756 : Vector<BC>& bc, int bccomp,
757 : const Vector<IntVect>& ratio,
758 : Interp* mapper,
759 : const Vector<BCRec>& bcr, int bcrcomp);
760 :
761 : }
762 :
763 : #include <AMReX_FillPatchUtil_I.H>
764 :
765 : #endif
|