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
46 bool isAllRegular () const noexcept { return m_allregular; }
48 bool isOK () const noexcept { return m_ok; }
59 void fillVolFrac (MultiFab& vfrac, const Geometry& geom) const;
63 void fillCentroid (MultiCutFab& centroid, const Geometry& geom) const;
65 void fillCentroid ( MultiFab& centroid, const Geometry& geom) const;
69 void fillBndryArea (MultiCutFab& bndryarea, const Geometry& geom) const;
71 void fillBndryArea ( MultiFab& bndryarea, const Geometry& geom) const;
75 void fillBndryCent (MultiCutFab& bndrycent, const Geometry& geom) const;
77 void fillBndryCent ( MultiFab& bndrycent, const Geometry& geom) const;
81 void fillBndryNorm (MultiCutFab& bndrynorm, const Geometry& geom) const;
83 void fillBndryNorm ( MultiFab& bndrynorm, const Geometry& geom) const;
87 void fillAreaFrac (Array<MultiCutFab*,AMREX_SPACEDIM> const& areafrac, const Geometry& geom) const;
89 void fillAreaFrac (Array< MultiFab*,AMREX_SPACEDIM> const& areafrac, const Geometry& geom) const;
93 void fillFaceCent (Array<MultiCutFab*,AMREX_SPACEDIM> const& facecent, const Geometry& geom) const;
95 void fillFaceCent (Array< MultiFab*,AMREX_SPACEDIM> const& facecent, const Geometry& geom) const;
99 void fillEdgeCent (Array<MultiCutFab*,AMREX_SPACEDIM> const& edgecent, const Geometry& geom) const;
101 void fillEdgeCent (Array< MultiFab*,AMREX_SPACEDIM> const& edgecent, const Geometry& geom) const;
105 void fillLevelSet (MultiFab& levelset, const Geometry& geom) const;
106
108 const BoxArray& boxArray () const noexcept { return m_grids; }
110 const DistributionMapping& DistributionMap () const noexcept { return m_dmap; }
111
113 Level (IndexSpace const* is, const Geometry& geom) : m_geom(geom), m_parent(is) {}
117 void prepareForCoarsening (const Level& rhs, int max_grid_size, IntVect const& ngrow);
118
120 const Geometry& Geom () const noexcept { return m_geom; }
122 IndexSpace const* getEBIndexSpace () const noexcept { return m_parent; }
123
125 IntVect const& nGrowVect () const noexcept { return m_ngrow; }
126
134 void write_to_chkpt_file (const std::string& fname, bool extend_domain_face, int max_grid_size) const;
135
137 bool hasEBInfo () const noexcept { return m_has_eb_info; }
139 void fillCutCellMask (iMultiFab& cutcellmask, const Geometry& geom) const;
140
144 void setShift (int direction, int ncells);
145
146// public: // for cuda
154 int coarsenFromFine (Level& fineLevel, bool fill_boundary);
156 void buildCellFlag ();
158 void buildCutCellMask (Level const& fine_level);
159
160protected:
161
168#if (AMREX_SPACEDIM == 3)
169 std::map<int,std::unique_ptr<MC::MCFab>> m_marching_cubes;
170#endif
183 bool m_allregular = false;
184 bool m_ok = false;
185 bool m_has_eb_info = true;
188
189private:
190 template <typename G> friend class GShopLevel;
191
192 // Need this function to work around a gcc bug.
193 void setRegularLevel () {
194 m_allregular = true;
195 m_ok = true;
196 m_has_eb_info = false;
197 }
198};
199
200template <typename G>
202 : public Level
203{
204public:
208 GShopLevel (IndexSpace const* is, G const& gshop, const Geometry& geom, int max_grid_size,
209 int ngrow, bool extend_domain_face, int num_crse_opt);
213 GShopLevel (IndexSpace const* is, int ilev, int max_grid_size, int ngrow,
214 const Geometry& geom, GShopLevel<G>& fineLevel);
216 GShopLevel (IndexSpace const* is, const Geometry& geom);
218 void define_fine (G const& gshop, const Geometry& geom,
219 int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt);
220
221#if (AMREX_SPACEDIM == 3)
225 template <typename GS = G>
226 requires (std::same_as<GS,STLtools>)
227 void define_fine_mvmc (GS const& gshop, const Geometry& geom,
228 int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt);
229#endif
230
234 static GShopLevel<G>
235 makeAllRegular(IndexSpace const* is, const Geometry& geom)
236 {
237 GShopLevel<G> r(is, geom);
238 r.setRegularLevel();
239 return r;
240 }
241
242private:
243
244 Box m_bounding_box;
245
246 void prepare_grids (G const& gshop, Geometry const& geom, int max_grid_size, int ngrow,
247 bool extend_domain_face, int num_crse_opt);
248
249};
250
251template <typename G>
253 : Level(is, geom)
254{}
255
256template <typename G>
257GShopLevel<G>::GShopLevel (IndexSpace const* is, G const& gshop, const Geometry& geom,
258 int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
259 : Level(is, geom)
260{
261 if (std::is_same_v<typename G::FunctionType, AllRegularIF>) {
262 m_allregular = true;
263 m_ok = true;
264 return;
265 }
266
267 define_fine(gshop, geom, max_grid_size, ngrow, extend_domain_face, num_crse_opt);
268}
269
270template <typename G>
271void
272GShopLevel<G>::prepare_grids (G const& gshop, const Geometry& geom,
273 int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
274{
275 if (amrex::Verbose() > 0 && extend_domain_face == false) {
276 amrex::Print() << "AMReX WARNING: extend_domain_face=false is not recommended!\n";
277 }
278
279 BL_PROFILE("EB2::GShopLevel()-prepare_grids");
280
281 // make sure ngrow is multiple of 16
282 m_ngrow = IntVect{static_cast<int>(std::ceil(ngrow/16.)) * 16};
283
284 Box const& domain = geom.Domain();
285 Box domain_grown = domain;
286 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
287 if (geom.isPeriodic(idim)) {
288 m_ngrow[idim] = 0;
289 } else {
290 m_ngrow[idim] = std::min(m_ngrow[idim], domain_grown.length(idim));
291 }
292 }
293 domain_grown.grow(m_ngrow);
294 m_bounding_box = (extend_domain_face) ? domain : domain_grown;
295 m_bounding_box.surroundingNodes();
296
297 BoxList cut_boxes;
298 BoxList covered_boxes;
299
300 const int nprocs = ParallelDescriptor::NProcs();
301 const int iproc = ParallelDescriptor::MyProc();
302
303 num_crse_opt = std::max(0,std::min(8,num_crse_opt));
304 for (int clev = num_crse_opt; clev >= 0; --clev) {
305 IntVect crse_ratio(1 << clev);
306 if (domain.coarsenable(crse_ratio)) {
307 Box const& crse_bounding_box = amrex::coarsen(m_bounding_box, crse_ratio);
308 Geometry const& crse_geom = amrex::coarsen(geom, crse_ratio);
309 BoxList test_boxes;
310 if (cut_boxes.isEmpty()) {
311 covered_boxes.clear();
312 test_boxes = BoxList(crse_geom.Domain());
313 test_boxes.maxSize(max_grid_size);
314 } else {
315 test_boxes.swap(cut_boxes);
316 test_boxes.coarsen(crse_ratio);
317 test_boxes.maxSize(max_grid_size);
318 }
319
320 const Long nboxes = test_boxes.size();
321 const auto& boxes = test_boxes.data();
322 for (Long i = iproc; i < nboxes; i += nprocs) {
323 const Box& vbx = boxes[i];
324 const Box& gbx = amrex::surroundingNodes(amrex::grow(vbx,1));
325 auto box_type = gshop.getBoxType(gbx&crse_bounding_box,crse_geom,RunOn::Gpu);
326 if (box_type == gshop.allcovered) {
327 covered_boxes.push_back(amrex::refine(vbx, crse_ratio));
328 } else if (box_type == gshop.mixedcells) {
329 cut_boxes.push_back(amrex::refine(vbx, crse_ratio));
330 }
331 }
332
333 amrex::AllGatherBoxes(cut_boxes.data());
334 }
335 }
336
337 amrex::AllGatherBoxes(covered_boxes.data());
338
339 if (m_ngrow != 0) {
340 auto grow_at_domain_boundary = [&] (BoxList& bl)
341 {
342 const IntVect& domlo = domain.smallEnd();
343 const IntVect& domhi = domain.bigEnd();
344 for (auto& b : bl) {
345 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
346 if (m_ngrow[idim] != 0) {
347 if (b.smallEnd(idim) == domlo[idim]) {
348 b.growLo(idim,m_ngrow[idim]);
349 }
350 if (b.bigEnd(idim) == domhi[idim]) {
351 b.growHi(idim,m_ngrow[idim]);
352 }
353 }
354 }
355 }
356 };
357 grow_at_domain_boundary(covered_boxes);
358 grow_at_domain_boundary(cut_boxes);
359 }
360
361 if ( cut_boxes.isEmpty() &&
362 !covered_boxes.isEmpty())
363 {
364 amrex::Abort("AMReX_EB2_Level.H: Domain is completely covered");
365 }
366
367 if (!covered_boxes.isEmpty()) {
368 if (num_crse_opt > 2) { // don't want the box too big
369 covered_boxes.maxSize(max_grid_size*4);
370 }
371 m_covered_grids = BoxArray(std::move(covered_boxes));
372 }
373
374 if (cut_boxes.isEmpty()) {
375 m_grids = BoxArray();
376 m_dmap = DistributionMapping();
377 m_allregular = true;
378 m_ok = true;
379 } else {
380 m_grids = BoxArray(std::move(cut_boxes));
381 m_dmap = DistributionMapping(m_grids);
382 }
383}
384
385template <typename G>
386void
387GShopLevel<G>::define_fine (G const& gshop, const Geometry& geom,
388 int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
389{
390 BL_PROFILE("EB2::GShopLevel()-fine");
391
392#ifdef AMREX_USE_FLOAT
393 Real small_volfrac = 1.e-5_rt;
394#else
395 Real small_volfrac = 1.e-14;
396#endif
397 bool cover_multiple_cuts = false;
398 int maxiter = 32;
399 {
400 ParmParse pp("eb2");
401 pp.queryAdd("small_volfrac", small_volfrac);
402 pp.queryAdd("cover_multiple_cuts", cover_multiple_cuts);
403 pp.queryAdd("maxiter", maxiter);
404 }
405 maxiter = std::min(100000, maxiter);
406
407 prepare_grids(gshop, geom, max_grid_size, ngrow, extend_domain_face, num_crse_opt);
408
409 if (m_allregular) { return; }
410
411 m_mgf.define(m_grids, m_dmap);
412 const int ng = GFab::ng;
413 MFInfo mf_info;
414 mf_info.SetTag("EB2::Level");
415 m_cellflag.define(m_grids, m_dmap, 1, ng, mf_info);
416 m_volfrac.define(m_grids, m_dmap, 1, ng, mf_info);
417 m_centroid.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
418 m_bndryarea.define(m_grids, m_dmap, 1, ng, mf_info);
419 m_bndrycent.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
420 m_bndrynorm.define(m_grids, m_dmap, AMREX_SPACEDIM, ng, mf_info);
421 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
422 m_areafrac[idim].define(amrex::convert(m_grids, IntVect::TheDimensionVector(idim)),
423 m_dmap, 1, ng, mf_info);
424 m_facecent[idim].define(amrex::convert(m_grids, IntVect::TheDimensionVector(idim)),
425 m_dmap, AMREX_SPACEDIM-1, ng, mf_info);
426 IntVect edge_type{1}; edge_type[idim] = 0;
427 m_edgecent[idim].define(amrex::convert(m_grids, edge_type), m_dmap, 1, ng, mf_info);
428 }
429
430 const auto dx = geom.CellSizeArray();
431 const auto problo = geom.ProbLoArray();
432
433 auto bounding_box = m_bounding_box;
434 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
435 if (!extend_domain_face || geom.isPeriodic(idim)) {
436 bounding_box.grow(idim,GFab::ng);
437 }
438 }
439
440 RunOn gshop_run_on = (Gpu::inLaunchRegion() && gshop.isGPUable())
442
443 bool hybrid = Gpu::inLaunchRegion() && (gshop_run_on == RunOn::Cpu);
444 amrex::ignore_unused(hybrid);
445
446 int iter = 0;
447 for (; iter < maxiter; ++iter)
448 {
449 int nsmallcells = 0;
450 int nmulticuts = 0;
451#ifdef AMREX_USE_OMP
452#pragma omp parallel if (Gpu::notInLaunchRegion()) reduction(+:nsmallcells,nmulticuts)
453#endif
454 {
455#if (AMREX_SPACEDIM == 3)
456 Array<BaseFab<Real>, AMREX_SPACEDIM> M2;
457 EBCellFlagFab cellflagtmp;
458#endif
459 for (MFIter mfi(m_mgf); mfi.isValid(); ++mfi)
460 {
461 auto& gfab = m_mgf[mfi];
462 const Box& vbx = gfab.validbox();
463
464 auto& levelset = gfab.getLevelSet();
465 if (iter == 0) {
466 gshop.fillFab(levelset, geom, gshop_run_on, bounding_box);
467#ifdef AMREX_USE_GPU
468 if (hybrid) {
469 levelset.prefetchToDevice();
470 }
471#endif
472 }
473
474 auto& cellflag = m_cellflag[mfi];
475
476 gfab.buildTypes(cellflag);
477
478 Array4<Real const> const& clst = levelset.const_array();
479 Array4<Real > const& lst = levelset.array();
480 Array4<EBCellFlag> const& cfg = m_cellflag.array(mfi);
481 Array4<Real> const& vfr = m_volfrac.array(mfi);
482 Array4<Real> const& ctr = m_centroid.array(mfi);
483 Array4<Real> const& bar = m_bndryarea.array(mfi);
484 Array4<Real> const& bct = m_bndrycent.array(mfi);
485 Array4<Real> const& bnm = m_bndrynorm.array(mfi);
486 AMREX_D_TERM(Array4<Real> const& apx = m_areafrac[0].array(mfi);,
487 Array4<Real> const& apy = m_areafrac[1].array(mfi);,
488 Array4<Real> const& apz = m_areafrac[2].array(mfi););
489 AMREX_D_TERM(Array4<Real> const& fcx = m_facecent[0].array(mfi);,
490 Array4<Real> const& fcy = m_facecent[1].array(mfi);,
491 Array4<Real> const& fcz = m_facecent[2].array(mfi););
492
493 auto& facetype = gfab.getFaceType();
494 AMREX_D_TERM(Array4<Type_t> const& ftx = facetype[0].array();,
495 Array4<Type_t> const& fty = facetype[1].array();,
496 Array4<Type_t> const& ftz = facetype[2].array(););
497
498 int nmc = 0;
499 int nsm = 0;
500
501#if (AMREX_SPACEDIM == 3)
502 auto& edgetype = gfab.getEdgeType();
503 Array4<Type_t const> const& xdg = edgetype[0].const_array();
504 Array4<Type_t const> const& ydg = edgetype[1].const_array();
505 Array4<Type_t const> const& zdg = edgetype[2].const_array();
506
507 Array4<Real> const& xip = m_edgecent[0].array(mfi);
508 Array4<Real> const& yip = m_edgecent[1].array(mfi);
509 Array4<Real> const& zip = m_edgecent[2].array(mfi);
510
511 if (iter == 0) {
512#ifdef AMREX_USE_GPU
513 if (hybrid) {
515 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
516 edgetype[idim].prefetchToHost();
517 m_edgecent[idim][mfi].prefetchToHost();
518 }
519 }
520#endif
521
522 gshop.getIntercept({xip,yip,zip}, {xdg,ydg,zdg}, clst,
523 geom, gshop_run_on, bounding_box);
524
525#ifdef AMREX_USE_GPU
526 if (hybrid) {
527 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
528 edgetype[idim].prefetchToDevice();
529 m_edgecent[idim][mfi].prefetchToDevice();
530 }
531 }
532#endif
533 } else {
534 gshop.updateIntercept({xip,yip,zip}, {xdg,ydg,zdg}, clst, geom);
535 }
536
537 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
538 const Box& b = facetype[idim].box();
539 M2[idim].resize(b,3);
540 }
541 Array4<Real> const& xm2 = M2[0].array();
542 Array4<Real> const& ym2 = M2[1].array();
543 Array4<Real> const& zm2 = M2[2].array();
544
545 nmc = build_faces(vbx, cfg, ftx, fty, ftz, xdg, ydg, zdg, lst,
546 xip, yip, zip, apx, apy, apz, fcx, fcy, fcz,
547 xm2, ym2, zm2, dx, problo, cover_multiple_cuts);
548
549 cellflagtmp.resize(m_cellflag[mfi].box());
550 Elixir cellflagtmp_eli = cellflagtmp.elixir();
551 Array4<EBCellFlag> const& cfgtmp = cellflagtmp.array();
552
553 build_cells(vbx, cfg, ftx, fty, ftz, apx, apy, apz,
554 fcx, fcy, fcz, xm2, ym2, zm2, dx, vfr, ctr,
555 bar, bct, bnm, cfgtmp, lst,
556 small_volfrac, geom, extend_domain_face, cover_multiple_cuts,
557 nsm, nmc);
558
559 // Because it is used in a synchronous reduction kernel in
560 // build_cells, we do not need to worry about M2's lifetime.
561 // But we still need to use Elixir to extend the life of
562 // cellflagtmp.
563
564#elif (AMREX_SPACEDIM == 2)
565 Array4<Real> const& xip = m_edgecent[0].array(mfi);
566 Array4<Real> const& yip = m_edgecent[1].array(mfi);
567
568 if (iter == 0) {
569#ifdef AMREX_USE_GPU
570 if (hybrid) {
572 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
573 facetype[idim].prefetchToHost();
574 m_edgecent[idim][mfi].prefetchToHost();
575 }
576 }
577#endif
578
579 // yes, factype[1] and then [0]
580 gshop.getIntercept({xip,yip},
581 {facetype[1].const_array(), facetype[0].const_array()},
582 clst, geom, gshop_run_on, bounding_box);
583
584#ifdef AMREX_USE_GPU
585 if (hybrid) {
586 for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
587 facetype[idim].prefetchToDevice();
588 m_edgecent[idim][mfi].prefetchToDevice();
589 }
590 }
591#endif
592 } else {
593 gshop.updateIntercept({xip,yip},
594 {facetype[1].const_array(), facetype[0].const_array()},
595 clst, geom);
596 }
597
598 nmc = build_faces(vbx, cfg, ftx, fty, lst, xip, yip, apx, apy, fcx, fcy,
599 dx, problo, cover_multiple_cuts, nsm);
600
601 build_cells(vbx, cfg, ftx, fty, apx, apy, dx, vfr, ctr,
602 bar, bct, bnm, lst, small_volfrac, geom, extend_domain_face,
603 nsm, nmc);
604#endif
605 nsmallcells += nsm;
606 nmulticuts += nmc;
607 }
608 }
609
610 ParallelAllReduce::Sum<int>({nsmallcells,nmulticuts}, ParallelContext::CommunicatorSub());
611 if (nsmallcells == 0 && nmulticuts == 0) {
612 break;
613 } else {
614 auto ls = m_mgf.getLevelSet();
615 // This is an alias MultiFab, therefore FillBoundary on it is fine.
616 ls.FillBoundary(geom.periodicity());
617 if (amrex::Verbose() > 0) {
618 if (nsmallcells) {
619 amrex::Print() << "AMReX EB: Iter. " << iter+1 << " fixed " << nsmallcells
620 << " small cells" << '\n';
621 }
622 if (nmulticuts) {
623 amrex::Print() << "AMReX EB: Iter. " << iter+1 << " fixed " << nmulticuts
624 << " multicuts" << '\n';
625 }
626 }
627 }
628 }
629
630 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(iter < maxiter, "EB: failed to fix small cells");
631
632#ifdef AMREX_USE_OMP
633#pragma omp parallel if (Gpu::notInLaunchRegion())
634#endif
635 for (MFIter mfi(m_mgf); mfi.isValid(); ++mfi)
636 {
637 auto& gfab = m_mgf[mfi];
638 auto const& levelset = gfab.getLevelSet();
639 Array4<Real const> const& clst = levelset.const_array();
640
641 AMREX_D_TERM(Array4<Real> const& xip = m_edgecent[0].array(mfi);,
642 Array4<Real> const& yip = m_edgecent[1].array(mfi);,
643 Array4<Real> const& zip = m_edgecent[2].array(mfi);)
644#if (AMREX_SPACEDIM == 3)
645 auto const& edgetype = gfab.getEdgeType();
646 Array4<Type_t const> const& xdg = edgetype[0].const_array();
647 Array4<Type_t const> const& ydg = edgetype[1].const_array();
648 Array4<Type_t const> const& zdg = edgetype[2].const_array();
649
650 intercept_to_edge_centroid(xip, yip, zip, xdg, ydg, zdg, clst, dx, problo);
651
652#elif (AMREX_SPACEDIM == 2)
653 auto& facetype = gfab.getFaceType();
654 Array4<Type_t const> const& ftx = facetype[0].const_array();
655 Array4<Type_t const> const& fty = facetype[1].const_array();
656 // fty then ftx
657 intercept_to_edge_centroid(xip, yip, fty, ftx, clst, dx, problo);
658#endif
659 }
660
661 m_levelset = m_mgf.getLevelSet();
662
663 m_ok = true;
664}
665
666#if (AMREX_SPACEDIM == 3)
667template <typename G>
668template <typename GS>
669requires (std::same_as<GS,STLtools>)
670void
671GShopLevel<G>::define_fine_mvmc (GS const& gshop, const Geometry& geom,
672 int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
673{
674 BL_PROFILE("EB2::GShopLevel()-fine-mvmc");
675
676 prepare_grids(gshop, geom, max_grid_size, ngrow, extend_domain_face, num_crse_opt);
677
678 if (m_allregular) { return; }
679
681
682 m_sdf.define(amrex::convert(m_grids,IntVect(1)), m_dmap, 1, IntVect(1));
683 gshop.fillSignedDistance(m_sdf, m_sdf.nGrowVect(), geom);
684
685 for (MFIter mfi(m_sdf,MFItInfo().DisableDeviceSync()); mfi.isValid(); ++mfi) {
686 m_marching_cubes[mfi.index()] = std::make_unique<MC::MCFab>();
687 }
688
689 MFItInfo info{};
690
691#if defined(AMREX_USE_OMP) && !defined(AMREX_USE_GPU)
692 info.SetDynamic(true);
693#pragma omp parallel
694#endif
695 for (MFIter mfi(m_sdf,info); mfi.isValid(); ++mfi) {
696 amrex::MC::marching_cubes(geom, m_sdf[mfi], *m_marching_cubes[mfi.index()]);
697 }
698
699 MC::write_stl("test.stl", m_marching_cubes);
700
701 m_ok = false; // xxxxx TODO
702}
703#endif
704
705template <typename G>
706GShopLevel<G>::GShopLevel (IndexSpace const* is, int /*ilev*/, int max_grid_size, int /*ngrow*/,
707 const Geometry& geom, GShopLevel<G>& fineLevel)
708 : Level(is, geom)
709{
710 if (fineLevel.isAllRegular()) {
711 m_allregular = true;
712 m_ok = true;
713 return;
714 }
715
716 BL_PROFILE("EB2::GShopLevel()-coarse");
717
718 const BoxArray& fine_grids = fineLevel.m_grids;
719 const BoxArray& fine_covered_grids = fineLevel.m_covered_grids;
720
721 const int coarse_ratio = 2;
722 const int min_width = 8;
723 bool coarsenable = fine_grids.coarsenable(coarse_ratio, min_width)
724 && (fine_covered_grids.empty() || fine_covered_grids.coarsenable(coarse_ratio));
725
726 m_ngrow = amrex::coarsen(fineLevel.m_ngrow,2);
727 if (amrex::scale(m_ngrow,2) != fineLevel.m_ngrow) {
729 }
730
731 if (coarsenable)
732 {
733 int ierr = coarsenFromFine(fineLevel, true);
734 m_ok = (ierr == 0);
735 }
736 else
737 {
738 Level fine_level_2(is, fineLevel.m_geom);
739 fine_level_2.prepareForCoarsening(fineLevel, max_grid_size, amrex::scale(m_ngrow,2));
740 int ierr = coarsenFromFine(fine_level_2, false);
741 m_ok = (ierr == 0);
742 }
743}
744
745}
746
747#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:1633
Array4< T const > array() const noexcept
Definition AMReX_BaseFab.H:387
Elixir elixir() noexcept
Definition AMReX_BaseFab.H:1676
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:564
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:617
__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:901
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:26
Definition AMReX_EB2_Level.H:203
GShopLevel(IndexSpace const *is, const Geometry &geom)
Construct an empty level (usually regular) bound to geom.
Definition AMReX_EB2_Level.H:252
void define_fine(G const &gshop, const Geometry &geom, int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
Define data for the finest level using the supplied GeometryShop.
Definition AMReX_EB2_Level.H:387
GShopLevel(IndexSpace const *is, int ilev, int max_grid_size, int ngrow, const Geometry &geom, GShopLevel< G > &fineLevel)
Construct a coarse level by coarsening fineLevel.
Definition AMReX_EB2_Level.H:706
GShopLevel(IndexSpace const *is, G const &gshop, const Geometry &geom, int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
Build a level directly from a GeometryShop object.
Definition AMReX_EB2_Level.H:257
static GShopLevel< G > makeAllRegular(IndexSpace const *is, const Geometry &geom)
Build a regular (no EB) level.
Definition AMReX_EB2_Level.H:235
void define_fine_mvmc(GS const &gshop, const Geometry &geom, int max_grid_size, int ngrow, bool extend_domain_face, int num_crse_opt)
Define a fine level for multi-valued multi-cut STL geometries (3-D only)
Definition AMReX_EB2_Level.H:671
Definition AMReX_EB2.H:28
Definition AMReX_EB2_Level.H:42
MultiFab m_bndryarea
Definition AMReX_EB2_Level.H:176
bool isOK() const noexcept
Definition AMReX_EB2_Level.H:48
void fillFaceCent(Array< MultiCutFab *, 3 > const &facecent, const Geometry &geom) const
Populate face centroids for each direction.
Definition AMReX_EB2_Level.cpp:794
MultiFab m_volfrac
Definition AMReX_EB2_Level.H:174
MultiGFab m_mgf
Definition AMReX_EB2_Level.H:171
bool m_allregular
Definition AMReX_EB2_Level.H:183
IntVect m_shift
Definition AMReX_EB2_Level.H:186
iMultiFab m_cutcellmask
Definition AMReX_EB2_Level.H:182
Array< MultiFab, 3 > m_facecent
Definition AMReX_EB2_Level.H:180
Array< MultiFab, 3 > m_edgecent
Definition AMReX_EB2_Level.H:181
void prepareForCoarsening(const Level &rhs, int max_grid_size, IntVect const &ngrow)
Prepare for coarsening by copying from rhs, respecting max grid size and grow vector.
Definition AMReX_EB2_Level.cpp:10
void fillCutCellMask(iMultiFab &cutcellmask, const Geometry &geom) const
Populate the cut-cell mask into cutcellmask.
Definition AMReX_EB2_Level.cpp:946
MultiFab m_centroid
Definition AMReX_EB2_Level.H:175
FabArray< EBCellFlagFab > m_cellflag
Definition AMReX_EB2_Level.H:173
const DistributionMapping & DistributionMap() const noexcept
Definition AMReX_EB2_Level.H:110
BoxArray m_covered_grids
Definition AMReX_EB2_Level.H:165
void fillLevelSet(MultiFab &levelset, const Geometry &geom) const
Write the implicit function values into levelset.
Definition AMReX_EB2_Level.cpp:907
Geometry m_geom
Definition AMReX_EB2_Level.H:162
Array< MultiFab, 3 > m_areafrac
Definition AMReX_EB2_Level.H:179
IntVect m_ngrow
Definition AMReX_EB2_Level.H:163
IndexSpace const * m_parent
Definition AMReX_EB2_Level.H:187
BoxArray m_grids
Definition AMReX_EB2_Level.H:164
const Geometry & Geom() const noexcept
Definition AMReX_EB2_Level.H:120
std::map< int, std::unique_ptr< MC::MCFab > > m_marching_cubes
Definition AMReX_EB2_Level.H:169
bool hasEBInfo() const noexcept
Definition AMReX_EB2_Level.H:137
void write_to_chkpt_file(const std::string &fname, bool extend_domain_face, int max_grid_size) const
Write this level’s EB data to a checkpoint file fname.
Definition AMReX_EB2_Level.cpp:956
void fillEBCellFlag(FabArray< EBCellFlagFab > &cellflag, const Geometry &geom) const
Fill cellflag with EBCellFlag data using geometry geom.
Definition AMReX_EB2_Level.cpp:429
void fillEdgeCent(Array< MultiCutFab *, 3 > const &edgecent, const Geometry &geom) const
Populate edge centroids for each direction.
Definition AMReX_EB2_Level.cpp:835
int coarsenFromFine(Level &fineLevel, bool fill_boundary)
Coarsen EB data from fineLevel into this level.
Definition AMReX_EB2_Level.cpp:84
MultiFab m_bndrynorm
Definition AMReX_EB2_Level.H:178
bool m_has_eb_info
Definition AMReX_EB2_Level.H:185
void fillAreaFrac(Array< MultiCutFab *, 3 > const &areafrac, const Geometry &geom) const
Populate face area fractions for each direction.
Definition AMReX_EB2_Level.cpp:640
MultiFab m_bndrycent
Definition AMReX_EB2_Level.H:177
MultiFab m_levelset
Definition AMReX_EB2_Level.H:172
void setShift(int direction, int ncells)
Shift this level by ncells cells along direction direction (for periodic tiling).
Definition AMReX_EB2_Level.cpp:1047
bool isAllRegular() const noexcept
Definition AMReX_EB2_Level.H:46
void fillBndryCent(MultiCutFab &bndrycent, const Geometry &geom) const
Populate boundary centroids into bndrycent.
Definition AMReX_EB2_Level.cpp:592
void fillBndryNorm(MultiCutFab &bndrynorm, const Geometry &geom) const
Populate boundary normals into bndrynorm.
Definition AMReX_EB2_Level.cpp:616
MultiFab m_sdf
Definition AMReX_EB2_Level.H:167
IntVect const & nGrowVect() const noexcept
Definition AMReX_EB2_Level.H:125
void fillVolFrac(MultiFab &vfrac, const Geometry &geom) const
Populate volume fractions into vfrac using geom.
Definition AMReX_EB2_Level.cpp:483
const BoxArray & boxArray() const noexcept
Definition AMReX_EB2_Level.H:108
Level(IndexSpace const *is, const Geometry &geom)
Construct a Level bound to EB index space is and Geometry geom.
Definition AMReX_EB2_Level.H:113
friend class GShopLevel
Definition AMReX_EB2_Level.H:190
bool m_ok
Definition AMReX_EB2_Level.H:184
void fillCentroid(MultiCutFab &centroid, const Geometry &geom) const
Populate cut-cell centroids into a MultiCutFab.
Definition AMReX_EB2_Level.cpp:544
DistributionMapping m_dmap
Definition AMReX_EB2_Level.H:166
void fillBndryArea(MultiCutFab &bndryarea, const Geometry &geom) const
Populate boundary-face areas into bndryarea.
Definition AMReX_EB2_Level.cpp:568
void buildCellFlag()
Build EBCellFlag data from the current implicit function field.
Definition AMReX_EB2_Level.cpp:403
void buildCutCellMask(Level const &fine_level)
Build the cut-cell mask by coarsening from fine_level.
Definition AMReX_EB2_Level.cpp:966
IndexSpace const * getEBIndexSpace() const noexcept
Definition AMReX_EB2_Level.H:122
LayoutData wrapper that allocates a GFab per FAB/Box.
Definition AMReX_EB2_MultiGFab.H:86
FAB storing EBCellFlag data per cell.
Definition AMReX_EBCellFlag.H:324
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:344
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:75
const Box & Domain() const noexcept
Returns our rectangular domain.
Definition AMReX_Geometry.H:216
Periodicity periodicity() const noexcept
Definition AMReX_Geometry.H:361
GpuArray< Real, 3 > ProbLoArray() const noexcept
Definition AMReX_Geometry.H:192
bool isPeriodic(int dir) const noexcept
Is the domain periodic in the specified direction?
Definition AMReX_Geometry.H:337
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:771
__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:790
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
FabArray façade that only allocates cut-cell Fabs (skipping regular cells).
Definition AMReX_MultiCutFab.H:87
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:351
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:1047
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:1418
__host__ __device__ BoxND< dim > grow(const BoxND< dim > &b, int i) noexcept
Grow BoxND in all directions by given amount.
Definition AMReX_Box.H:1289
__host__ __device__ BoxND< dim > refine(const BoxND< dim > &b, int ref_ratio) noexcept
Definition AMReX_Box.H:1468
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:310
bool inLaunchRegion() noexcept
Definition AMReX_GpuControl.H:88
void write_stl(std::string const &filename, std::map< int, std::unique_ptr< MCFab > > const &mc_fabs)
Write the collected marching-cubes output to an STL file.
Definition AMReX_MarchingCubes.cpp:924
void marching_cubes(Geometry const &geom, FArrayBox &sdf_fab, MCFab &mc_fab)
Run marching cubes on signed-distance field sdf_fab.
Definition AMReX_MarchingCubes.cpp:693
void Initialize()
Initialize internal lookup tables and device buffers for marching cubes.
Definition AMReX_MarchingCubes.cpp:651
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
Definition AMReX_Amr.cpp:50
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1567
__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:1531
RunOn
Definition AMReX_GpuControl.H:65
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:1102
int Verbose() noexcept
Definition AMReX.cpp:181
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241
void AllGatherBoxes(Vector< Box > &bxs, int n_extra_reserve)
Definition AMReX_Box.cpp:124
A multidimensional array accessor.
Definition AMReX_Array4.H:285
FabArray memory allocation information.
Definition AMReX_FabArray.H:68
MFInfo & SetTag() noexcept
Definition AMReX_FabArray.H:81
Definition AMReX_MFIter.H:20
MFItInfo & SetDynamic(bool f) noexcept
Definition AMReX_MFIter.H:43