예제 #1
0
def main(write_output=True):
    cl_ctx = cl.create_some_context()
    queue = cl.CommandQueue(cl_ctx)
    actx = PyOpenCLArrayContext(
        queue,
        allocator=cl_tools.MemoryPool(cl_tools.ImmediateAllocator(queue)),
        force_device_scalars=True,
    )

    from meshmode.mesh import BTAG_ALL
    from meshmode.mesh.generation import generate_warped_rect_mesh
    mesh = generate_warped_rect_mesh(dim=2, order=4, nelements_side=6)

    dcoll = DiscretizationCollection(actx, mesh, order=4)

    nodes = thaw(dcoll.nodes(), actx)
    bdry_nodes = thaw(dcoll.nodes(dd=BTAG_ALL), actx)
    bdry_normals = thaw(dcoll.normal(dd=BTAG_ALL), actx)

    if write_output:
        vis = shortcuts.make_visualizer(dcoll)
        vis.write_vtk_file("geo.vtu", [("nodes", nodes)])

        bvis = shortcuts.make_boundary_visualizer(dcoll)
        bvis.write_vtk_file("bgeo.vtu", [("bdry normals", bdry_normals),
                                         ("bdry nodes", bdry_nodes)])
예제 #2
0
def test_empty_boundary(actx_factory):
    # https://github.com/inducer/grudge/issues/54

    from meshmode.mesh import BTAG_NONE

    actx = actx_factory()

    dim = 2
    mesh = mgen.generate_regular_rect_mesh(a=(-0.5, ) * dim,
                                           b=(0.5, ) * dim,
                                           nelements_per_axis=(8, ) * dim,
                                           order=4)
    dcoll = DiscretizationCollection(actx, mesh, order=4)
    normal = dcoll.normal(BTAG_NONE)
    from meshmode.dof_array import DOFArray
    for component in normal:
        assert isinstance(component, DOFArray)
        assert len(component) == len(dcoll.discr_from_dd(BTAG_NONE).groups)
예제 #3
0
def test_2d_gauss_theorem(actx_factory):
    """Verify Gauss's theorem explicitly on a mesh"""

    pytest.importorskip("meshpy")

    from meshpy.geometry import make_circle, GeometryBuilder
    from meshpy.triangle import MeshInfo, build

    geob = GeometryBuilder()
    geob.add_geometry(*make_circle(1))
    mesh_info = MeshInfo()
    geob.set(mesh_info)

    mesh_info = build(mesh_info)

    from meshmode.mesh.io import from_meshpy
    from meshmode.mesh import BTAG_ALL

    mesh = from_meshpy(mesh_info, order=1)

    actx = actx_factory()

    dcoll = DiscretizationCollection(actx, mesh, order=2)
    volm_disc = dcoll.discr_from_dd(dof_desc.DD_VOLUME)
    x_volm = thaw(volm_disc.nodes(), actx)

    def f(x):
        return flat_obj_array(
            actx.np.sin(3 * x[0]) + actx.np.cos(3 * x[1]),
            actx.np.sin(2 * x[0]) + actx.np.cos(x[1]))

    f_volm = f(x_volm)
    int_1 = op.integral(dcoll, "vol", op.local_div(dcoll, f_volm))

    prj_f = op.project(dcoll, "vol", BTAG_ALL, f_volm)
    normal = thaw(dcoll.normal(BTAG_ALL), actx)
    int_2 = op.integral(dcoll, BTAG_ALL, prj_f.dot(normal))

    assert abs(int_1 - int_2) < 1e-13
