def test_chained_full_resample_matrix(actx_factory, ndim, visualize=False): from meshmode.discretization.connection.chained import \ make_full_resample_matrix actx = actx_factory() discr = create_discretization(actx, ndim, order=2, nelements=12) connections = [] conn = create_refined_connection(actx, discr) connections.append(conn) conn = create_refined_connection(actx, conn.to_discr) connections.append(conn) from meshmode.discretization.connection import \ ChainedDiscretizationConnection chained = ChainedDiscretizationConnection(connections) def f(x): from functools import reduce return 0.1 * reduce(lambda x, y: x * actx.np.sin(5 * y), x) resample_mat = actx.to_numpy(make_full_resample_matrix(actx, chained)) x = thaw(connections[0].from_discr.nodes(), actx) fx = f(x) f1 = resample_mat @ flatten_to_numpy(actx, fx) f2 = flatten_to_numpy(actx, chained(fx)) f3 = flatten_to_numpy(actx, connections[1](connections[0](fx))) assert np.allclose(f1, f2) assert np.allclose(f2, f3)
def test_flatten_unflatten(actx_factory): actx = actx_factory() ambient_dim = 2 from meshmode.mesh.generation import generate_regular_rect_mesh mesh = generate_regular_rect_mesh(a=(-0.5, ) * ambient_dim, b=(+0.5, ) * ambient_dim, n=(3, ) * ambient_dim, order=1) discr = Discretization(actx, mesh, PolynomialWarpAndBlendGroupFactory(3)) a = np.random.randn(discr.ndofs) from meshmode.dof_array import flatten, unflatten a_round_trip = actx.to_numpy( flatten(unflatten(actx, discr, actx.from_numpy(a)))) assert np.array_equal(a, a_round_trip) from meshmode.dof_array import flatten_to_numpy, unflatten_from_numpy a_round_trip = flatten_to_numpy(actx, unflatten_from_numpy(actx, discr, a)) assert np.array_equal(a, a_round_trip) x = thaw(discr.nodes(), actx) avg_mass = DOFArray( actx, tuple([(np.pi + actx.zeros((grp.nelements, 1), a.dtype)) for grp in discr.groups])) c = MyContainer(name="flatten", mass=avg_mass, momentum=make_obj_array([x, x, x]), enthalpy=x) from meshmode.dof_array import unflatten_like c_round_trip = unflatten_like(actx, flatten(c), c) assert flat_norm(c - c_round_trip) < 1.0e-8
def test_orientation_3d(actx_factory, what, mesh_gen_func, visualize=False): pytest.importorskip("pytential") logging.basicConfig(level=logging.INFO) actx = actx_factory() mesh = mesh_gen_func() logger.info("%d elements", mesh.nelements) from meshmode.discretization import Discretization discr = Discretization(actx, mesh, PolynomialWarpAndBlendGroupFactory(3)) from pytential import bind, sym # {{{ check normals point outward if what == "torus": nodes = sym.nodes(mesh.ambient_dim).as_vector() angle = sym.atan2(nodes[1], nodes[0]) center_nodes = sym.make_obj_array( [5 * sym.cos(angle), 5 * sym.sin(angle), 0 * angle]) normal_outward_expr = (sym.normal(mesh.ambient_dim) | (nodes - center_nodes)) else: normal_outward_expr = (sym.normal(mesh.ambient_dim) | sym.nodes(mesh.ambient_dim)) normal_outward_check = flatten_to_numpy( actx, bind(discr, normal_outward_expr)(actx).as_scalar()) > 0 assert normal_outward_check.all(), normal_outward_check # }}} normals = bind(discr, sym.normal(mesh.ambient_dim).xproject(1))(actx) if visualize: from meshmode.discretization.visualization import make_visualizer vis = make_visualizer(actx, discr, 3) vis.write_vtk_file("orientation_3d_%s_normals.vtu" % what, [ ("normals", normals), ])
def test_boundary_interpolation(actx_factory, group_factory, boundary_tag, mesh_name, dim, mesh_pars, per_face_groups): if (group_factory is LegendreGaussLobattoTensorProductGroupFactory and mesh_name == "blob"): pytest.skip("tensor products not implemented on blobs") actx = actx_factory() if group_factory is LegendreGaussLobattoTensorProductGroupFactory: group_cls = TensorProductElementGroup else: group_cls = SimplexElementGroup from meshmode.discretization import Discretization from meshmode.discretization.connection import (make_face_restriction, check_connection) from pytools.convergence import EOCRecorder eoc_rec = EOCRecorder() order = 4 def f(x): return 0.1 * actx.np.sin(30 * x) for mesh_par in mesh_pars: # {{{ get mesh if mesh_name == "blob": assert dim == 2 h = float(mesh_par) #from meshmode.mesh.io import generate_gmsh, FileSource # print("BEGIN GEN") # mesh = generate_gmsh( # FileSource("blob-2d.step"), 2, order=order, # force_ambient_dim=2, # other_options=[ # "-string", "Mesh.CharacteristicLengthMax = %s;" % h] # ) # print("END GEN") from meshmode.mesh.io import read_gmsh mesh = read_gmsh("blob2d-order%d-h%s.msh" % (order, mesh_par), force_ambient_dim=2) elif mesh_name == "warp": mesh = mgen.generate_warped_rect_mesh(dim, order=order, nelements_side=mesh_par, group_cls=group_cls) h = 1 / mesh_par elif mesh_name == "rect": mesh = mgen.generate_regular_rect_mesh( a=(0, ) * dim, b=(1, ) * dim, order=order, nelements_per_axis=(mesh_par, ) * dim, group_cls=group_cls) h = 1 / mesh_par else: raise ValueError("mesh_name not recognized") # }}} vol_discr = Discretization(actx, mesh, group_factory(order)) print("h=%s -> %d elements" % (h, sum(mgrp.nelements for mgrp in mesh.groups))) x = thaw(vol_discr.nodes()[0], actx) vol_f = f(x) bdry_connection = make_face_restriction( actx, vol_discr, group_factory(order), boundary_tag, per_face_groups=per_face_groups) check_connection(actx, bdry_connection) bdry_discr = bdry_connection.to_discr bdry_x = thaw(bdry_discr.nodes()[0], actx) bdry_f = f(bdry_x) bdry_f_2 = bdry_connection(vol_f) if mesh_name == "blob" and dim == 2 and mesh.nelements < 500: from meshmode.discretization.connection.direct import \ make_direct_full_resample_matrix mat = actx.to_numpy( make_direct_full_resample_matrix(actx, bdry_connection)) bdry_f_2_by_mat = mat.dot(flatten_to_numpy(actx, vol_f)) mat_error = la.norm( flatten_to_numpy(actx, bdry_f_2) - bdry_f_2_by_mat) assert mat_error < 1e-14, mat_error err = flat_norm(bdry_f - bdry_f_2, np.inf) eoc_rec.add_data_point(h, err) order_slack = 0.75 if mesh_name == "blob" else 0.5 print(eoc_rec) assert (eoc_rec.order_estimate() >= order - order_slack or eoc_rec.max_error() < 3.6e-13)
def test_sanity_balls(actx_factory, src_file, dim, mesh_order, visualize=False): pytest.importorskip("pytential") logging.basicConfig(level=logging.INFO) actx = actx_factory() from pytools.convergence import EOCRecorder vol_eoc_rec = EOCRecorder() surf_eoc_rec = EOCRecorder() # overkill quad_order = mesh_order from pytential import bind, sym for h in [0.2, 0.1, 0.05]: from meshmode.mesh.io import generate_gmsh, FileSource mesh = generate_gmsh(FileSource(src_file), dim, order=mesh_order, other_options=[ "-string", "Mesh.CharacteristicLengthMax = %g;" % h ], force_ambient_dim=dim, target_unit="MM") logger.info("%d elements", mesh.nelements) # {{{ discretizations and connections from meshmode.discretization import Discretization vol_discr = Discretization( actx, mesh, InterpolatoryQuadratureSimplexGroupFactory(quad_order)) from meshmode.discretization.connection import make_face_restriction bdry_connection = make_face_restriction( actx, vol_discr, InterpolatoryQuadratureSimplexGroupFactory(quad_order), BTAG_ALL) bdry_discr = bdry_connection.to_discr # }}} from math import gamma true_surf = 2 * np.pi**(dim / 2) / gamma(dim / 2) true_vol = true_surf / dim vol_x = thaw(vol_discr.nodes(), actx) vol_one = vol_x[0] * 0 + 1 from pytential import norm, integral # noqa comp_vol = integral(vol_discr, vol_one) rel_vol_err = abs(true_vol - comp_vol) / true_vol vol_eoc_rec.add_data_point(h, rel_vol_err) print("VOL", true_vol, comp_vol) bdry_x = thaw(bdry_discr.nodes(), actx) bdry_one_exact = bdry_x[0] * 0 + 1 bdry_one = bdry_connection(vol_one) intp_err = norm(bdry_discr, bdry_one - bdry_one_exact) assert intp_err < 1e-14 comp_surf = integral(bdry_discr, bdry_one) rel_surf_err = abs(true_surf - comp_surf) / true_surf surf_eoc_rec.add_data_point(h, rel_surf_err) print("SURF", true_surf, comp_surf) if visualize: from meshmode.discretization.visualization import make_visualizer vol_vis = make_visualizer(actx, vol_discr, 7) bdry_vis = make_visualizer(actx, bdry_discr, 7) name = src_file.split("-")[0] vol_vis.write_vtk_file(f"sanity_balls_volume_{name}_{h:g}.vtu", [ ("f", vol_one), ("area_el", bind(vol_discr, sym.area_element(mesh.ambient_dim, mesh.ambient_dim))(actx)), ]) bdry_vis.write_vtk_file(f"sanity_balls_boundary_{name}_{h:g}.vtu", [("f", bdry_one)]) # {{{ check normals point outward normal_outward_check = bind( bdry_discr, sym.normal(mesh.ambient_dim) | sym.nodes(mesh.ambient_dim), )(actx).as_scalar() normal_outward_check = flatten_to_numpy(actx, normal_outward_check > 0) assert normal_outward_check.all(), normal_outward_check # }}} print("---------------------------------") print("VOLUME") print("---------------------------------") print(vol_eoc_rec) assert vol_eoc_rec.order_estimate() >= mesh_order print("---------------------------------") print("SURFACE") print("---------------------------------") print(surf_eoc_rec) assert surf_eoc_rec.order_estimate() >= mesh_order
def test_sanity_single_element(actx_factory, dim, mesh_order, group_cls, visualize=False): pytest.importorskip("pytential") actx = actx_factory() if group_cls is SimplexElementGroup: group_factory = PolynomialWarpAndBlendGroupFactory(mesh_order + 3) elif group_cls is TensorProductElementGroup: group_factory = LegendreGaussLobattoTensorProductGroupFactory( mesh_order + 3) else: raise TypeError import modepy as mp shape = group_cls._modepy_shape_cls(dim) space = mp.space_for_shape(shape, mesh_order) vertices = mp.unit_vertices_for_shape(shape) nodes = mp.edge_clustered_nodes_for_space(space, shape).reshape(dim, 1, -1) vertex_indices = np.arange(shape.nvertices, dtype=np.int32).reshape(1, -1) center = np.empty(dim, np.float64) center.fill(-0.5) mg = group_cls(mesh_order, vertex_indices, nodes, dim=dim) mesh = Mesh(vertices, [mg], is_conforming=True) from meshmode.discretization import Discretization vol_discr = Discretization(actx, mesh, group_factory) # {{{ volume calculation check if isinstance(mg, SimplexElementGroup): from pytools import factorial true_vol = 1 / factorial(dim) * 2**dim elif isinstance(mg, TensorProductElementGroup): true_vol = 2**dim else: raise TypeError nodes = thaw(vol_discr.nodes(), actx) vol_one = 1 + 0 * nodes[0] from pytential import norm, integral # noqa comp_vol = integral(vol_discr, vol_one) rel_vol_err = abs(true_vol - comp_vol) / true_vol assert rel_vol_err < 1e-12 # }}} # {{{ boundary discretization from meshmode.discretization.connection import make_face_restriction bdry_connection = make_face_restriction(actx, vol_discr, group_factory, BTAG_ALL) bdry_discr = bdry_connection.to_discr # }}} from pytential import bind, sym bdry_normals = bind(bdry_discr, sym.normal(dim).as_vector())(actx) if visualize: from meshmode.discretization.visualization import make_visualizer bdry_vis = make_visualizer(actx, bdry_discr, 4) bdry_vis.write_vtk_file("sanity_single_element_boundary.vtu", [("normals", bdry_normals)]) normal_outward_check = bind( bdry_discr, sym.normal(dim) | (sym.nodes(dim) + 0.5 * sym.ones_vec(dim)), )(actx).as_scalar() normal_outward_check = flatten_to_numpy(actx, normal_outward_check > 0) assert normal_outward_check.all(), normal_outward_check
def test_chained_to_direct(actx_factory, ndim, chain_type, nelements=128, visualize=False): import time from meshmode.discretization.connection.chained import \ flatten_chained_connection actx = actx_factory() discr = create_discretization(actx, ndim, nelements=nelements) connections = [] if chain_type == 1: conn = create_refined_connection(actx, discr) connections.append(conn) conn = create_refined_connection(actx, conn.to_discr) connections.append(conn) elif chain_type == 2: conn = create_refined_connection(actx, discr) connections.append(conn) conn = create_refined_connection(actx, conn.to_discr) connections.append(conn) conn = create_refined_connection(actx, conn.to_discr) connections.append(conn) elif chain_type == 3 and ndim == 3: conn = create_refined_connection(actx, discr, threshold=np.inf) connections.append(conn) conn = create_face_connection(actx, conn.to_discr) connections.append(conn) else: raise ValueError("unknown test case") from meshmode.discretization.connection import \ ChainedDiscretizationConnection chained = ChainedDiscretizationConnection(connections) t_start = time.time() direct = flatten_chained_connection(actx, chained) t_end = time.time() if visualize: print("[TIME] Flatten: {:.5e}".format(t_end - t_start)) if chain_type < 3: to_element_indices = np.full(direct.to_discr.mesh.nelements, 0, dtype=np.int64) for grp in direct.groups: for batch in grp.batches: for i in batch.to_element_indices.get(actx.queue): to_element_indices[i] += 1 assert np.min(to_element_indices) > 0 def f(x): from functools import reduce return 0.1 * reduce(lambda x, y: x * actx.np.sin(5 * y), x) x = thaw(connections[0].from_discr.nodes(), actx) fx = f(x) t_start = time.time() f1 = flatten_to_numpy(actx, direct(fx)) t_end = time.time() if visualize: print("[TIME] Direct: {:.5e}".format(t_end - t_start)) t_start = time.time() f2 = flatten_to_numpy(actx, chained(fx)) t_end = time.time() if visualize: print("[TIME] Chained: {:.5e}".format(t_end - t_start)) if visualize and ndim == 2: import matplotlib.pyplot as pt pt.figure(figsize=(10, 8), dpi=300) pt.plot(f1, label="Direct") pt.plot(f2, label="Chained") pt.ylim([np.min(f2) - 0.1, np.max(f2) + 0.1]) pt.legend() pt.savefig("test_chained_to_direct.png") pt.clf() assert np.allclose(f1, f2)