예제 #1
0
def _visualize_refinement(actx: PyOpenCLArrayContext,
                          discr,
                          niter,
                          stage_nr,
                          stage_name,
                          flags,
                          visualize=False):
    if not visualize:
        return

    if stage_nr not in (1, 2):
        raise ValueError("unexpected stage number")

    flags = actx.to_numpy(flags)
    logger.info("for stage %s: splitting %d/%d stage-%d elements", stage_name,
                np.sum(flags), discr.mesh.nelements, stage_nr)

    from meshmode.discretization.visualization import make_visualizer
    vis = make_visualizer(actx, discr, 3)

    assert len(flags) == discr.mesh.nelements

    flags = flags.astype(bool)
    nodes_flags_template = discr.zeros(actx)
    nodes_flags = []
    for grp in discr.groups:
        meg = grp.mesh_el_group
        nodes_flags_grp = actx.to_numpy(nodes_flags_template[grp.index])
        nodes_flags_grp[flags[meg.element_nr_base:meg.nelements +
                              meg.element_nr_base]] = 1
        nodes_flags.append(actx.from_numpy(nodes_flags_grp))

    nodes_flags = DOFArray(actx, tuple(nodes_flags))

    vis_data = [
        ("refine_flags", nodes_flags),
    ]

    if 0:
        from pytential import sym, bind
        bdry_normals = bind(discr, sym.normal(
            discr.ambient_dim))(actx).as_vector(dtype=object)
        vis_data.append(("bdry_normals", bdry_normals), )

    vis.write_vtk_file(f"refinement-{stage_name}-{niter:03d}.vtu", vis_data)
예제 #2
0
def test_to_fd_idempotency(ctx_factory, mm_mesh, fspace_degree):
    """
    Make sure mm->fd->mm and (mm->)->fd->mm->fd are identity
    """
    cl_ctx = ctx_factory()
    queue = cl.CommandQueue(cl_ctx)
    actx = PyOpenCLArrayContext(queue)

    # make sure degree is higher order than mesh
    fspace_degree += mm_mesh.groups[0].order

    # Make a function space and a function with unique values at each node
    factory = InterpolatoryQuadratureSimplexGroupFactory(fspace_degree)
    discr = Discretization(actx, mm_mesh, factory)
    fdrake_connection = build_connection_to_firedrake(discr)
    fdrake_mesh = fdrake_connection.firedrake_fspace().mesh()
    dtype = fdrake_mesh.coordinates.dat.data.dtype

    mm_unique = discr.zeros(actx, dtype=dtype)
    unique_vals = np.arange(np.size(mm_unique[0]), dtype=dtype)
    mm_unique[0].set(unique_vals.reshape(mm_unique[0].shape))
    mm_unique_copy = DOFArray(actx, (mm_unique[0].copy(), ))

    # Test for idempotency mm->fd->mm
    fdrake_unique = fdrake_connection.from_meshmode(mm_unique)
    fdrake_connection.from_firedrake(fdrake_unique, out=mm_unique_copy)

    np.testing.assert_allclose(actx.to_numpy(mm_unique_copy[0]),
                               actx.to_numpy(mm_unique[0]),
                               atol=CLOSE_ATOL)

    # Test for idempotency (mm->)fd->mm->fd
    fdrake_unique_copy = fdrake_connection.from_meshmode(mm_unique_copy)
    np.testing.assert_allclose(fdrake_unique_copy.dat.data,
                               fdrake_unique.dat.data,
                               atol=CLOSE_ATOL)