예제 #4
0
def test_surface_divergence_theorem(actx_factory, mesh_name, visualize=False):
    r"""Check the surface divergence theorem.

        .. math::

            \int_Sigma \phi \nabla_i f_i =
            \int_\Sigma \nabla_i \phi f_i +
            \int_\Sigma \kappa \phi f_i n_i +
            \int_{\partial \Sigma} \phi f_i m_i

        where :math:`n_i` is the surface normal and :class:`m_i` is the
        face normal (which should be orthogonal to both the surface normal
        and the face tangent).
    """
    actx = actx_factory()

    # {{{ cases

    if mesh_name == "2-1-ellipse":
        from mesh_data import EllipseMeshBuilder
        builder = EllipseMeshBuilder(radius=3.1, aspect_ratio=2.0)
    elif mesh_name == "spheroid":
        from mesh_data import SpheroidMeshBuilder
        builder = SpheroidMeshBuilder()
    elif mesh_name == "circle":
        from mesh_data import EllipseMeshBuilder
        builder = EllipseMeshBuilder(radius=1.0, aspect_ratio=1.0)
    elif mesh_name == "starfish":
        from mesh_data import StarfishMeshBuilder
        builder = StarfishMeshBuilder()
    elif mesh_name == "sphere":
        from mesh_data import SphereMeshBuilder
        builder = SphereMeshBuilder(radius=1.0, mesh_order=16)
    else:
        raise ValueError("unknown mesh name: %s" % mesh_name)

    # }}}

    # {{{ convergence

    def f(x):
        return flat_obj_array(
            actx.np.sin(3 * x[1]) + actx.np.cos(3 * x[0]) + 1.0,
            actx.np.sin(2 * x[0]) + actx.np.cos(x[1]),
            3.0 * actx.np.cos(x[0] / 2) + actx.np.cos(x[1]),
        )[:ambient_dim]

    from pytools.convergence import EOCRecorder
    eoc_global = EOCRecorder()
    eoc_local = EOCRecorder()

    theta = np.pi / 3.33
    ambient_dim = builder.ambient_dim
    if ambient_dim == 2:
        mesh_rotation = np.array([
            [np.cos(theta), -np.sin(theta)],
            [np.sin(theta), np.cos(theta)],
        ])
    else:
        mesh_rotation = np.array([
            [1.0, 0.0, 0.0],
            [0.0, np.cos(theta), -np.sin(theta)],
            [0.0, np.sin(theta), np.cos(theta)],
        ])

    mesh_offset = np.array([0.33, -0.21, 0.0])[:ambient_dim]

    for i, resolution in enumerate(builder.resolutions):
        from meshmode.mesh.processing import affine_map
        from meshmode.discretization.connection import FACE_RESTR_ALL

        mesh = builder.get_mesh(resolution, builder.mesh_order)
        mesh = affine_map(mesh, A=mesh_rotation, b=mesh_offset)

        from meshmode.discretization.poly_element import \
                QuadratureSimplexGroupFactory

        qtag = dof_desc.DISCR_TAG_QUAD
        dcoll = DiscretizationCollection(actx,
                                         mesh,
                                         order=builder.order,
                                         discr_tag_to_group_factory={
                                             qtag:
                                             QuadratureSimplexGroupFactory(
                                                 2 * builder.order)
                                         })

        volume = dcoll.discr_from_dd(dof_desc.DD_VOLUME)
        logger.info("ndofs:     %d", volume.ndofs)
        logger.info("nelements: %d", volume.mesh.nelements)

        dd = dof_desc.DD_VOLUME
        dq = dd.with_discr_tag(qtag)
        df = dof_desc.as_dofdesc(FACE_RESTR_ALL)
        ambient_dim = dcoll.ambient_dim

        # variables
        f_num = f(thaw(dcoll.nodes(dd=dd), actx))
        f_quad_num = f(thaw(dcoll.nodes(dd=dq), actx))

        from grudge.geometry import normal, summed_curvature

        kappa = summed_curvature(actx, dcoll, dd=dq)
        normal = normal(actx, dcoll, dd=dq)
        face_normal = thaw(dcoll.normal(df), actx)
        face_f = op.project(dcoll, dd, df, f_num)

        # operators
        stiff = op.mass(
            dcoll,
            sum(
                op.local_d_dx(dcoll, i, f_num_i)
                for i, f_num_i in enumerate(f_num)))
        stiff_t = sum(
            op.weak_local_d_dx(dcoll, i, f_num_i)
            for i, f_num_i in enumerate(f_num))
        kterm = op.mass(dcoll, dq, kappa * f_quad_num.dot(normal))
        flux = op.face_mass(dcoll, face_f.dot(face_normal))

        # sum everything up
        op_global = op.nodal_sum(dcoll, dd, stiff - (stiff_t + kterm))
        op_local = op.elementwise_sum(dcoll, dd,
                                      stiff - (stiff_t + kterm + flux))

        err_global = abs(op_global)
        err_local = op.norm(dcoll, op_local, np.inf)
        logger.info("errors: global %.5e local %.5e", err_global, err_local)

        # compute max element size
        from grudge.dt_utils import h_max_from_volume

        h_max = h_max_from_volume(dcoll)

        eoc_global.add_data_point(h_max, actx.to_numpy(err_global))
        eoc_local.add_data_point(h_max, err_local)

        if visualize:
            from grudge.shortcuts import make_visualizer
            vis = make_visualizer(dcoll)

            filename = f"surface_divergence_theorem_{mesh_name}_{i:04d}.vtu"
            vis.write_vtk_file(filename, [("r", actx.np.log10(op_local))],
                               overwrite=True)

    # }}}

    order = min(builder.order, builder.mesh_order) - 0.5
    logger.info("\n%s", str(eoc_global))
    logger.info("\n%s", str(eoc_local))

    assert eoc_global.max_error() < 1.0e-12 \
            or eoc_global.order_estimate() > order - 0.5

    assert eoc_local.max_error() < 1.0e-12 \
            or eoc_local.order_estimate() > order - 0.5
