Line data Source code
1 :
2 : #ifndef BL_IMULTIFAB_H
3 : #define BL_IMULTIFAB_H
4 : #include <AMReX_Config.H>
5 :
6 : #include <AMReX_BLassert.H>
7 : #include <AMReX_IArrayBox.H>
8 : #include <AMReX_FabArray.H>
9 : #include <AMReX_FabArrayUtility.H>
10 : #include <AMReX_Geometry.H>
11 : #include <memory>
12 :
13 : namespace amrex {
14 :
15 : /*
16 : * A Collection of IArrayBoxes
17 : *
18 : * The iMultiFab class is publicly derived from the
19 : * FabArray<int,IArrayBox> class. It is a collection (stored as an array) of
20 : * IArrayBoxes useful for storing integer data on a domain defined by
21 : * a union of rectangular regions embedded in a uniform index space. The
22 : * iMultiFab class extends the function of the underlying FabArray class just
23 : * as the IArrayBox class extends the function of BaseFab<int>. Additional
24 : * member functions are defined for I/O and simple arithmetic operations on
25 : * these aggregate objects.
26 : *
27 : * This class does NOT provide a copy constructor or assignment operator.
28 : */
29 : class iMultiFab
30 : :
31 : public FabArray<IArrayBox>
32 : {
33 : public:
34 :
35 : /**
36 : * \brief Constructs an empty iMultiFab. Data can be defined at a later
37 : * time using the define member functions inherited
38 : * from FabArray.
39 : */
40 : iMultiFab () noexcept = default;
41 :
42 : /**
43 : * \brief Constructs an empty iMultiFab. Data can be defined at a later
44 : * time using the define member functions inherited from FabArray. If
45 : * `define` is called later with a nullptr as MFInfo's arena, the default
46 : * Arena `a` will be used. If the arena in MFInfo is not a nullptr, the
47 : * MFInfo's arena will be used.
48 : */
49 : explicit iMultiFab (Arena* a) noexcept;
50 :
51 : /**
52 : * \brief Constructs a iMultiFab with a valid region defined by bxs and
53 : * a region of definition defined by the grow factor ngrow.
54 : */
55 : iMultiFab (const BoxArray& bxs,
56 : const DistributionMapping& dm,
57 : int ncomp,
58 : int ngrow,
59 : #ifdef AMREX_STRICT_MODE
60 : const MFInfo& info,
61 : const FabFactory<IArrayBox>& factory);
62 : #else
63 : const MFInfo& info = MFInfo(),
64 : const FabFactory<IArrayBox>& factory = DefaultFabFactory<IArrayBox>());
65 : #endif
66 :
67 : iMultiFab (const BoxArray& bxs,
68 : const DistributionMapping& dm,
69 : int ncomp,
70 : const IntVect& ngrow,
71 : #ifdef AMREX_STRICT_MODE
72 : const MFInfo& info,
73 : const FabFactory<IArrayBox>& factory);
74 : #else
75 : const MFInfo& info = MFInfo(),
76 : const FabFactory<IArrayBox>& factory = DefaultFabFactory<IArrayBox>());
77 : #endif
78 :
79 : /**
80 : * \brief Make an alias iMultiFab. maketype must be
81 : * amrex::make_alias. scomp is the starting component of the
82 : * alias and ncomp is the number of components in the new aliasing
83 : * iMultiFab.
84 : *
85 : * \param rhs
86 : * \param maketype
87 : * \param scomp
88 : * \param ncomp
89 : */
90 : iMultiFab (const iMultiFab& rhs, MakeType maketype, int scomp, int ncomp);
91 :
92 1587 : ~iMultiFab () = default;
93 :
94 : iMultiFab (iMultiFab&& rhs) noexcept = default;
95 : iMultiFab& operator= (iMultiFab&& rhs) noexcept = default;
96 :
97 : iMultiFab (const iMultiFab& rhs) = delete;
98 : iMultiFab& operator= (const iMultiFab& rhs) = delete;
99 :
100 : iMultiFab& operator= (int r);
101 :
102 : /**
103 : * \brief Returns the minimum value contained in component comp of the
104 : * iMultiFab. The parameter nghost determines the number of
105 : * boundary cells to search for the minimum. The default is to
106 : * search only the valid regions of the IArrayBoxes.
107 : *
108 : * \param comp
109 : * \param nghost
110 : * \param local
111 : */
112 : [[nodiscard]] int min (int comp,
113 : int nghost = 0,
114 : bool local = false) const;
115 :
116 : /**
117 : * \brief Identical to the previous min() function, but confines its
118 : * search to intersection of Box b and the iMultiFab.
119 : *
120 : * \param b
121 : * \param comp
122 : * \param nghost
123 : * \param local
124 : */
125 : [[nodiscard]] int min (const Box& region,
126 : int comp,
127 : int nghost = 0,
128 : bool local = false) const;
129 :
130 : /**
131 : * \brief Returns the maximum value contained in component comp of the
132 : * iMultiFab. The parameter nghost determines the number of
133 : * boundary cells to search for the maximum. The default is to
134 : * search only the valid regions of the IArrayBoxes.
135 : *
136 : * \param comp
137 : * \param nghost
138 : * \param local
139 : */
140 : [[nodiscard]] int max (int comp,
141 : int nghost = 0,
142 : bool local = false) const;
143 :
144 : /**
145 : * \brief Identical to the previous max() function, but confines its
146 : * search to intersection of Box b and the iMultiFab.
147 : *
148 : * \param b
149 : * \param comp
150 : * \param nghost
151 : * \param local
152 : */
153 : [[nodiscard]] int max (const Box& region,
154 : int comp,
155 : int nghost = 0,
156 : bool local = false) const;
157 :
158 : /**
159 : * \brief Returns the sum in component comp
160 : *
161 : * \param comp
162 : * \param nghost
163 : * \param local
164 : */
165 : [[nodiscard]] Long sum (int comp, int nghost = 0, bool local = false) const;
166 :
167 : /**
168 : * \brief Adds the scalar value val to the value of each cell in the
169 : * specified subregion of the iMultiFab. The subregion consists
170 : * of the num_comp components starting at component comp.
171 : * The value of nghost specifies the number of cells in the
172 : * boundary region of each IArrayBox in the subregion that should
173 : * be modified.
174 : *
175 : * \param val
176 : * \param comp
177 : * \param num_comp
178 : * \param nghost
179 : */
180 : void plus (int val,
181 : int comp,
182 : int num_comp,
183 : int nghost = 0);
184 :
185 : /**
186 : * \brief Identical to the previous version of plus(), with the
187 : * restriction that the subregion is further constrained to
188 : * the intersection with Box region.
189 : *
190 : * \param val
191 : * \param region
192 : * \param comp
193 : * \param num_comp
194 : * \param nghost
195 : */
196 : void plus (int val,
197 : const Box& region,
198 : int comp,
199 : int num_comp,
200 : int nghost = 0);
201 :
202 : /**
203 : * \brief Adds the scalar value val to the value of each cell in the
204 : * valid region of each component of the iMultiFab. The value
205 : * of nghost specifies the number of cells in the boundary
206 : * region that should be modified.
207 : *
208 : * \param val
209 : * \param nghost
210 : */
211 : void plus (int val,
212 : int nghost);
213 :
214 : /**
215 : * \brief Adds the scalar value val to the value of each cell in the
216 : * valid region of each component of the iMultiFab, that also
217 : * intersects the Box region. The value of nghost specifies the
218 : * number of cells in the boundary region of each IArrayBox in
219 : * the subregion that should be modified.
220 : *
221 : * \param val
222 : * \param region
223 : * \param nghost
224 : */
225 : void plus (int val,
226 : const Box& region,
227 : int nghost);
228 :
229 : /**
230 : * \brief Scales the value of each cell in the specified subregion of the
231 : * iMultiFab by the scalar val (a[i] \<- a[i]*val). The subregion
232 : * consists of the num_comp components starting at component comp.
233 : * The value of nghost specifies the number of cells in the
234 : * boundary region of each IArrayBox in the subregion that should
235 : * be modified.
236 : *
237 : * \param val
238 : * \param comp
239 : * \param num_comp
240 : * \param nghost
241 : */
242 : void mult (int val,
243 : int comp,
244 : int num_comp,
245 : int nghost = 0);
246 :
247 : /**
248 : * \brief Identical to the previous version of mult(), with the
249 : * restriction that the subregion is further constrained to the
250 : * intersection with Box region. The value of nghost specifies the
251 : * number of cells in the boundary region of each IArrayBox in
252 : * the subregion that should be modified.
253 : *
254 : * \param val
255 : * \param region
256 : * \param comp
257 : * \param num_comp
258 : * \param nghost
259 : */
260 : void mult (int val,
261 : const Box& region,
262 : int comp,
263 : int num_comp,
264 : int nghost = 0);
265 :
266 : /**
267 : * \brief Scales the value of each cell in the valid region of each
268 : * component of the iMultiFab by the scalar val (a[i] \<- a[i]*val).
269 : * The value of nghost specifies the number of cells in the
270 : * boundary region that should be modified.
271 : *
272 : * \param val
273 : * \param nghost
274 : */
275 : void mult (int val,
276 : int nghost = 0);
277 :
278 : /**
279 : * \brief Scales the value of each cell in the valid region of each
280 : * component of the iMultiFab by the scalar val (a[i] \<- a[i]*val),
281 : * that also intersects the Box region. The value of nghost
282 : * specifies the number of cells in the boundary region of each
283 : * IArrayBox in the subregion that should be modified.
284 : *
285 : * \param val
286 : * \param region
287 : * \param nghost
288 : */
289 : void mult (int val,
290 : const Box& region,
291 : int nghost = 0);
292 :
293 : /**
294 : * \brief Negates the value of each cell in the specified subregion of
295 : * the iMultiFab. The subregion consists of the num_comp
296 : * components starting at component comp. The value of nghost
297 : * specifies the number of cells in the boundary region of each
298 : * IArrayBox in the subregion that should be modified.
299 : *
300 : * \param comp
301 : * \param num_comp
302 : * \param nghost
303 : */
304 : void negate (int comp,
305 : int num_comp,
306 : int nghost = 0);
307 :
308 : /**
309 : * \brief Identical to the previous version of negate(), with the
310 : * restriction that the subregion is further constrained to
311 : * the intersection with Box region.
312 : *
313 : * \param region
314 : * \param comp
315 : * \param num_comp
316 : * \param nghost
317 : */
318 : void negate (const Box& region,
319 : int comp,
320 : int num_comp,
321 : int nghost = 0);
322 :
323 : /**
324 : * \brief Negates the value of each cell in the valid region of
325 : * the iMultiFab. The value of nghost specifies the number of
326 : * cells in the boundary region that should be modified.
327 : *
328 : * \param nghost
329 : */
330 : void negate (int nghost = 0);
331 :
332 : /**
333 : * \brief Negates the value of each cell in the valid region of
334 : * the iMultiFab that also intersects the Box region. The value
335 : * of nghost specifies the number of cells in the boundary region
336 : * that should be modified.
337 : *
338 : * \param region
339 : * \param nghost
340 : */
341 : void negate (const Box& region,
342 : int nghost = 0);
343 :
344 : [[nodiscard]] IntVect minIndex (int comp,
345 : int nghost = 0) const;
346 :
347 : [[nodiscard]] IntVect maxIndex (int comp,
348 : int nghost = 0) const;
349 :
350 : /**
351 : * \brief This function adds the values of the cells in mf to the corresponding
352 : * cells of this iMultiFab. mf is required to have the same BoxArray or
353 : * "valid region" as this iMultiFab. The addition is done only to num_comp
354 : * components, starting with component number strt_comp. The parameter
355 : * nghost specifies the number of boundary cells that will be modified.
356 : * If nghost == 0, only the valid region of each IArrayBox will be
357 : * modified.
358 : *
359 : * \param mf
360 : * \param strt_comp
361 : * \param num_comp
362 : * \param nghost
363 : */
364 : void plus (const iMultiFab& mf,
365 : int strt_comp,
366 : int num_comp,
367 : int nghost);
368 :
369 : /**
370 : * \brief This function subtracts the values of the cells in mf from the
371 : * corresponding cells of this iMultiFab. mf is required to have the
372 : * same BoxArray or "valid region" as this iMultiFab. The subtraction is
373 : * done only to num_comp components, starting with component number
374 : * strt_comp. The parameter nghost specifies the number of boundary
375 : * cells that will be modified. If nghost == 0, only the valid region of
376 : * each IArrayBox will be modified.
377 : *
378 : * \param mf
379 : * \param strt_comp
380 : * \param num_comp
381 : * \param nghost
382 : */
383 : void minus (const iMultiFab& mf,
384 : int strt_comp,
385 : int num_comp,
386 : int nghost);
387 :
388 : /**
389 : * \brief This function divides the values of the cells in mf from the
390 : * corresponding cells of this iMultiFab. mf is required to have the
391 : * same BoxArray or "valid region" as this iMultiFab. The division is
392 : * done only to num_comp components, starting with component number
393 : * strt_comp. The parameter nghost specifies the number of boundary
394 : * cells that will be modified. If nghost == 0, only the valid region of
395 : * each IArrayBox will be modified. Note, nothing is done to protect
396 : * against divide by zero.
397 : *
398 : * \param mf
399 : * \param strt_comp
400 : * \param num_comp
401 : * \param nghost
402 : */
403 : void divide (const iMultiFab& mf,
404 : int strt_comp,
405 : int num_comp,
406 : int nghost);
407 :
408 : /**
409 : * \brief Add src to dst including nghost ghost cells.
410 : * The two iMultiFabs MUST have the same underlying BoxArray.
411 : *
412 : * \param dst
413 : * \param src
414 : * \param srccomp
415 : * \param dstcomp
416 : * \param numcomp
417 : * \param nghost
418 : */
419 : static void Add (iMultiFab& dst,
420 : const iMultiFab& src,
421 : int srccomp,
422 : int dstcomp,
423 : int numcomp,
424 : int nghost);
425 :
426 : /**
427 : * \brief Copy from src to dst including nghost ghost cells.
428 : * The two iMultiFabs MUST have the same underlying BoxArray.
429 : *
430 : * \param dst
431 : * \param src
432 : * \param srccomp
433 : * \param dstcomp
434 : * \param numcomp
435 : * \param nghost
436 : */
437 : static void Copy (iMultiFab& dst,
438 : const iMultiFab& src,
439 : int srccomp,
440 : int dstcomp,
441 : int numcomp,
442 : int nghost);
443 :
444 : static void Copy (iMultiFab& dst,
445 : const iMultiFab& src,
446 : int srccomp,
447 : int dstcomp,
448 : int numcomp,
449 : const IntVect& nghost);
450 :
451 : /**
452 : * \brief Subtract src from dst including nghost ghost cells.
453 : * The two iMultiFabs MUST have the same underlying BoxArray.
454 : *
455 : * \param dst
456 : * \param src
457 : * \param srccomp
458 : * \param dstcomp
459 : * \param numcomp
460 : * \param nghost
461 : */
462 : static void Subtract (iMultiFab& dst,
463 : const iMultiFab& src,
464 : int srccomp,
465 : int dstcomp,
466 : int numcomp,
467 : int nghost);
468 :
469 : /**
470 : * \brief Multiply dst by src including nghost ghost cells.
471 : * The two iMultiFabs MUST have the same underlying BoxArray.
472 : *
473 : * \param dst
474 : * \param src
475 : * \param srccomp
476 : * \param dstcomp
477 : * \param numcomp
478 : * \param nghost
479 : */
480 : static void Multiply (iMultiFab& dst,
481 : const iMultiFab& src,
482 : int srccomp,
483 : int dstcomp,
484 : int numcomp,
485 : int nghost);
486 :
487 : /**
488 : * \brief Divide dst by src including nghost ghost cells.
489 : * The two iMultiFabs MUST have the same underlying BoxArray.
490 : *
491 : * \param dst
492 : * \param src
493 : * \param srccomp
494 : * \param dstcomp
495 : * \param numcomp
496 : * \param nghost
497 : */
498 : static void Divide (iMultiFab& dst,
499 : const iMultiFab& src,
500 : int srccomp,
501 : int dstcomp,
502 : int numcomp,
503 : int nghost);
504 :
505 : void define (const BoxArray& bxs,
506 : const DistributionMapping& dm,
507 : int nvar,
508 : const IntVect& ngrow,
509 : #ifdef AMREX_STRICT_MODE
510 : const MFInfo& info,
511 : const FabFactory<IArrayBox>& factory);
512 : #else
513 : const MFInfo& info = MFInfo(),
514 : const FabFactory<IArrayBox>& factory = DefaultFabFactory<IArrayBox>());
515 : #endif
516 :
517 : void define (const BoxArray& bxs,
518 : const DistributionMapping& dm,
519 : int nvar,
520 : int ngrow,
521 : #ifdef AMREX_STRICT_MODE
522 : const MFInfo& info,
523 : const FabFactory<IArrayBox>& factory);
524 : #else
525 : const MFInfo& info = MFInfo(),
526 : const FabFactory<IArrayBox>& factory = DefaultFabFactory<IArrayBox>());
527 : #endif
528 :
529 : static void Initialize ();
530 : static void Finalize ();
531 : };
532 :
533 : // ngrow != IntVect{0} is a special case that should not be used in normal cases,
534 : // because it may mark valid cells as non-owner and ghost cells as owners.
535 : std::unique_ptr<iMultiFab>
536 : OwnerMask (FabArrayBase const& mf, const Periodicity& period, const IntVect& ngrow=IntVect{0});
537 :
538 : }
539 :
540 : #endif /*BL_IMULTIFAB_H*/
|