예제 #3
0
    def exec_compute_potential_insn_fmm(self, actx: PyOpenCLArrayContext,
            insn, bound_expr, evaluate, fmm_driver):
        """
        :arg fmm_driver: A function that accepts four arguments:
            *wrangler*, *strength*, *geo_data*, *kernel*, *kernel_arguments*
        :returns: a tuple ``(assignments, extra_outputs)``, where *assignments*
            is a list of tuples containing pairs ``(name, value)`` representing
            assignments to be performed in the evaluation context.
            *extra_outputs* is data that *fmm_driver* may return
            (such as timing data), passed through unmodified.
        """
        target_name_and_side_to_number, target_discrs_and_qbx_sides = (
                self.get_target_discrs_and_qbx_sides(insn, bound_expr))

        geo_data = self.qbx_fmm_geometry_data(
                bound_expr.places,
                insn.source.geometry,
                target_discrs_and_qbx_sides)

        # FIXME Exert more positive control over geo_data attribute lifetimes using
        # geo_data.<method>.clear_cache(geo_data).

        # FIXME Synthesize "bad centers" around corners and edges that have
        # inadequate QBX coverage.

        # FIXME don't compute *all* output kernels on all targets--respect that
        # some target discretizations may only be asking for derivatives (e.g.)

        flat_strengths = get_flat_strengths_from_densities(
                actx, bound_expr.places, evaluate, insn.densities,
                dofdesc=insn.source)

        base_kernel = single_valued(knl.get_base_kernel() for
            knl in insn.source_kernels)

        output_and_expansion_dtype = (
                self.get_fmm_output_and_expansion_dtype(insn.source_kernels,
                    flat_strengths[0]))
        kernel_extra_kwargs, source_extra_kwargs = (
                self.get_fmm_expansion_wrangler_extra_kwargs(
                    actx, insn.target_kernels + insn.source_kernels,
                    geo_data.tree().user_source_ids,
                    insn.kernel_arguments, evaluate))

        tree_indep = self._tree_indep_data_for_wrangler(
                target_kernels=insn.target_kernels,
                source_kernels=insn.source_kernels)

        wrangler = tree_indep.wrangler_cls(
                        tree_indep, geo_data, output_and_expansion_dtype,
                        self.qbx_order,
                        self.fmm_level_to_order,
                        source_extra_kwargs=source_extra_kwargs,
                        kernel_extra_kwargs=kernel_extra_kwargs)

        from pytential.qbx.geometry import target_state
        if actx.to_numpy(actx.np.any(
                actx.thaw(geo_data.user_target_to_center())
                == target_state.FAILED)):
            raise RuntimeError("geometry has failed targets")

        # {{{ geometry data inspection hook

        if self.geometry_data_inspector is not None:
            perform_fmm = self.geometry_data_inspector(insn, bound_expr, geo_data)
            if not perform_fmm:
                return [(o.name, 0) for o in insn.outputs]

        # }}}

        # Execute global QBX.
        all_potentials_on_every_target, extra_outputs = (
                fmm_driver(
                    wrangler, flat_strengths, geo_data,
                    base_kernel, kernel_extra_kwargs))

        results = []

        for o in insn.outputs:
            target_side_number = target_name_and_side_to_number[
                    o.target_name, o.qbx_forced_limit]
            target_discr, _ = target_discrs_and_qbx_sides[target_side_number]
            target_slice = slice(*geo_data.target_info().target_discr_starts[
                    target_side_number:target_side_number+2])

            result = \
                all_potentials_on_every_target[o.target_kernel_index][target_slice]

            from meshmode.discretization import Discretization
            if isinstance(target_discr, Discretization):
                template_ary = thaw(target_discr.nodes()[0], actx)
                result = unflatten(template_ary, result, actx, strict=False)

            results.append((o.name, result))

        return results, extra_outputs