예제 #5
0
def test_face_normal_surface(actx_factory, mesh_name):
    """Check that face normals are orthogonal to the surface normal"""
    actx = actx_factory()

    # {{{ geometry

    if mesh_name == "2-1-ellipse":
        from mesh_data import EllipseMeshBuilder
        builder = EllipseMeshBuilder(radius=3.1, aspect_ratio=2.0)
    elif mesh_name == "spheroid":
        from mesh_data import SpheroidMeshBuilder
        builder = SpheroidMeshBuilder()
    else:
        raise ValueError("unknown mesh name: %s" % mesh_name)

    mesh = builder.get_mesh(builder.resolutions[0], builder.mesh_order)
    dcoll = DiscretizationCollection(actx, mesh, order=builder.order)

    volume_discr = dcoll.discr_from_dd(dof_desc.DD_VOLUME)
    logger.info("ndofs:    %d", volume_discr.ndofs)
    logger.info("nelements: %d", volume_discr.mesh.nelements)

    # }}}

    # {{{ Compute surface and face normals
    from meshmode.discretization.connection import FACE_RESTR_INTERIOR
    from grudge.geometry import normal

    dv = dof_desc.DD_VOLUME
    df = dof_desc.as_dofdesc(FACE_RESTR_INTERIOR)

    ambient_dim = mesh.ambient_dim

    surf_normal = op.project(dcoll, dv, df, normal(actx, dcoll, dd=dv))
    surf_normal = surf_normal / actx.np.sqrt(sum(surf_normal**2))

    face_normal_i = thaw(dcoll.normal(df), actx)
    face_normal_e = dcoll.opposite_face_connection()(face_normal_i)

    if mesh.ambient_dim == 3:
        from grudge.geometry import pseudoscalar, area_element
        # NOTE: there's only one face tangent in 3d
        face_tangent = (pseudoscalar(actx, dcoll, dd=df) /
                        area_element(actx, dcoll, dd=df)).as_vector(
                            dtype=object)

    # }}}

    # {{{ checks

    def _eval_error(x):
        return op.norm(dcoll, x, np.inf, dd=df)

    rtol = 1.0e-14

    # check interpolated surface normal is orthogonal to face normal
    error = _eval_error(surf_normal.dot(face_normal_i))
    logger.info("error[n_dot_i]:    %.5e", error)
    assert error < rtol

    # check angle between two neighboring elements
    error = _eval_error(face_normal_i.dot(face_normal_e) + 1.0)
    logger.info("error[i_dot_e]:    %.5e", error)
    assert error > rtol

    # check orthogonality with face tangent
    if ambient_dim == 3:
        error = _eval_error(face_tangent.dot(face_normal_i))
        logger.info("error[t_dot_i]:  %.5e", error)
        assert error < 5 * rtol
