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