예제 #4
0
def test_from_fd_idempotency(ctx_factory, fdrake_mesh, fspace_degree,
                             fspace_type, only_convert_bdy):
    """
    Make sure fd->mm->fd and (fd->)->mm->fd->mm are identity
    """
    # Make a function space and a function with unique values at each node
    if fspace_type == "scalar":
        fdrake_fspace = FunctionSpace(fdrake_mesh, "DG", fspace_degree)
        # Just use the node nr
        fdrake_unique = Function(fdrake_fspace)
        fdrake_unique.dat.data[:] = np.arange(fdrake_unique.dat.data.shape[0])
    elif fspace_type == "vector":
        fdrake_fspace = VectorFunctionSpace(fdrake_mesh, "DG", fspace_degree)
        # use the coordinates
        xx = SpatialCoordinate(fdrake_fspace.mesh())
        fdrake_unique = Function(fdrake_fspace).interpolate(xx)
    elif fspace_type == "tensor":
        fdrake_fspace = TensorFunctionSpace(fdrake_mesh, "DG", fspace_degree)
        # use the coordinates, duplicated into the right tensor shape
        xx = SpatialCoordinate(fdrake_fspace.mesh())
        dim = fdrake_fspace.mesh().geometric_dimension()
        unique_expr = as_tensor([xx for _ in range(dim)])
        fdrake_unique = Function(fdrake_fspace).interpolate(unique_expr)

    # Make connection
    cl_ctx = ctx_factory()
    queue = cl.CommandQueue(cl_ctx)
    actx = PyOpenCLArrayContext(queue)

    # If only converting boundary, first go ahead and do one round of
    # fd->mm->fd. This will zero out any degrees of freedom absent in
    # the meshmode mesh (because they are not associated to cells
    #                    with >= 1 node on the boundary)
    #
    # Otherwise, just continue as normal
    if only_convert_bdy:
        fdrake_connection = \
            build_connection_from_firedrake(actx,
                                            fdrake_fspace,
                                            restrict_to_boundary="on_boundary")
        temp = fdrake_connection.from_firedrake(fdrake_unique, actx=actx)
        fdrake_unique = fdrake_connection.from_meshmode(temp)
    else:
        fdrake_connection = build_connection_from_firedrake(
            actx, fdrake_fspace)

    # Test for idempotency fd->mm->fd
    mm_field = fdrake_connection.from_firedrake(fdrake_unique, actx=actx)
    fdrake_unique_copy = Function(fdrake_fspace)
    fdrake_connection.from_meshmode(mm_field, out=fdrake_unique_copy)

    np.testing.assert_allclose(fdrake_unique_copy.dat.data,
                               fdrake_unique.dat.data,
                               atol=CLOSE_ATOL)

    # Test for idempotency (fd->)mm->fd->mm
    mm_field_copy = fdrake_connection.from_firedrake(fdrake_unique_copy,
                                                     actx=actx)
    if fspace_type == "scalar":
        np.testing.assert_allclose(actx.to_numpy(mm_field_copy[0]),
                                   actx.to_numpy(mm_field[0]),
                                   atol=CLOSE_ATOL)
    else:
        for dof_arr_cp, dof_arr in zip(mm_field_copy.flatten(),
                                       mm_field.flatten()):
            np.testing.assert_allclose(actx.to_numpy(dof_arr_cp[0]),
                                       actx.to_numpy(dof_arr[0]),
                                       atol=CLOSE_ATOL)
예제 #5
0
def test_to_fd_transfer(ctx_factory, fspace_degree, mesh_name, mesh_pars, dim):
    """
    Make sure creating a function which projects onto
    one dimension then transports it is the same
    (up to resampling error) as projecting to one
    dimension on the transported mesh
    """
    # build estimate-of-convergence recorder
    from pytools.convergence import EOCRecorder
    # dimension projecting onto -> EOCRecorder
    eoc_recorders = {d: EOCRecorder() for d in range(dim)}

    # make a computing context
    cl_ctx = ctx_factory()
    queue = cl.CommandQueue(cl_ctx)
    actx = PyOpenCLArrayContext(queue)

    # Get each of the refinements of the meshmeshes and record
    # conversions errors
    for mesh_par in mesh_pars:
        if mesh_name in ("blob2d-order1", "blob2d-order4"):
            assert dim == 2
            from meshmode.mesh.io import read_gmsh
            mm_mesh = read_gmsh(f"{mesh_name}-h{mesh_par}.msh",
                                force_ambient_dim=dim)
            h = float(mesh_par)
        elif mesh_name == "warp":
            from meshmode.mesh.generation import generate_warped_rect_mesh
            mm_mesh = generate_warped_rect_mesh(dim,
                                                order=4,
                                                nelements_side=mesh_par)
            h = 1 / mesh_par
        else:
            raise ValueError("mesh_name not recognized")

        # Make discr and connect it to firedrake
        factory = InterpolatoryQuadratureSimplexGroupFactory(fspace_degree)
        discr = Discretization(actx, mm_mesh, factory)

        fdrake_connection = build_connection_to_firedrake(discr)
        fdrake_fspace = fdrake_connection.firedrake_fspace()
        spatial_coord = SpatialCoordinate(fdrake_fspace.mesh())

        # get the group's nodes in a numpy array
        nodes = discr.nodes()
        group_nodes = np.array(
            [actx.to_numpy(dof_arr[0]) for dof_arr in nodes])

        for d in range(dim):
            meshmode_f = discr.zeros(actx)
            meshmode_f[0][:] = group_nodes[d, :, :]

            # connect to firedrake and evaluate expr in firedrake
            fdrake_f = Function(fdrake_fspace).interpolate(spatial_coord[d])

            # transport to firedrake and record error
            mm2fd_f = fdrake_connection.from_meshmode(meshmode_f)

            err = np.max(np.abs(fdrake_f.dat.data - mm2fd_f.dat.data))
            eoc_recorders[d].add_data_point(h, err)

    # assert that order is correct or error is "low enough"
    for d, eoc_rec in eoc_recorders.items():
        print("\nvector *x* -> *x[%s]*\n" % d, eoc_rec)
        assert (eoc_rec.order_estimate() >= fspace_degree
                or eoc_rec.max_error() < 2e-14)