예제 #6
0
def main(ctx_factory, dim=2, order=4, use_quad=False, visualize=False):
    cl_ctx = ctx_factory()
    queue = cl.CommandQueue(cl_ctx)
    actx = PyOpenCLArrayContext(
        queue,
        allocator=cl_tools.MemoryPool(cl_tools.ImmediateAllocator(queue)),
        force_device_scalars=True,
    )

    # {{{ parameters

    # sphere radius
    radius = 1.0
    # sphere resolution
    resolution = 64 if dim == 2 else 1

    # final time
    final_time = np.pi

    # flux
    flux_type = "lf"

    # }}}

    # {{{ discretization

    if dim == 2:
        from meshmode.mesh.generation import make_curve_mesh, ellipse
        mesh = make_curve_mesh(
                lambda t: radius * ellipse(1.0, t),
                np.linspace(0.0, 1.0, resolution + 1),
                order)
    elif dim == 3:
        from meshmode.mesh.generation import generate_icosphere
        mesh = generate_icosphere(radius, order=4 * order,
                uniform_refinement_rounds=resolution)
    else:
        raise ValueError("unsupported dimension")

    discr_tag_to_group_factory = {}
    if use_quad:
        qtag = dof_desc.DISCR_TAG_QUAD
    else:
        qtag = None

    from meshmode.discretization.poly_element import \
            default_simplex_group_factory, \
            QuadratureSimplexGroupFactory

    discr_tag_to_group_factory[dof_desc.DISCR_TAG_BASE] = \
        default_simplex_group_factory(base_dim=dim-1, order=order)

    if use_quad:
        discr_tag_to_group_factory[qtag] = \
            QuadratureSimplexGroupFactory(order=4*order)

    from grudge import DiscretizationCollection

    dcoll = DiscretizationCollection(
        actx, mesh,
        discr_tag_to_group_factory=discr_tag_to_group_factory
    )

    volume_discr = dcoll.discr_from_dd(dof_desc.DD_VOLUME)
    logger.info("ndofs:     %d", volume_discr.ndofs)
    logger.info("nelements: %d", volume_discr.mesh.nelements)

    # }}}

    # {{{ Surface advection operator

    # velocity field
    x = thaw(dcoll.nodes(), actx)
    c = make_obj_array([-x[1], x[0], 0.0])[:dim]

    def f_initial_condition(x):
        return x[0]

    from grudge.models.advection import SurfaceAdvectionOperator
    adv_operator = SurfaceAdvectionOperator(
        dcoll,
        c,
        flux_type=flux_type,
        quad_tag=qtag
    )

    u0 = f_initial_condition(x)

    def rhs(t, u):
        return adv_operator.operator(t, u)

    # check velocity is tangential
    from grudge.geometry import normal

    surf_normal = normal(actx, dcoll, dd=dof_desc.DD_VOLUME)

    error = op.norm(dcoll, c.dot(surf_normal), 2)
    logger.info("u_dot_n:   %.5e", error)

    # }}}

    # {{{ time stepping

    # FIXME: dt estimate is not necessarily valid for surfaces
    dt = actx.to_numpy(
        0.45 * adv_operator.estimate_rk4_timestep(actx, dcoll, fields=u0))
    nsteps = int(final_time // dt) + 1

    logger.info("dt:        %.5e", dt)
    logger.info("nsteps:    %d", nsteps)

    from grudge.shortcuts import set_up_rk4
    dt_stepper = set_up_rk4("u", dt, u0, rhs)
    plot = Plotter(actx, dcoll, order, visualize=visualize)

    norm_u = actx.to_numpy(op.norm(dcoll, u0, 2))

    step = 0

    event = dt_stepper.StateComputed(0.0, 0, 0, u0)
    plot(event, "fld-surface-%04d" % 0)

    if visualize and dim == 3:
        from grudge.shortcuts import make_visualizer
        vis = make_visualizer(dcoll)
        vis.write_vtk_file(
            "fld-surface-velocity.vtu",
            [
                ("u", c),
                ("n", surf_normal)
            ],
            overwrite=True
        )

        df = dof_desc.DOFDesc(FACE_RESTR_INTERIOR)
        face_discr = dcoll.discr_from_dd(df)
        face_normal = thaw(dcoll.normal(dd=df), actx)

        from meshmode.discretization.visualization import make_visualizer
        vis = make_visualizer(actx, face_discr)
        vis.write_vtk_file("fld-surface-face-normals.vtu", [
            ("n", face_normal)
            ], overwrite=True)

    for event in dt_stepper.run(t_end=final_time, max_steps=nsteps + 1):
        if not isinstance(event, dt_stepper.StateComputed):
            continue

        step += 1
        if step % 10 == 0:
            norm_u = actx.to_numpy(op.norm(dcoll, event.state_component, 2))
            plot(event, "fld-surface-%04d" % step)

        logger.info("[%04d] t = %.5f |u| = %.5e", step, event.t, norm_u)

        # NOTE: These are here to ensure the solution is bounded for the
        # time interval specified
        assert norm_u < 3