Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_EB2_Level.H
Go to the documentation of this file.
1#ifndef AMREX_EB2_LEVEL_H_
2#define AMREX_EB2_LEVEL_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_ParmParse.H>
6#include <AMReX_Geometry.H>
7#include <AMReX_MultiFab.H>
8#include <AMReX_iMultiFab.H>
9#include <AMReX_LayoutData.H>
10#include <AMReX_VisMF.H>
11#include <AMReX_Array.H>
12#include <AMReX_EBCellFlag.H>
13#include <AMReX_MultiCutFab.H>
14#if (AMREX_SPACEDIM == 3)
15# include <AMReX_MarchingCubes.H>
16#endif
17#include <AMReX_EB2_MultiGFab.H>
18#include <AMReX_EB2_C.H>
20
21#ifdef AMREX_USE_OMP
22#include <omp.h>
23#endif
24
25#include <cmath>
26#include <limits>
27#include <map>
28#include <type_traits>
29#include <unordered_map>
30
31namespace amrex
32{
33 class STLtools;
34}
35
36namespace amrex::EB2 {
37
38class IndexSpace;
39template <typename G> class GShopLevel;
40
41class Level
42{
43public:
44
45 bool isAllRegular () const noexcept { return m_allregular; }
46 bool isOK () const noexcept { return m_ok; }
48 void fillVolFrac (MultiFab& vfrac, const Geometry& geom) const;
49 void fillCentroid (MultiCutFab& centroid, const Geometry& geom) const;
50 void fillCentroid ( MultiFab& centroid, const Geometry& geom) const;
51 void fillBndryArea (MultiCutFab& bndryarea, const Geometry& geom) const;
52 void fillBndryArea ( MultiFab& bndryarea, const Geometry& geom) const;
53 void fillBndryCent (MultiCutFab& bndrycent, const Geometry& geom) const;
54 void fillBndryCent ( MultiFab& bndrycent, const Geometry& geom) const;
55 void fillBndryNorm (MultiCutFab& bndrynorm, const Geometry& geom) const;
56 void fillBndryNorm ( MultiFab& bndrynorm, const Geometry& geom) const;
57 void fillAreaFrac (Array<MultiCutFab*,AMREX_SPACEDIM> const& areafrac, const Geometry& geom) const;
58 void fillAreaFrac (Array< MultiFab*,AMREX_SPACEDIM> const& areafrac, const Geometry& geom) const;
59 void fillFaceCent (Array<MultiCutFab*,AMREX_SPACEDIM> const& facecent, const Geometry& geom) const;
60 void fillFaceCent (Array< MultiFab*,AMREX_SPACEDIM> const& facecent, const Geometry& geom) const;
61 void fillEdgeCent (Array<MultiCutFab*,AMREX_SPACEDIM> const& edgecent, const Geometry& geom) const;
62 void fillEdgeCent (Array< MultiFab*,AMREX_SPACEDIM> const& edgecent, const Geometry& geom) const;
63 void fillLevelSet (MultiFab& levelset, const Geometry& geom) const;
64
65 const BoxArray& boxArray () const noexcept { return m_grids; }
66 const DistributionMapping& DistributionMap () const noexcept { return m_dmap; }
67
68 Level (IndexSpace const* is, const Geometry& geom) : m_geom(geom), m_parent(is) {}
69 void prepareForCoarsening (const Level& rhs, int max_grid_size, IntVect const& ngrow);
70
71 const Geometry& Geom () const noexcept { return m_geom; }
72 IndexSpace const* getEBIndexSpace () const noexcept { return m_parent; }
73
74 IntVect const& nGrowVect () const noexcept { return m_ngrow; }
75
76 void write_to_chkpt_file (const std::string& fname, bool extend_domain_face, int max_grid_size) const;
77
78 bool hasEBInfo () const noexcept { return m_has_eb_info; }
79 void fillCutCellMask (iMultiFab& cutcellmask, const Geometry& geom) const;
80
81 void setShift (int direction, int ncells);
82
83// public: // for cuda
84 int coarsenFromFine (Level& fineLevel, bool fill_boundary);
85 void buildCellFlag ();
86 void buildCutCellMask (Level const& fine_level);
87
88protected:
89
96#if (AMREX_SPACEDIM == 3)
97 std::map<int,std::unique_ptr<MC::MCFab>> m_marching_cubes;
98#endif
111 bool m_allregular = false;
112 bool m_ok = false;
113 bool m_has_eb_info = true;
116
117private:
118 template <typename G> friend class GShopLevel;
119
120 // Need this function to work around a gcc bug.
121 void setRegularLevel () {
122 m_allregular = true;
123 m_ok = true;
124 m_has_eb_info = false;
125 }
126};
127
128template <typename G>
130 : public Level
131{
132public:
133 GShopLevel (IndexSpace const* is, G const& gshop, const Geometry& geom, int max_grid_size,
134 int ngrow, bool extend_domain_face, int num_crse_opt);
135 GShopLevel (IndexSpace const* is, int ilev, int max_grid_size, int ngrow,
136 const Geometry& geom, GShopLevel<G>& fineLevel);
137 GShopLevel (IndexSpace const* is, const Geometry& geom);
138 void define_fine (G const& gshop, const Geometry& geom,
139 int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt);
140
141#if (AMREX_SPACEDIM == 3)
142 template <typename GS = G, std::enable_if_t<std::is_same_v<GS,STLtools>,int> FOO = 0>
143 void define_fine_mvmc (GS const& gshop, const Geometry& geom,
144 int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt);
145#endif
146
147 static GShopLevel<G>
148 makeAllRegular(IndexSpace const* is, const Geometry& geom)
149 {
150 GShopLevel<G> r(is, geom);
151 r.setRegularLevel();
152 return r;
153 }
154
155private:
156
157 Box m_bounding_box;
158
159 void prepare_grids (G const& gshop, Geometry const& geom, int max_grid_size, int ngrow,
160 bool extend_domain_face, int num_crse_opt);
161
162};
163
164template <typename G>
166 : Level(is, geom)
167{}
168
169template <typename G>
170GShopLevel<G>::GShopLevel (IndexSpace const* is, G const& gshop, const Geometry& geom,
171 int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
172 : Level(is, geom)
173{
174 if (std::is_same_v<typename G::FunctionType, AllRegularIF>) {
175 m_allregular = true;
176 m_ok = true;
177 return;
178 }
179
180 define_fine(gshop, geom, max_grid_size, ngrow, extend_domain_face, num_crse_opt);
181}
182
183template <typename G>
184void
185GShopLevel<G>::prepare_grids (G const& gshop, const Geometry& geom,
186 int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
187{
188 if (amrex::Verbose() > 0 && extend_domain_face == false) {
189 amrex::Print() << "AMReX WARNING: extend_domain_face=false is not recommended!\n";
190 }
191
192 BL_PROFILE("EB2::GShopLevel()-prepare_grids");
193
194 // make sure ngrow is multiple of 16
195 m_ngrow = IntVect{static_cast<int>(std::ceil(ngrow/16.)) * 16};
196
197 Box const& domain = geom.Domain();
198 Box domain_grown = domain;
199 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
200 if (geom.isPeriodic(idim)) {
201 m_ngrow[idim] = 0;
202 } else {
203 m_ngrow[idim] = std::min(m_ngrow[idim], domain_grown.length(idim));
204 }
205 }
206 domain_grown.grow(m_ngrow);
207 m_bounding_box = (extend_domain_face) ? domain : domain_grown;
208 m_bounding_box.surroundingNodes();
209
210 BoxList cut_boxes;
211 BoxList covered_boxes;
212
213 const int nprocs = ParallelDescriptor::NProcs();
214 const int iproc = ParallelDescriptor::MyProc();
215
216 num_crse_opt = std::max(0,std::min(8,num_crse_opt));
217 for (int clev = num_crse_opt; clev >= 0; --clev) {
218 IntVect crse_ratio(1 << clev);
219 if (domain.coarsenable(crse_ratio)) {
220 Box const& crse_bounding_box = amrex::coarsen(m_bounding_box, crse_ratio);
221 Geometry const& crse_geom = amrex::coarsen(geom, crse_ratio);
222 BoxList test_boxes;
223 if (cut_boxes.isEmpty()) {
224 covered_boxes.clear();
225 test_boxes = BoxList(crse_geom.Domain());
226 test_boxes.maxSize(max_grid_size);
227 } else {
228 test_boxes.swap(cut_boxes);
229 test_boxes.coarsen(crse_ratio);
230 test_boxes.maxSize(max_grid_size);
231 }
232
233 const Long nboxes = test_boxes.size();
234 const auto& boxes = test_boxes.data();
235 for (Long i = iproc; i < nboxes; i += nprocs) {
236 const Box& vbx = boxes[i];
237 const Box& gbx = amrex::surroundingNodes(amrex::grow(vbx,1));
238 auto box_type = gshop.getBoxType(gbx&crse_bounding_box,crse_geom,RunOn::Gpu);
239 if (box_type == gshop.allcovered) {
240 covered_boxes.push_back(amrex::refine(vbx, crse_ratio));
241 } else if (box_type == gshop.mixedcells) {
242 cut_boxes.push_back(amrex::refine(vbx, crse_ratio));
243 }
244 }
245
246 amrex::AllGatherBoxes(cut_boxes.data());
247 }
248 }
249
250 amrex::AllGatherBoxes(covered_boxes.data());
251
252 if (m_ngrow != 0) {
253 auto grow_at_domain_boundary = [&] (BoxList& bl)
254 {
255 const IntVect& domlo = domain.smallEnd();
256 const IntVect& domhi = domain.bigEnd();
257 for (auto& b : bl) {
258 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
259 if (m_ngrow[idim] != 0) {
260 if (b.smallEnd(idim) == domlo[idim]) {
261 b.growLo(idim,m_ngrow[idim]);
262 }
263 if (b.bigEnd(idim) == domhi[idim]) {
264 b.growHi(idim,m_ngrow[idim]);
265 }
266 }
267 }
268 }
269 };
270 grow_at_domain_boundary(covered_boxes);
271 grow_at_domain_boundary(cut_boxes);
272 }
273
274 if ( cut_boxes.isEmpty() &&
275 !covered_boxes.isEmpty())
276 {
277 amrex::Abort("AMReX_EB2_Level.H: Domain is completely covered");
278 }
279
280 if (!covered_boxes.isEmpty()) {
281 if (num_crse_opt > 2) { // don't want the box too big
282 covered_boxes.maxSize(max_grid_size*4);
283 }
284 m_covered_grids = BoxArray(std::move(covered_boxes));
285 }
286
287 if (cut_boxes.isEmpty()) {
288 m_grids = BoxArray();
289 m_dmap = DistributionMapping();
290 m_allregular = true;
291 m_ok = true;
292 } else {
293 m_grids = BoxArray(std::move(cut_boxes));
294 m_dmap = DistributionMapping(m_grids);
295 }
296}
297
298template <typename G>
299void
300GShopLevel<G>::define_fine (G const& gshop, const Geometry& geom,
301 int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
302{
303 BL_PROFILE("EB2::GShopLevel()-fine");
304
305#ifdef AMREX_USE_FLOAT
306 Real small_volfrac = 1.e-5_rt;
307#else
308 Real small_volfrac = 1.e-14;
309#endif
310 bool cover_multiple_cuts = false;
311 int maxiter = 32;
312 {
313 ParmParse pp("eb2");
314 pp.queryAdd("small_volfrac", small_volfrac);
315 pp.queryAdd("cover_multiple_cuts", cover_multiple_cuts);
316 pp.queryAdd("maxiter", maxiter);
317 }
318 maxiter = std::min(100000, maxiter);
319
320 prepare_grids(gshop, geom, max_grid_size, ngrow, extend_domain_face, num_crse_opt);
321
322 if (m_allregular) { return; }
323
324 m_mgf.define(m_grids, m_dmap);
325 const int ng = GFab::ng;
326 MFInfo mf_info;
327 mf_info.SetTag("EB2::Level");
328 m_cellflag.define(m_grids, m_dmap, 1, ng, mf_info);
329 m_volfrac.define(m_grids, m_dmap, 1, ng, mf_info);
330 m_centroid.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
331 m_bndryarea.define(m_grids, m_dmap, 1, ng, mf_info);
332 m_bndrycent.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
333 m_bndrynorm.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
334 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
335 m_areafrac[idim].define(amrex::convert(m_grids, IntVect::TheDimensionVector(idim)),
336 m_dmap, 1, ng, mf_info);
337 m_facecent[idim].define(amrex::convert(m_grids, IntVect::TheDimensionVector(idim)),
338 m_dmap, AMREX_SPACEDIM-1, ng, mf_info);
339 IntVect edge_type{1}; edge_type[idim] = 0;
340 m_edgecent[idim].define(amrex::convert(m_grids, edge_type), m_dmap, 1, ng, mf_info);
341 }
342
343 const auto dx = geom.CellSizeArray();
344 const auto problo = geom.ProbLoArray();
345
346 auto bounding_box = m_bounding_box;
347 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
348 if (!extend_domain_face || geom.isPeriodic(idim)) {
349 bounding_box.grow(idim,GFab::ng);
350 }
351 }
352
353 RunOn gshop_run_on = (Gpu::inLaunchRegion() && gshop.isGPUable())
355
356 bool hybrid = Gpu::inLaunchRegion() && (gshop_run_on == RunOn::Cpu);
357 amrex::ignore_unused(hybrid);
358
359 int iter = 0;
360 for (; iter < maxiter; ++iter)
361 {
362 int nsmallcells = 0;
363 int nmulticuts = 0;
364#ifdef AMREX_USE_OMP
365#pragma omp parallel if (Gpu::notInLaunchRegion()) reduction(+:nsmallcells,nmulticuts)
366#endif
367 {
368#if (AMREX_SPACEDIM == 3)
369 Array<BaseFab<Real>, AMREX_SPACEDIM> M2;
370 EBCellFlagFab cellflagtmp;
371#endif
372 for (MFIter mfi(m_mgf); mfi.isValid(); ++mfi)
373 {
374 auto& gfab = m_mgf[mfi];
375 const Box& vbx = gfab.validbox();
376
377 auto& levelset = gfab.getLevelSet();
378 if (iter == 0) {
379 gshop.fillFab(levelset, geom, gshop_run_on, bounding_box);
380#ifdef AMREX_USE_GPU
381 if (hybrid) {
382 levelset.prefetchToDevice();
383 }
384#endif
385 }
386
387 auto& cellflag = m_cellflag[mfi];
388
389 gfab.buildTypes(cellflag);
390
391 Array4<Real const> const& clst = levelset.const_array();
392 Array4<Real > const& lst = levelset.array();
393 Array4<EBCellFlag> const& cfg = m_cellflag.array(mfi);
394 Array4<Real> const& vfr = m_volfrac.array(mfi);
395 Array4<Real> const& ctr = m_centroid.array(mfi);
396 Array4<Real> const& bar = m_bndryarea.array(mfi);
397 Array4<Real> const& bct = m_bndrycent.array(mfi);
398 Array4<Real> const& bnm = m_bndrynorm.array(mfi);
399 AMREX_D_TERM(Array4<Real> const& apx = m_areafrac[0].array(mfi);,
400 Array4<Real> const& apy = m_areafrac[1].array(mfi);,
401 Array4<Real> const& apz = m_areafrac[2].array(mfi););
402 AMREX_D_TERM(Array4<Real> const& fcx = m_facecent[0].array(mfi);,
403 Array4<Real> const& fcy = m_facecent[1].array(mfi);,
404 Array4<Real> const& fcz = m_facecent[2].array(mfi););
405
406 auto& facetype = gfab.getFaceType();
407 AMREX_D_TERM(Array4<Type_t> const& ftx = facetype[0].array();,
408 Array4<Type_t> const& fty = facetype[1].array();,
409 Array4<Type_t> const& ftz = facetype[2].array(););
410
411 int nmc = 0;
412 int nsm = 0;
413
414#if (AMREX_SPACEDIM == 3)
415 auto& edgetype = gfab.getEdgeType();
416 Array4<Type_t const> const& xdg = edgetype[0].const_array();
417 Array4<Type_t const> const& ydg = edgetype[1].const_array();
418 Array4<Type_t const> const& zdg = edgetype[2].const_array();
419
420 Array4<Real> const& xip = m_edgecent[0].array(mfi);
421 Array4<Real> const& yip = m_edgecent[1].array(mfi);
422 Array4<Real> const& zip = m_edgecent[2].array(mfi);
423
424 if (iter == 0) {
425#ifdef AMREX_USE_GPU
426 if (hybrid) {
428 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
429 edgetype[idim].prefetchToHost();
430 m_edgecent[idim][mfi].prefetchToHost();
431 }
432 }
433#endif
434
435 gshop.getIntercept({xip,yip,zip}, {xdg,ydg,zdg}, clst,
436 geom, gshop_run_on, bounding_box);
437
438#ifdef AMREX_USE_GPU
439 if (hybrid) {
440 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
441 edgetype[idim].prefetchToDevice();
442 m_edgecent[idim][mfi].prefetchToDevice();
443 }
444 }
445#endif
446 } else {
447 gshop.updateIntercept({xip,yip,zip}, {xdg,ydg,zdg}, clst, geom);
448 }
449
450 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
451 const Box& b = facetype[idim].box();
452 M2[idim].resize(b,3);
453 }
454 Array4<Real> const& xm2 = M2[0].array();
455 Array4<Real> const& ym2 = M2[1].array();
456 Array4<Real> const& zm2 = M2[2].array();
457
458 nmc = build_faces(vbx, cfg, ftx, fty, ftz, xdg, ydg, zdg, lst,
459 xip, yip, zip, apx, apy, apz, fcx, fcy, fcz,
460 xm2, ym2, zm2, dx, problo, cover_multiple_cuts);
461
462 cellflagtmp.resize(m_cellflag[mfi].box());
463 Elixir cellflagtmp_eli = cellflagtmp.elixir();
464 Array4<EBCellFlag> const& cfgtmp = cellflagtmp.array();
465
466 build_cells(vbx, cfg, ftx, fty, ftz, apx, apy, apz,
467 fcx, fcy, fcz, xm2, ym2, zm2, dx, vfr, ctr,
468 bar, bct, bnm, cfgtmp, lst,
469 small_volfrac, geom, extend_domain_face, cover_multiple_cuts,
470 nsm, nmc);
471
472 // Because it is used in a synchronous reduction kernel in
473 // build_cells, we do not need to worry about M2's lifetime.
474 // But we still need to use Elixir to extend the life of
475 // cellflagtmp.
476
477#elif (AMREX_SPACEDIM == 2)
478 Array4<Real> const& xip = m_edgecent[0].array(mfi);
479 Array4<Real> const& yip = m_edgecent[1].array(mfi);
480
481 if (iter == 0) {
482#ifdef AMREX_USE_GPU
483 if (hybrid) {
485 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
486 facetype[idim].prefetchToHost();
487 m_edgecent[idim][mfi].prefetchToHost();
488 }
489 }
490#endif
491
492 // yes, factype[1] and then [0]
493 gshop.getIntercept({xip,yip},
494 {facetype[1].const_array(), facetype[0].const_array()},
495 clst, geom, gshop_run_on, bounding_box);
496
497#ifdef AMREX_USE_GPU
498 if (hybrid) {
499 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
500 facetype[idim].prefetchToDevice();
501 m_edgecent[idim][mfi].prefetchToDevice();
502 }
503 }
504#endif
505 } else {
506 gshop.updateIntercept({xip,yip},
507 {facetype[1].const_array(), facetype[0].const_array()},
508 clst, geom);
509 }
510
511 nmc = build_faces(vbx, cfg, ftx, fty, lst, xip, yip, apx, apy, fcx, fcy,
512 dx, problo, cover_multiple_cuts, nsm);
513
514 build_cells(vbx, cfg, ftx, fty, apx, apy, dx, vfr, ctr,
515 bar, bct, bnm, lst, small_volfrac, geom, extend_domain_face,
516 nsm, nmc);
517#endif
518 nsmallcells += nsm;
519 nmulticuts += nmc;
520 }
521 }
522
523 ParallelAllReduce::Sum<int>({nsmallcells,nmulticuts}, ParallelContext::CommunicatorSub());
524 if (nsmallcells == 0 && nmulticuts == 0) {
525 break;
526 } else {
527 auto ls = m_mgf.getLevelSet();
528 // This is an alias MultiFab, therefore FillBoundary on it is fine.
529 ls.FillBoundary(geom.periodicity());
530 if (amrex::Verbose() > 0) {
531 if (nsmallcells) {
532 amrex::Print() << "AMReX EB: Iter. " << iter+1 << " fixed " << nsmallcells
533 << " small cells" << '\n';
534 }
535 if (nmulticuts) {
536 amrex::Print() << "AMReX EB: Iter. " << iter+1 << " fixed " << nmulticuts
537 << " multicuts" << '\n';
538 }
539 }
540 }
541 }
542
543 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < maxiter, "EB: failed to fix small cells");
544
545#ifdef AMREX_USE_OMP
546#pragma omp parallel if (Gpu::notInLaunchRegion())
547#endif
548 for (MFIter mfi(m_mgf); mfi.isValid(); ++mfi)
549 {
550 auto& gfab = m_mgf[mfi];
551 auto const& levelset = gfab.getLevelSet();
552 Array4<Real const> const& clst = levelset.const_array();
553
554 AMREX_D_TERM(Array4<Real> const& xip = m_edgecent[0].array(mfi);,
555 Array4<Real> const& yip = m_edgecent[1].array(mfi);,
556 Array4<Real> const& zip = m_edgecent[2].array(mfi);)
557#if (AMREX_SPACEDIM == 3)
558 auto const& edgetype = gfab.getEdgeType();
559 Array4<Type_t const> const& xdg = edgetype[0].const_array();
560 Array4<Type_t const> const& ydg = edgetype[1].const_array();
561 Array4<Type_t const> const& zdg = edgetype[2].const_array();
562
563 intercept_to_edge_centroid(xip, yip, zip, xdg, ydg, zdg, clst, dx, problo);
564
565#elif (AMREX_SPACEDIM == 2)
566 auto& facetype = gfab.getFaceType();
567 Array4<Type_t const> const& ftx = facetype[0].const_array();
568 Array4<Type_t const> const& fty = facetype[1].const_array();
569 // fty then ftx
570 intercept_to_edge_centroid(xip, yip, fty, ftx, clst, dx, problo);
571#endif
572 }
573
574 m_levelset = m_mgf.getLevelSet();
575
576 m_ok = true;
577}
578
579#if (AMREX_SPACEDIM == 3)
580template <typename G>
581template <typename GS, std::enable_if_t<std::is_same_v<GS,STLtools>,int> FOO>
582void
583GShopLevel<G>::define_fine_mvmc (GS const& gshop, const Geometry& geom,
584 int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
585{
586 BL_PROFILE("EB2::GShopLevel()-fine-mvmc");
587
588 prepare_grids(gshop, geom, max_grid_size, ngrow, extend_domain_face, num_crse_opt);
589
590 if (m_allregular) { return; }
591
593
594 m_sdf.define(amrex::convert(m_grids,IntVect(1)), m_dmap, 1, IntVect(1));
595 gshop.fillSignedDistance(m_sdf, m_sdf.nGrowVect(), geom);
596
597 for (MFIter mfi(m_sdf,MFItInfo().DisableDeviceSync()); mfi.isValid(); ++mfi) {
598 m_marching_cubes[mfi.index()] = std::make_unique<MC::MCFab>();
599 }
600
601 MFItInfo info{};
602
603#if defined(AMREX_USE_OMP) && !defined(AMREX_USE_GPU)
604 info.SetDynamic(true);
605#pragma omp parallel
606#endif
607 for (MFIter mfi(m_sdf,info); mfi.isValid(); ++mfi) {
608 amrex::MC::marching_cubes(geom, m_sdf[mfi], *m_marching_cubes[mfi.index()]);
609 }
610
611 MC::write_stl("test.stl", m_marching_cubes);
612
613 m_ok = false; // xxxxx TODO
614}
615#endif
616
617template <typename G>
618GShopLevel<G>::GShopLevel (IndexSpace const* is, int /*ilev*/, int max_grid_size, int /*ngrow*/,
619 const Geometry& geom, GShopLevel<G>& fineLevel)
620 : Level(is, geom)
621{
622 if (fineLevel.isAllRegular()) {
623 m_allregular = true;
624 m_ok = true;
625 return;
626 }
627
628 BL_PROFILE("EB2::GShopLevel()-coarse");
629
630 const BoxArray& fine_grids = fineLevel.m_grids;
631 const BoxArray& fine_covered_grids = fineLevel.m_covered_grids;
632
633 const int coarse_ratio = 2;
634 const int min_width = 8;
635 bool coarsenable = fine_grids.coarsenable(coarse_ratio, min_width)
636 && (fine_covered_grids.empty() || fine_covered_grids.coarsenable(coarse_ratio));
637
638 m_ngrow = amrex::coarsen(fineLevel.m_ngrow,2);
639 if (amrex::scale(m_ngrow,2) != fineLevel.m_ngrow) {
641 }
642
643 if (coarsenable)
644 {
645 int ierr = coarsenFromFine(fineLevel, true);
646 m_ok = (ierr == 0);
647 }
648 else
649 {
650 Level fine_level_2(is, fineLevel.m_geom);
651 fine_level_2.prepareForCoarsening(fineLevel, max_grid_size, amrex::scale(m_ngrow,2));
652 int ierr = coarsenFromFine(fine_level_2, false);
653 m_ok = (ierr == 0);
654 }
655}
656
657}
658
659#endif
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:49
amrex::ParmParse pp
Input file parser instance for the given namespace.
Definition AMReX_HypreIJIface.cpp:15
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
void resize(const Box &b, int N=1, Arena *ar=nullptr)
This function resizes a BaseFab so it covers the Box b with N components.
Definition AMReX_BaseFab.H:1629
Array4< T const > array() const noexcept
Definition AMReX_BaseFab.H:382
Elixir elixir() noexcept
Definition AMReX_BaseFab.H:1671
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
bool coarsenable(int refinement_ratio, int min_width=1) const
Coarsen each Box in the BoxArray to the specified ratio.
Definition AMReX_BoxArray.cpp:615
bool empty() const noexcept
Return whether the BoxArray is empty.
Definition AMReX_BoxArray.H:621
__host__ __device__ BoxND< new_dim > resize() const noexcept
Return a new BoxND of size new_dim by either shrinking or expanding this BoxND.
Definition AMReX_Box.H:893
GpuArray< Real, 3 > CellSizeArray() const noexcept
Definition AMReX_CoordSys.H:76
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
static constexpr int ng
Definition AMReX_EB2_MultiGFab.H:19
Definition AMReX_EB2_Level.H:131
GShopLevel(IndexSpace const *is, const Geometry &geom)
Definition AMReX_EB2_Level.H:165
void define_fine(G const &gshop, const Geometry &geom, int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
Definition AMReX_EB2_Level.H:300
GShopLevel(IndexSpace const *is, int ilev, int max_grid_size, int ngrow, const Geometry &geom, GShopLevel< G > &fineLevel)
Definition AMReX_EB2_Level.H:618
void define_fine_mvmc(GS const &gshop, const Geometry &geom, int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
Definition AMReX_EB2_Level.H:583
GShopLevel(IndexSpace const *is, G const &gshop, const Geometry &geom, int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
Definition AMReX_EB2_Level.H:170
static GShopLevel< G > makeAllRegular(IndexSpace const *is, const Geometry &geom)
Definition AMReX_EB2_Level.H:148
Definition AMReX_EB2.H:28
Definition AMReX_EB2_Level.H:42
MultiFab m_bndryarea
Definition AMReX_EB2_Level.H:104
bool isOK() const noexcept
Definition AMReX_EB2_Level.H:46
void fillFaceCent(Array< MultiCutFab *, 3 > const &facecent, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:794
MultiFab m_volfrac
Definition AMReX_EB2_Level.H:102
MultiGFab m_mgf
Definition AMReX_EB2_Level.H:99
bool m_allregular
Definition AMReX_EB2_Level.H:111
IntVect m_shift
Definition AMReX_EB2_Level.H:114
iMultiFab m_cutcellmask
Definition AMReX_EB2_Level.H:110
Array< MultiFab, 3 > m_facecent
Definition AMReX_EB2_Level.H:108
Array< MultiFab, 3 > m_edgecent
Definition AMReX_EB2_Level.H:109
void prepareForCoarsening(const Level &rhs, int max_grid_size, IntVect const &ngrow)
Definition AMReX_EB2_Level.cpp:10
void fillCutCellMask(iMultiFab &cutcellmask, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:946
MultiFab m_centroid
Definition AMReX_EB2_Level.H:103
FabArray< EBCellFlagFab > m_cellflag
Definition AMReX_EB2_Level.H:101
const DistributionMapping & DistributionMap() const noexcept
Definition AMReX_EB2_Level.H:66
BoxArray m_covered_grids
Definition AMReX_EB2_Level.H:93
void fillLevelSet(MultiFab &levelset, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:907
Geometry m_geom
Definition AMReX_EB2_Level.H:90
Array< MultiFab, 3 > m_areafrac
Definition AMReX_EB2_Level.H:107
IntVect m_ngrow
Definition AMReX_EB2_Level.H:91
IndexSpace const * m_parent
Definition AMReX_EB2_Level.H:115
BoxArray m_grids
Definition AMReX_EB2_Level.H:92
const Geometry & Geom() const noexcept
Definition AMReX_EB2_Level.H:71
std::map< int, std::unique_ptr< MC::MCFab > > m_marching_cubes
Definition AMReX_EB2_Level.H:97
bool hasEBInfo() const noexcept
Definition AMReX_EB2_Level.H:78
void write_to_chkpt_file(const std::string &fname, bool extend_domain_face, int max_grid_size) const
Definition AMReX_EB2_Level.cpp:956
void fillEBCellFlag(FabArray< EBCellFlagFab > &cellflag, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:429
void fillEdgeCent(Array< MultiCutFab *, 3 > const &edgecent, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:835
int coarsenFromFine(Level &fineLevel, bool fill_boundary)
Definition AMReX_EB2_Level.cpp:84
MultiFab m_bndrynorm
Definition AMReX_EB2_Level.H:106
bool m_has_eb_info
Definition AMReX_EB2_Level.H:113
void fillAreaFrac(Array< MultiCutFab *, 3 > const &areafrac, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:640
MultiFab m_bndrycent
Definition AMReX_EB2_Level.H:105
MultiFab m_levelset
Definition AMReX_EB2_Level.H:100
void setShift(int direction, int ncells)
Definition AMReX_EB2_Level.cpp:1047
bool isAllRegular() const noexcept
Definition AMReX_EB2_Level.H:45
void fillBndryCent(MultiCutFab &bndrycent, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:592
void fillBndryNorm(MultiCutFab &bndrynorm, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:616
MultiFab m_sdf
Definition AMReX_EB2_Level.H:95
IntVect const & nGrowVect() const noexcept
Definition AMReX_EB2_Level.H:74
void fillVolFrac(MultiFab &vfrac, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:483
const BoxArray & boxArray() const noexcept
Definition AMReX_EB2_Level.H:65
Level(IndexSpace const *is, const Geometry &geom)
Definition AMReX_EB2_Level.H:68
friend class GShopLevel
Definition AMReX_EB2_Level.H:118
bool m_ok
Definition AMReX_EB2_Level.H:112
void fillCentroid(MultiCutFab &centroid, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:544
DistributionMapping m_dmap
Definition AMReX_EB2_Level.H:94
void fillBndryArea(MultiCutFab &bndryarea, const Geometry &geom) const
Definition AMReX_EB2_Level.cpp:568
void buildCellFlag()
Definition AMReX_EB2_Level.cpp:403
void buildCutCellMask(Level const &fine_level)
Definition AMReX_EB2_Level.cpp:966
IndexSpace const * getEBIndexSpace() const noexcept
Definition AMReX_EB2_Level.H:72
Definition AMReX_EB2_MultiGFab.H:65
Definition AMReX_EBCellFlag.H:287
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:349
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:211
Periodicity periodicity() const noexcept
Definition AMReX_Geometry.H:356
GpuArray< Real, 3 > ProbLoArray() const noexcept
Definition AMReX_Geometry.H:187
bool isPeriodic(int dir) const noexcept
Is the domain periodic in the specified direction?
Definition AMReX_Geometry.H:332
Definition AMReX_GpuElixir.H:13
__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:679
__host__ static __device__ constexpr IntVectND< dim > TheDimensionVector(int d) noexcept
This static member function returns a reference to a constant IntVectND object, all of whose dim argu...
Definition AMReX_IntVect.H:698
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:85
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:169
Definition AMReX_MultiCutFab.H:81
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
Parse Parameters From Command Line and Input Files.
Definition AMReX_ParmParse.H:348
int queryAdd(std::string_view name, T &ref)
If name is found, the value in the ParmParse database will be stored in the ref argument....
Definition AMReX_ParmParse.H:1044
This class provides the user with a few print options.
Definition AMReX_Print.H:35
A Collection of IArrayBoxes.
Definition AMReX_iMultiFab.H:34
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
amrex_long Long
Definition AMReX_INT.H:30
__host__ __device__ BoxND< dim > coarsen(const BoxND< dim > &b, int ref_ratio) noexcept
Coarsen BoxND by given (positive) coarsening ratio.
Definition AMReX_Box.H:1409
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1280
__host__ __device__ BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition AMReX_Box.H:1459
std::array< T, N > Array
Definition AMReX_Array.H:26
int MyProc() noexcept
Definition AMReX_ParallelDescriptor.H:128
int NProcs() noexcept
Definition AMReX_ParallelDescriptor.H:255
Definition AMReX_FabArrayBase.H:33
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:263
bool inLaunchRegion() noexcept
Definition AMReX_GpuControl.H:92
void write_stl(std::string const &filename, std::map< int, std::unique_ptr< MCFab > > const &mc_fabs)
Definition AMReX_MarchingCubes.cpp:924
void marching_cubes(Geometry const &geom, FArrayBox &sdf_fab, MCFab &mc_fab)
Definition AMReX_MarchingCubes.cpp:693
void Initialize()
Definition AMReX_MarchingCubes.cpp:651
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
Definition AMReX_Amr.cpp:49
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1558
__host__ __device__ BoxND< dim > surroundingNodes(const BoxND< dim > &b, int dir) noexcept
Return a BoxND with NODE based coordinates in direction dir that encloses BoxND b....
Definition AMReX_Box.H:1522
RunOn
Definition AMReX_GpuControl.H:69
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
__host__ __device__ constexpr IntVectND< dim > scale(const IntVectND< dim > &p, int s) noexcept
Returns a IntVectND obtained by multiplying each of the components of this IntVectND by s.
Definition AMReX_IntVect.H:1013
int Verbose() noexcept
Definition AMReX.cpp:169
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:230
void AllGatherBoxes(Vector< Box > &bxs, int n_extra_reserve)
Definition AMReX_Box.cpp:124
A multidimensional array accessor.
Definition AMReX_Array4.H:283
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
MFInfo & SetTag() noexcept
Definition AMReX_FabArray.H:79
Definition AMReX_MFIter.H:20
MFItInfo & SetDynamic(bool f) noexcept
Definition AMReX_MFIter.H:40