예제 #6
0
def test_from_fd_transfer(ctx_factory, fspace_degree, fdrake_mesh_name,
                          fdrake_mesh_pars, dim, only_convert_bdy):
    """
    Make sure creating a function which projects onto
    one dimension then transports it is the same
    (up to resampling error) as projecting to one
    dimension on the transported mesh
    """
    # build estimate-of-convergence recorder
    from pytools.convergence import EOCRecorder
    # (fd -> mm ? True : False, dimension projecting onto)
    eoc_recorders = {(True, d): EOCRecorder() for d in range(dim)}
    if not only_convert_bdy:
        for d in range(dim):
            eoc_recorders[(False, d)] = EOCRecorder()

    # make a computing context
    cl_ctx = ctx_factory()
    queue = cl.CommandQueue(cl_ctx)
    actx = PyOpenCLArrayContext(queue)

    def get_fdrake_mesh_and_h_from_par(mesh_par):
        if fdrake_mesh_name == "UnitInterval":
            assert dim == 1
            n = mesh_par
            fdrake_mesh = UnitIntervalMesh(n)
            h = 1 / n
        elif fdrake_mesh_name == "UnitSquare":
            assert dim == 2
            n = mesh_par
            fdrake_mesh = UnitSquareMesh(n, n)
            h = 1 / n
        elif fdrake_mesh_name == "UnitCube":
            assert dim == 3
            n = mesh_par
            fdrake_mesh = UnitCubeMesh(n, n, n)
            h = 1 / n
        elif fdrake_mesh_name in ("blob2d-order1", "blob2d-order4"):
            assert dim == 2
            if fdrake_mesh_name == "blob2d-order1":
                from firedrake import Mesh
                fdrake_mesh = Mesh(f"{fdrake_mesh_name}-h{mesh_par}.msh",
                                   dim=dim)
            else:
                from meshmode.mesh.io import read_gmsh
                from meshmode.interop.firedrake import export_mesh_to_firedrake
                mm_mesh = read_gmsh(f"{fdrake_mesh_name}-h{mesh_par}.msh",
                                    force_ambient_dim=dim)
                fdrake_mesh, _, _ = export_mesh_to_firedrake(mm_mesh)
            h = float(mesh_par)
        elif fdrake_mesh_name == "warp":
            from meshmode.mesh.generation import generate_warped_rect_mesh
            from meshmode.interop.firedrake import export_mesh_to_firedrake
            mm_mesh = generate_warped_rect_mesh(dim,
                                                order=4,
                                                nelements_side=mesh_par)
            fdrake_mesh, _, _ = export_mesh_to_firedrake(mm_mesh)
            h = 1 / mesh_par
        else:
            raise ValueError("fdrake_mesh_name not recognized")

        return (fdrake_mesh, h)

    # Record error for each refinement of each mesh
    for mesh_par in fdrake_mesh_pars:
        fdrake_mesh, h = get_fdrake_mesh_and_h_from_par(mesh_par)
        # make function space and build connection
        fdrake_fspace = FunctionSpace(fdrake_mesh, "DG", fspace_degree)
        if only_convert_bdy:
            fdrake_connection = \
                build_connection_from_firedrake(actx,
                                                fdrake_fspace,
                                                restrict_to_boundary="on_boundary")
        else:
            fdrake_connection = build_connection_from_firedrake(
                actx, fdrake_fspace)
        # get this for making functions in firedrake
        spatial_coord = SpatialCoordinate(fdrake_mesh)

        # get nodes in handier format for making meshmode functions
        discr = fdrake_connection.discr
        # nodes is np array (ambient_dim,) of DOFArray (ngroups,)
        # of arrays (nelements, nunit_dofs), we want a single np array
        # of shape (ambient_dim, nelements, nunit_dofs)
        nodes = discr.nodes()
        group_nodes = np.array(
            [actx.to_numpy(dof_arr[0]) for dof_arr in nodes])

        # Now, for each coordinate d, test transferring the function
        # x -> sin(dth component of x)
        for d in range(dim):
            fdrake_f = Function(fdrake_fspace).interpolate(
                sin(spatial_coord[d]))
            # transport fdrake function and put in numpy
            fd2mm_f = fdrake_connection.from_firedrake(fdrake_f, actx=actx)
            fd2mm_f = actx.to_numpy(fd2mm_f[0])
            meshmode_f = np.sin(group_nodes[d, :, :])

            # record fd -> mm error
            err = np.max(np.abs(fd2mm_f - meshmode_f))
            eoc_recorders[(True, d)].add_data_point(h, err)

            if not only_convert_bdy:
                # now transport mm -> fd
                meshmode_f_dofarr = discr.zeros(actx)
                meshmode_f_dofarr[0][:] = meshmode_f
                mm2fd_f = fdrake_connection.from_meshmode(meshmode_f_dofarr)
                # record mm -> fd error
                err = np.max(np.abs(fdrake_f.dat.data - mm2fd_f.dat.data))
                eoc_recorders[(False, d)].add_data_point(h, err)

    # assert that order is correct or error is "low enough"
    for ((fd2mm, d), eoc_rec) in eoc_recorders.items():
        print(
            "\nfiredrake -> meshmode: %s\nvector *x* -> *sin(x[%s])*\n" %
            (fd2mm, d), eoc_rec)
        assert (eoc_rec.order_estimate() >= fspace_degree
                or eoc_rec.max_error() < 2e-14)
