Line data Source code
1 :
2 : #ifndef AMREX_INTERPBNDRYDATA_H_
3 : #define AMREX_INTERPBNDRYDATA_H_
4 : #include <AMReX_Config.H>
5 :
6 : #include <AMReX_BCRec.H>
7 : #include <AMReX_BndryData.H>
8 : #include <AMReX_InterpBndryData_K.H>
9 :
10 : namespace amrex {
11 :
12 : /**
13 : * \brief An InterpBndryData object adds to a BndryData object the ability to
14 : manipulate and set the data stored in the boundary cells.
15 :
16 : The InterpBndryData class is a class derived from
17 : BndryData. It is intended to provide a more physical method for
18 : filling boundary-related data. Boundary values in a BndryData object
19 : are stored in FabSets around each grid in the domain, and the
20 : InterpBndryData class provides a mechanism for filling these FabSets,
21 : consistent with AMR-like numerical discretizations. When asked to
22 : set its boundary values, an InterpBndryData object:
23 :
24 : Fills with physical boundary values if the FAB is on the
25 : domain boundary - the corresponding values are presumed to be
26 : stored in the ghost cells of a MultiFab given to the boundary filling
27 : routine
28 :
29 : Fills on intersection with data from the VALID region of the
30 : input MultiFab, ensuring that adjacent FABs use consistent data at
31 : their intersection, and otherwise,
32 :
33 : Fills with values interpolated from a coarser FAB that
34 : bounds the cells that do not meet the above two criteria
35 : */
36 : template <typename MF>
37 : class InterpBndryDataT
38 : :
39 : public BndryDataT<MF>
40 : {
41 : public:
42 :
43 : using value_type = typename MF::value_type;
44 :
45 : /**
46 : * \brief default constructor
47 : */
48 : InterpBndryDataT () noexcept = default;
49 :
50 : /**
51 : * \brief constructor for given BoxArray, etc
52 : *
53 : * \param _grids
54 : * \param _dmap
55 : * \param _ncomp
56 : * \param geom
57 : */
58 : InterpBndryDataT (const BoxArray& _grids,
59 : const DistributionMapping& _dmap,
60 : int _ncomp,
61 : const Geometry& geom);
62 :
63 0 : ~InterpBndryDataT () = default;
64 :
65 : InterpBndryDataT (const InterpBndryDataT<MF>& rhs) = delete;
66 : InterpBndryDataT (InterpBndryDataT<MF>&& rhs) = delete;
67 : InterpBndryDataT<MF>& operator= (const InterpBndryDataT<MF>& rhs) = delete;
68 : InterpBndryDataT<MF>& operator= (InterpBndryDataT<MF>&& rhs) = delete;
69 :
70 : /**
71 : * \brief Set bndry values at physical boundaries
72 : *
73 : * \param mf MF containing physical boundary values
74 : * \param mf_start starting component of the data in MF
75 : * \param bnd_start starting component in this boundary register
76 : * \param num_comp number of components
77 : */
78 : void setPhysBndryValues (const MF& mf, int mf_start, int bnd_start, int num_comp);
79 :
80 : /**
81 : * \briefSset bndry values at coarse/fine and physical boundaries
82 : *
83 : * \param crse BndryRegister storing coarse level data
84 : * \param c_start starting component of the coarse data
85 : * \param fine MF containing physical boundary values
86 : * \param f_start starting component of the fine data
87 : * \param bnd_start starting component in this InterpBndryData
88 : * \param num_comp number of component
89 : * \param ratio refinement ratio
90 : * \param max_order interpolation order
91 : */
92 : void setBndryValues (BndryRegisterT<MF> const& crse, int c_start, const MF& fine, int f_start,
93 : int bnd_start, int num_comp, const IntVect& ratio,
94 : int max_order = IBD_max_order_DEF, int max_width = 2);
95 :
96 : /**
97 : * \brief Update boundary values at coarse/fine boundaries
98 : *
99 : * \param crse BndryRegister storing coarse level data
100 : * \param c_start starting component of the coarse data
101 : * \param bnd_start starting component in this InterpBndryData
102 : * \param num_comp number of component
103 : * \param ratio refinement ratio
104 : * \param max_order interpolation order
105 : */
106 : void updateBndryValues (BndryRegisterT<MF>& crse, int c_start, int bnd_start, int num_comp,
107 : const IntVect& ratio, int max_order = IBD_max_order_DEF, int max_width = 2);
108 :
109 : //! Set boundary values to zero
110 : void setHomogValues ();
111 :
112 : // For sliding parabolic interp in bdfuncs
113 : static constexpr int IBD_max_order_DEF = 3;
114 :
115 : };
116 :
117 : template <typename MF>
118 : InterpBndryDataT<MF>::InterpBndryDataT (const BoxArray& _grids,
119 : const DistributionMapping& _dmap,
120 : int _ncomp,
121 : const Geometry& _geom)
122 : :
123 : BndryDataT<MF>(_grids,_dmap,_ncomp,_geom)
124 : {}
125 :
126 : template <typename MF>
127 : void
128 : InterpBndryDataT<MF>::setPhysBndryValues (const MF& mf, int mf_start, int bnd_start, int num_comp)
129 : {
130 : AMREX_ASSERT(this->grids == mf.boxArray());
131 :
132 : const Box& fine_domain = this->geom.Domain();
133 :
134 : #ifdef AMREX_USE_OMP
135 : #pragma omp parallel if (Gpu::notInLaunchRegion())
136 : #endif
137 : for (MFIter mfi(mf,MFItInfo().SetDynamic(true)); mfi.isValid(); ++mfi)
138 : {
139 : const Box& bx = mfi.validbox();
140 : for (OrientationIter fi; fi; ++fi) {
141 : const Orientation face = fi();
142 : if (bx[face] == fine_domain[face] && !(this->geom.isPeriodic(face.coordDir())))
143 : {
144 : // Physical bndry, copy from grid.
145 : auto & bnd_fab = this->bndry[face][mfi];
146 : auto const& src_fab = mf[mfi];
147 : auto const& bnd_array = bnd_fab.array();
148 : auto const& src_array = src_fab.const_array();
149 : const Box& b = src_fab.box() & bnd_fab.box();
150 : AMREX_HOST_DEVICE_PARALLEL_FOR_4D ( b, num_comp, i, j, k, n,
151 : {
152 : bnd_array(i,j,k,n+bnd_start) = src_array(i,j,k,n+mf_start);
153 : });
154 : }
155 : }
156 : }
157 : }
158 :
159 : template <typename MF>
160 : void
161 : InterpBndryDataT<MF>::setBndryValues (BndryRegisterT<MF> const& crse, int c_start, const MF& fine,
162 : int f_start, int bnd_start, int num_comp, const IntVect& ratio,
163 : int max_order, int max_width)
164 : {
165 : AMREX_ASSERT(this->grids == fine.boxArray());
166 :
167 : const Box& fine_domain = this->geom.Domain();
168 :
169 : if (max_order==3 || max_order==1)
170 : {
171 : MFItInfo info;
172 : if (Gpu::notInLaunchRegion()) { info.SetDynamic(true); }
173 : #ifdef AMREX_USE_OMP
174 : #pragma omp parallel if (Gpu::notInLaunchRegion())
175 : #endif
176 : for (MFIter mfi(fine,info); mfi.isValid(); ++mfi)
177 : {
178 : const Box& fine_bx = mfi.validbox();
179 : for (OrientationIter fi; fi; ++fi) {
180 : const Orientation face = fi();
181 : const int dir = face.coordDir();
182 : if (fine_bx[face] != fine_domain[face] || this->geom.isPeriodic(dir))
183 : { // coarse/fine boundary: interpolate from coarse data stored in BndryRegister crse
184 : auto const& crse_array = crse[face].const_array(mfi);
185 : auto const& bdry_array = this->bndry[face].array(mfi);
186 : Box const& b = this->bndry[face][mfi].box();
187 : const auto rr = ratio.dim3();
188 :
189 : if (max_order == 1)
190 : {
191 : AMREX_HOST_DEVICE_FOR_4D ( b, num_comp, it, jt, kt, n,
192 : {
193 : interpbndrydata_o1(it,jt,kt,n,bdry_array,bnd_start,
194 : crse_array,c_start,rr);
195 : });
196 : }
197 : else
198 : {
199 : auto const& mask_array = this->masks[face].const_array(mfi);
200 : int is_not_covered = BndryData::not_covered;
201 : switch (dir) {
202 : case 0:
203 : {
204 : AMREX_HOST_DEVICE_FOR_4D ( b, num_comp, it, jt, kt, n,
205 : {
206 : interpbndrydata_x_o3(it,jt,kt,n,bdry_array,bnd_start,
207 : crse_array,c_start,rr,
208 : mask_array, is_not_covered, max_width);
209 : });
210 : break;
211 : }
212 : #if (AMREX_SPACEDIM >= 2)
213 : case 1:
214 : {
215 : AMREX_HOST_DEVICE_FOR_4D ( b, num_comp, it, jt, kt, n,
216 : {
217 : interpbndrydata_y_o3(it,jt,kt,n,bdry_array,bnd_start,
218 : crse_array,c_start,rr,
219 : mask_array, is_not_covered, max_width);
220 : });
221 : break;
222 : }
223 : #if (AMREX_SPACEDIM == 3)
224 : case 2:
225 : {
226 : AMREX_HOST_DEVICE_FOR_4D ( b, num_comp, it, jt, kt, n,
227 : {
228 : interpbndrydata_z_o3(it,jt,kt,n,bdry_array,bnd_start,
229 : crse_array,c_start,rr,
230 : mask_array, is_not_covered, max_width);
231 : });
232 : break;
233 : }
234 : #endif
235 : #endif
236 : default: {}
237 : }
238 : }
239 : }
240 : else if (fine.defined(mfi)) {
241 : // Physical bndry
242 : auto & bnd_fab = this->bndry[face][mfi];
243 : auto const& src_fab = fine[mfi];
244 : auto const& bnd_array = bnd_fab.array();
245 : auto const& src_array = src_fab.const_array();
246 : const Box& b = bnd_fab.box() & src_fab.box();
247 : AMREX_HOST_DEVICE_FOR_4D ( b, num_comp, ii, jj, kk, nn,
248 : {
249 : bnd_array(ii,jj,kk,nn+bnd_start) = src_array(ii,jj,kk,nn+f_start);
250 : });
251 : }
252 : }
253 : }
254 : }
255 : else
256 : {
257 : amrex::Abort("InterpBndryDataT<MF>::setBndryValues supports only max_order=1 or 3");
258 : }
259 : }
260 :
261 : template <typename MF>
262 : void
263 : InterpBndryDataT<MF>::updateBndryValues (BndryRegisterT<MF>& crse, int c_start, int bnd_start,
264 : int num_comp, const IntVect& ratio, int max_order, int max_width)
265 : {
266 : MF foo(this->grids, this->bndry[0].DistributionMap(), 1, num_comp, MFInfo().SetAlloc(false));
267 : setBndryValues(crse, c_start, foo, 0, bnd_start, num_comp, ratio, max_order, max_width);
268 : }
269 :
270 : template <typename MF>
271 : void
272 : InterpBndryDataT<MF>::setHomogValues ()
273 : {
274 : this->setVal(0.);
275 : }
276 :
277 : using InterpBndryData = InterpBndryDataT<MultiFab>;
278 : using fInterpBndryData = InterpBndryDataT<fMultiFab>;
279 :
280 : }
281 :
282 : #endif /*AMREX_INTERPBNDRYDATA_H_*/
|