|
1 | 1 | # -*- coding: utf-8 -*- |
2 | 2 |
|
| 3 | +import numpy |
3 | 4 | import pytest |
4 | 5 |
|
5 | 6 | import amrex.space3d as amr |
6 | 7 |
|
7 | 8 |
|
8 | 9 | @pytest.mark.skipif(not amr.Config.have_eb, reason="Requires -DAMReX_EB=ON") |
9 | | -def test_makeEBFab(): |
10 | | - pass |
| 10 | +def test_makeEBFabFactory(): |
| 11 | + n_cell = 64 |
| 12 | + max_grid_size = 16 |
11 | 13 |
|
12 | | - # TODO: |
13 | | - # EB2_Build(...) |
14 | | - # makeEBFabFactory(...) |
| 14 | + # Build Geometry |
| 15 | + domain = amr.Box(amr.IntVect(0,0,0), amr.IntVect(n_cell-1, n_cell-1, n_cell-1)) |
| 16 | + real_box = amr.RealBox([0., 0., 0.], [1. , 1., 1.]) |
| 17 | + coord = 0 # Cartesian |
| 18 | + is_per = [1, 1, 1] # is periodic? |
| 19 | + geom = amr.Geometry(domain, real_box, coord, is_per) |
| 20 | + |
| 21 | + # EB parameters |
| 22 | + pp = amr.ParmParse("eb2") |
| 23 | + pp.add("geom_type", "sphere") |
| 24 | + pp.addarr("sphere_center", [0.5, 0.5, 0.5]) |
| 25 | + rsphere = 0.25 |
| 26 | + pp.add("sphere_radius", rsphere) |
| 27 | + pp.add("sphere_has_fluid_inside", 1) |
| 28 | + |
| 29 | + # EB generation |
| 30 | + eb_requried_level = 0 |
| 31 | + eb_max_level = 2 |
| 32 | + amr.EB2_Build(geom, eb_requried_level, eb_max_level) |
| 33 | + |
| 34 | + # Build BoxArray |
| 35 | + ba = amr.BoxArray(domain) |
| 36 | + ba.max_size(max_grid_size) |
| 37 | + |
| 38 | + # Build DistributionMapping |
| 39 | + dm = amr.DistributionMapping(ba) |
| 40 | + |
| 41 | + # Make EB Factory |
| 42 | + ng = amr.Vector_int([1,1,1]) |
| 43 | + factory = amr.makeEBFabFactory(geom, ba, dm, ng, amr.EBSupport.full) |
| 44 | + |
| 45 | + # Get EB data |
| 46 | + vfrac = factory.getVolFrac() |
| 47 | + |
| 48 | + dx = geom.data().CellSize() |
| 49 | + total_vol = vfrac.sum()*dx[0]*dx[1]*dx[2] |
| 50 | + sphere_vol = 4./3.*numpy.pi*rsphere**3 |
| 51 | + assert abs(sphere_vol-total_vol)/sphere_vol < 2.e-3 |
0 commit comments