예제 #7
0
def main():
    # If can't import firedrake, do nothing
    #
    # filename MUST include "firedrake" (i.e. match *firedrake*.py) in order
    # to be run during CI
    try:
        import firedrake  # noqa : F401
    except ImportError:
        return 0

    # For this example, imagine we wish to solve the Laplace equation
    # on a meshmode mesh with some given Dirichlet boundary conditions,
    # and decide to use firedrake.
    #
    # To verify this is working, we use a solution to the wave equation
    # to get our boundary conditions

    # {{{ First we make a discretization in meshmode and get our bcs

    cl_ctx = cl.create_some_context()
    queue = cl.CommandQueue(cl_ctx)
    actx = PyOpenCLArrayContext(queue)

    nel_1d = 16
    from meshmode.mesh.generation import generate_regular_rect_mesh
    mesh = generate_regular_rect_mesh(
            a=(-0.5, -0.5),
            b=(0.5, 0.5),
            nelements_per_axis=(nel_1d, nel_1d))

    order = 3

    from meshmode.discretization import Discretization
    from meshmode.discretization.poly_element import \
        InterpolatoryQuadratureSimplexGroupFactory
    group_factory = InterpolatoryQuadratureSimplexGroupFactory(order=order)
    discr = Discretization(actx, mesh, group_factory)

    # Get our solution: we will use
    # Real(e^z) = Real(e^{x+iy})
    #           = e^x Real(e^{iy})
    #           = e^x cos(y)
    nodes = discr.nodes()
    for i in range(len(nodes)):
        nodes[i] = thaw(nodes[i], actx)
    # First index is dimension
    candidate_sol = actx.np.exp(nodes[0]) * actx.np.cos(nodes[1])

    # }}}

    # {{{ Now send candidate_sol into firedrake and use it for boundary conditions

    from meshmode.interop.firedrake import build_connection_to_firedrake
    fd_connection = build_connection_to_firedrake(discr, group_nr=0)
    # convert candidate_sol to firedrake
    fd_candidate_sol = fd_connection.from_meshmode(candidate_sol)
    # get the firedrake function space
    fd_fspace = fd_connection.firedrake_fspace()

    # set up dirichlet laplace problem in fd and solve
    from firedrake import (
        FunctionSpace, TrialFunction, TestFunction, Function, inner, grad, dx,
        Constant, project, DirichletBC, solve)

    # because it's easier to write down the variational problem,
    # we're going to project from our "DG" space
    # into a continuous one.
    cfd_fspace = FunctionSpace(fd_fspace.mesh(), "CG", order)
    u = TrialFunction(cfd_fspace)
    v = TestFunction(cfd_fspace)
    sol = Function(cfd_fspace)

    a = inner(grad(u), grad(v)) * dx
    rhs = Constant(0.0) * v * dx
    bc_value = project(fd_candidate_sol, cfd_fspace)
    bc = DirichletBC(cfd_fspace, bc_value, "on_boundary")
    params = {"ksp_monitor": None}
    solve(a == rhs, sol, bcs=[bc], solver_parameters=params)

    # project back into our "DG" space
    sol = project(sol, fd_fspace)

    # }}}

    # {{{ Take the solution from firedrake and compare it to candidate_sol

    true_sol = fd_connection.from_firedrake(sol, actx=actx)
    # pull back into numpy
    true_sol = actx.to_numpy(true_sol[0])
    candidate_sol = actx.to_numpy(candidate_sol[0])
    print("l^2 difference between candidate solution and firedrake solution=",
          np.linalg.norm(true_sol - candidate_sol))