Esempio n. 1
0
def tangential_onb(ambient_dim, dim=None, where=None):
    """Return a matrix of shape ``(ambient_dim, dim)`` with orthogonal columns
    spanning the tangential space of the surface of *where*.
    """

    if dim is None:
        dim = ambient_dim - 1

    pd_mat = parametrization_derivative_matrix(ambient_dim, dim, where)

    # {{{ Gram-Schmidt

    orth_pd_mat = np.zeros_like(pd_mat)
    for k in range(pd_mat.shape[1]):
        avec = pd_mat[:, k]
        q = avec
        for j in range(k):
            q = q - np.dot(avec, orth_pd_mat[:, j])*orth_pd_mat[:, j]
        q = cse(q, "q%d" % k)

        orth_pd_mat[:, k] = cse(q/sqrt(np.sum(q**2)), "orth_pd_vec%d_" % k)

    # }}}

    return orth_pd_mat
Esempio n. 2
0
def mv_normal(dd, ambient_dim, dim=None):
    """Exterior unit normal as a :class:`~pymbolic.geometric_algebra.MultiVector`."""
    dd = as_dofdesc(dd)
    if not dd.is_trace():
        raise ValueError("may only request normals on boundaries")

    if dim is None:
        dim = ambient_dim - 1

    if dim == ambient_dim - 1:
        return surface_normal(ambient_dim, dim=dim, dd=dd)

    # NOTE: In the case of (d - 2)-dimensional curves, we don't really have
    # enough information on the face to decide what an "exterior face normal"
    # is (e.g the "normal" to a 1D curve in 3D space is actually a
    # "normal plane")
    #
    # The trick done here is that we take the surface normal, move it to the
    # face and then take a cross product with the face tangent to get the
    # correct exterior face normal vector.
    assert dim == ambient_dim - 2

    from grudge.symbolic.operators import project
    volm_normal = (
            surface_normal(ambient_dim, dim=dim + 1, dd=DD_VOLUME)
            .map(project(DD_VOLUME, dd)))
    pder = pseudoscalar(ambient_dim, dim, dd=dd)

    mv = cse(-(volm_normal ^ pder) << volm_normal.I.inv(),
            "face_normal",
            cse_scope.DISCRETIZATION)

    return cse(mv / sqrt(mv.norm_squared()),
            "unit_face_normal",
            cse_scope.DISCRETIZATION)
Esempio n. 3
0
def int_tpair(expression, qtag=None):
    i = _sym().interp("vol", "int_faces")(expression)
    e = cse(_sym().OppositeInteriorFaceSwap()(i))

    if qtag is not None and qtag != _sym().QTAG_NONE:
        q_dd = _sym().DOFDesc("int_faces", qtag)
        i = cse(_sym().interp("int_faces", q_dd)(i))
        e = cse(_sym().interp("int_faces", q_dd)(e))
    else:
        q_dd = "int_faces"

    return TracePair(q_dd, i, e)
Esempio n. 4
0
def int_tpair(expression, qtag=None):
    from grudge.symbolic.operators import project, OppositeInteriorFaceSwap

    i = project("vol", "int_faces")(expression)
    e = cse(OppositeInteriorFaceSwap()(i))

    if qtag is not None and qtag != QTAG_NONE:
        q_dd = DOFDesc("int_faces", qtag)
        i = cse(project("int_faces", q_dd)(i))
        e = cse(project("int_faces", q_dd)(e))
    else:
        q_dd = "int_faces"

    return TracePair(q_dd, i, e)
Esempio n. 5
0
def forward_metric_derivative(xyz_axis, rst_axis, dd=None):
    r"""
    Pointwise metric derivatives representing

    .. math::

        \frac{\partial x_{\mathrm{xyz\_axis}} }{\partial r_{\mathrm{rst\_axis}} }
    """
    if dd is None:
        dd = DD_VOLUME

    inner_dd = dd.with_qtag(QTAG_NONE)

    from grudge.symbolic.operators import RefDiffOperator
    diff_op = RefDiffOperator(rst_axis, inner_dd)

    result = diff_op(NodeCoordinateComponent(xyz_axis, inner_dd))

    if dd.uses_quadrature():
        from grudge.symbolic.operators import project
        result = project(inner_dd, dd)(result)

    return cse(result,
               prefix="dx%d_dr%d" % (xyz_axis, rst_axis),
               scope=cse_scope.DISCRETIZATION)
Esempio n. 6
0
def S_surface_laplacian_S(u, dim, invertibility_scale=0, qbx_fix_scale=0):
    """
    :arg u: The field to which the surface Laplacian is applied.
    """

    # This converges, but appears to work quite poorly compared to the above.

    tgrad_Su = cse(
            project_to_tangential(grad_S(0, u, dim)),
            "tgrad_su_from_surflap")

    return (
            - IntGdSource(0, Ones(), ds_direction=real(tgrad_Su))
            - 1j*IntGdSource(0, Ones(), ds_direction=imag(tgrad_Su))
            - invertibility_scale * S(0, Ones()*mean(S(0, u)))
            - qbx_fix_scale * (
                u
                # D+ - D- = identity (but QBX will compute the
                # 'compact part of the identity' -- call that I*)
                - (
                    D(0, u, qbx_forced_limit=+1)
                    - D(0, u, qbx_forced_limit=-1))
                )
                # The above is I - I*, which means only the high-frequency
                # bits of the identity are left.
            )
Esempio n. 7
0
def mv_normal(dd, ambient_dim, dim=None):
    """Exterior unit normal as a :class:`~pymbolic.geometric_algebra.MultiVector`."""

    dd = as_dofdesc(dd)

    if not dd.is_trace():
        raise ValueError("may only request normals on boundaries")

    if dim is None:
        dim = ambient_dim - 1

    # NOTE: Don't be tempted to add a sign here. As it is, it produces
    # exterior normals for positively oriented curves.

    pder = pseudoscalar(ambient_dim, dim, dd=dd) \
            / area_element(ambient_dim, dim, dd=dd)

    # Dorst Section 3.7.2
    mv = pder << pder.I.inv()

    if dim == 0:
        # NOTE: when the mesh is 0D, we do not have a clear notion of
        # `exterior normal`, so we just take the tangent of the 1D element,
        # project it to the element faces and make it signed
        tangent = parametrization_derivative(ambient_dim, dim=1, dd=DD_VOLUME)
        tangent = tangent / sqrt(tangent.norm_squared())

        from grudge.symbolic.operators import project
        project = project(DD_VOLUME, dd)
        mv = MultiVector(
            np.array(
                [mv.as_scalar() * project(t) for t in tangent.as_vector()]))

    return cse(mv, "normal", cse_scope.DISCRETIZATION)
Esempio n. 8
0
def make_sym_surface_mv(name, ambient_dim, dim, where=None):
    par_grad = parametrization_derivative_matrix(ambient_dim, dim, where)

    return sum(
            var("%s%d" % (name, i))
            * cse(MultiVector(vec), "tangent%d" % i, cse_scope.DISCRETIZATION)
            for i, vec in enumerate(par_grad.T))
Esempio n. 9
0
def shape_operator(ambient_dim, dim=None, dd=None):
    if dim is None:
        dim = ambient_dim - 1

    inv_form1 = inverse_first_fundamental_form(ambient_dim, dim=dim, dd=dd)
    form2 = second_fundamental_form(ambient_dim, dim=dim, dd=dd)

    return cse(-form2.dot(inv_form1), "shape_operator", cse_scope.DISCRETIZATION)
Esempio n. 10
0
def pseudoscalar(ambient_dim, dim=None, dd=None):
    if dim is None:
        dim = ambient_dim

    return cse(
        parametrization_derivative(ambient_dim, dim,
                                   dd=dd).project_max_grade(), "pseudoscalar",
        cse_scope.DISCRETIZATION)
Esempio n. 11
0
def _simplex_mapping_max_stretch_factor(ambient_dim, dim=None, where=None,
        with_elementwise_max=True):
    """Return the largest factor by which the reference-to-global
    mapping stretches the bi-unit (i.e. :math:`[-1,1]`) reference
    element along any axis.

    Returns a DOF vector that is elementwise constant.
    """

    if dim is None:
        dim = ambient_dim - 1

    # The 'technique' here is ad-hoc, but I'm fairly confident it's better than
    # what we had. The idea is that singular values of the mapping Jacobian
    # yield "stretch factors" of the mapping Why? because it maps a right
    # singular vector $`v_1`$ (of unit length) to $`\sigma_1 u_1`$, where
    # $`u_1`$ is the corresponding left singular vector (also of unit length).
    # And so the biggest one tells us about the direction with the 'biggest'
    # stretching, where 'stretching' (*2 to remove bi-unit reference element)
    # reflects available quadrature resolution in that direction.

    equi_pder_mat = _equilateral_parametrization_derivative_matrix(
            ambient_dim, dim, where)

    # Compute eigenvalues of J^T to compute SVD.
    equi_pder_mat_jtj = cse(
            np.dot(equi_pder_mat.T, equi_pder_mat),
            "pd_mat_jtj")

    stretch_factors = [
            cse(sqrt(s), "mapping_singval_%d" % i)
            for i, s in enumerate(
                _small_mat_eigenvalues(
                    # Multiply by 4 to compensate for equilateral reference
                    # elements of side length 2. (J^T J contains two factors of
                    # two.)
                    4 * equi_pder_mat_jtj))]

    from pymbolic.primitives import Max
    result = Max(tuple(stretch_factors))

    if with_elementwise_max:
        result = ElementwiseMax(result, where=where)

    return cse(result, "mapping_max_stretch", cse_scope.DISCRETIZATION)
Esempio n. 12
0
def normal(where=None):
    """Exterior unit normals."""

    # Don't be tempted to add a sign here. As it is, it produces
    # exterior normals for positively oriented curves.

    pder = pseudoscalar(where) / area_element(where)
    return cse(pder.a.I | pder, "normal",
            cse_scope.DISCRETIZATION)
Esempio n. 13
0
def int_tpair(expression, qtag=None, from_dd=None):
    from grudge.symbolic.operators import project, OppositeInteriorFaceSwap

    if from_dd is None:
        from_dd = DD_VOLUME
    assert not from_dd.uses_quadrature()

    trace_dd = DOFDesc(FACE_RESTR_INTERIOR, qtag)
    if from_dd.domain_tag == trace_dd.domain_tag:
        i = expression
    else:
        i = project(from_dd, trace_dd.with_qtag(None))(expression)
    e = cse(OppositeInteriorFaceSwap()(i))

    if trace_dd.uses_quadrature():
        i = cse(project(trace_dd.with_qtag(None), trace_dd)(i))
        e = cse(project(trace_dd.with_qtag(None), trace_dd)(e))

    return TracePair(trace_dd, interior=i, exterior=e)
Esempio n. 14
0
def parametrization_derivative_matrix(ambient_dim, dim, where=None):
    """Return a :class:`np.array` representing the derivative of the
    reference-to-global parametrization.
    """

    return cse(
            reference_jacobian(
                [NodeCoordinateComponent(i, where) for i in range(ambient_dim)],
                ambient_dim, dim, where=where),
            "pd_matrix", cse_scope.DISCRETIZATION)
Esempio n. 15
0
def forward_metric_nth_derivative(xyz_axis, ref_axes, dd=None):
    r"""
    Pointwise metric derivatives representing repeated derivatives

    .. math::

        \frac{\partial^n x_{\mathrm{xyz\_axis}} }{\partial r_{\mathrm{ref\_axes}}}

    where *ref_axes* is a multi-index description.

    :arg ref_axes: a :class:`tuple` of tuples indicating indices of
        coordinate axes of the reference element to the number of derivatives
        which will be taken. For example, the value ``((0, 2), (1, 1))``
        indicates taking the second derivative with respect to the first
        axis and the first derivative with respect to the second
        axis. Each axis must occur only once and the tuple must be sorted
        by the axis index.

        May also be a singile integer *i*, which is viewed as equivalent
        to ``((i, 1),)``.
    """

    if isinstance(ref_axes, int):
        ref_axes = ((ref_axes, 1),)

    if not isinstance(ref_axes, tuple):
        raise ValueError("ref_axes must be a tuple")

    if tuple(sorted(ref_axes)) != ref_axes:
        raise ValueError("ref_axes must be sorted")

    if len(dict(ref_axes)) != len(ref_axes):
        raise ValueError("ref_axes must not contain an axis more than once")

    if dd is None:
        dd = DD_VOLUME
    inner_dd = dd.with_qtag(QTAG_NONE)

    from pytools import flatten
    flat_ref_axes = flatten([rst_axis] * n for rst_axis, n in ref_axes)

    from grudge.symbolic.operators import RefDiffOperator
    result = NodeCoordinateComponent(xyz_axis, inner_dd)
    for rst_axis in flat_ref_axes:
        result = RefDiffOperator(rst_axis, inner_dd)(result)

    if dd.uses_quadrature():
        from grudge.symbolic.operators import project
        result = project(inner_dd, dd)(result)

    prefix = "dx%d_%s" % (
            xyz_axis,
            "_".join("%sr%d" % ("d" * n, rst_axis) for rst_axis, n in ref_axes))

    return cse(result, prefix, cse_scope.DISCRETIZATION)
    def flux(self, w, k):
        F2 = self.F2(w) #get F'' from state vector w
        q = self.q(w)
        P = self.P(q)
        v = self.v(q)
        v_null = Field('state_null')

        dim = self.dimensions

        # One entry for each flux direction
        if dim == 1:
            raise NotImplementedError
            #return [cse(join_fields(
            #           v_null,
            #            (P[0]+F2[0])/k[0],  # flux rho_v
            #            (v[0]+F2[1])/k[0]   # flux F
            #            ), "x_flux")]
        elif dim == 2:
            return [cse(join_fields(
                        v_null, (P[0]+F2[0])/k[0],(P[3]+F2[1])/k[0],      # flux rho_v
                        (v[0]+F2[2])/k[0],v_null,v_null,(v[1]+F2[3])/k[0] # flux F
                        ), "x_flux"),
                    cse(join_fields(
                        v_null, (P[2]+F2[4])/k[1],(P[1]+F2[5])/k[1],      # flux rho_v
                        v_null,(v[1]+F2[6])/k[1],(v[0]+F2[7])/k[1],v_null # flux F
                        ), "y_flux")]
        elif dim == 3:
            raise NotImplementedError
            #return [cse(join_fields(
            #            v_null, (P[0]+F2[0])/k[0],(P[8]+F2[1])/k[0],(P[7]+F2[2])/k[0],                                     # flux rho_v
            #            (v[0]+F2[3])/k[0],v_null,v_null,v_null,v_null,v_null,v_null,(v[2]+F2[4])/k[0],(v[1]+F2[5])/k[0]    # flux F
            #            ), "x_flux"),
            #        cse(join_fields(
            #            v_null, (P[5]+F2[6])/k[1],(P[1]+F2[7])/k[1],(P[6]+F2[8])/k[1],                                     # flux rho_v
            #            v_null,(v[1]+F2[9])/k[1],v_null,v_null,v_null,(v[0]+F2[10])/k[1],(v[2]+F2[11])/k[1],v_null,v_null  # flux F
            #            ), "y_flux"),
            #        cse(join_fields(
            #            v_null, (P[4]+F2[12])/k[2],(P[3]+F2[13])/k[2],(P[2]+F2[14])/k[2],                                  # flux rho_v
            #            v_null,v_null,(v[2]+F2[15])/k[2],(v[1]+F2[16])/k[2],(v[0]+F2[17])/k[2],v_null,v_null,v_null,v_null # flux F
            #            ), "z_flux")]
        else:
            raise ValueError("Invalid dimension")
Esempio n. 17
0
def pseudoscalar(ambient_dim, dim=None, where=None):
    """
    Same as the outer product of all parametrization derivative columns.
    """
    if dim is None:
        dim = ambient_dim - 1

    return cse(
            parametrization_derivative(ambient_dim, dim, where)
            .project_max_grade(),
            "pseudoscalar", cse_scope.DISCRETIZATION)
Esempio n. 18
0
def bv_tpair(dd, interior, exterior):
    """
    :arg interior: an expression that lives in the volume
        and will be restricted to the boundary identified by
        *tag* before use.
    :arg exterior: an expression that already lives on the boundary
        representing the exterior value to be used
        for the flux.
    """
    from grudge.symbolic.operators import project
    interior = cse(project("vol", dd)(interior))
    return TracePair(dd, interior, exterior)
Esempio n. 19
0
def surface_laplacian_S_squared(u, invertibility_scale=0):
    """
    :arg u: The field to which the surface Laplacian is applied.
    """
    # http://wiki.tiker.net/HellsKitchen/SurfaceLaplacian

    Su = cse(S(0, u), "su_from_surflap")

    return (
            - 2*mean_curvature()*Sp(0, Su)
            - ((Spp(0, Su)+Dp(0, Su))-(-1/4*u+Sp(0, Sp(0, u))))
            - invertibility_scale * mean(S(0, Su))*Ones())
Esempio n. 20
0
def first_fundamental_form(ambient_dim, dim=None, where=None):
    if dim is None:
        dim = ambient_dim - 1

    if ambient_dim != 3 and dim != 2:
        raise NotImplementedError("only available for surfaces in 3D")

    pd_mat = parametrization_derivative_matrix(ambient_dim, dim, where)

    return cse(
            np.dot(pd_mat.T, pd_mat),
            "fundform1")
Esempio n. 21
0
def inverse_metric_derivative(rst_axis,
                              xyz_axis,
                              ambient_dim,
                              dim=None,
                              dd=None):
    if dim is None:
        dim = ambient_dim

    if dim != ambient_dim:
        raise ValueError("not clear what inverse_metric_derivative means if "
                         "the derivative matrix is not square")

    par_vecs = [
        forward_metric_derivative_mv(ambient_dim, rst, dd)
        for rst in range(dim)
    ]

    # Yay Cramer's rule! (o rly?)
    from functools import reduce, partial
    from operator import xor as outerprod_op
    outerprod = partial(reduce, outerprod_op)

    def outprod_with_unit(i, at):
        unit_vec = np.zeros(dim)
        unit_vec[i] = 1

        vecs = par_vecs[:]
        vecs[at] = MultiVector(unit_vec)

        return outerprod(vecs)

    volume_pseudoscalar_inv = cse(
        outerprod(
            forward_metric_derivative_mv(ambient_dim, rst_axis, dd=dd)
            for rst_axis in range(dim)).inv())

    return cse((outprod_with_unit(xyz_axis, rst_axis) *
                volume_pseudoscalar_inv).as_scalar(),
               prefix="dr%d_dx%d" % (rst_axis, xyz_axis),
               scope=cse_scope.DISCRETIZATION)
Esempio n. 22
0
def normal(ambient_dim, dim=None, where=None):
    """Exterior unit normals."""

    # Don't be tempted to add a sign here. As it is, it produces
    # exterior normals for positively oriented curves and surfaces.

    pder = (
            pseudoscalar(ambient_dim, dim, where)
            / area_element(ambient_dim, dim, where))
    return cse(
            # Dorst Section 3.7.2
            pder << pder.I.inv(),
            "normal",
            scope=cse_scope.DISCRETIZATION)
Esempio n. 23
0
def summed_curvature(ambient_dim, dim=None, dd=None):
    """Sum of the principal curvatures"""

    if dim is None:
        dim = ambient_dim - 1

    if ambient_dim == 1:
        return 0.0

    if ambient_dim == dim:
        return 0.0

    op = shape_operator(ambient_dim, dim=dim, dd=dd)
    return cse(np.trace(op), "summed_curvature", cse_scope.DISCRETIZATION)
Esempio n. 24
0
def second_fundamental_form(ambient_dim, dim=None, where=None):
    """Compute the second fundamental form of a surface. This is in reference
    to the reference-to-global mapping in use for each element.

    .. note::

        Some references assume that the second fundamental form is computed
        with respect to an orthonormal basis, which this is not.
    """
    if dim is None:
        dim = ambient_dim - 1

    if ambient_dim != 3 and dim != 2:
        raise NotImplementedError("only available for surfaces in 3D")

    r = nodes(ambient_dim, where=where).as_vector()

    # https://en.wikipedia.org/w/index.php?title=Second_fundamental_form&oldid=821047433#Classical_notation

    from functools import partial
    d = partial(NumReferenceDerivative, where=where)
    ruu = d(((0, 2),), r)
    ruv = d(((0, 1), (1, 1)), r)
    rvv = d(((1, 2),), r)

    nrml = normal(ambient_dim, dim, where).as_vector()

    ff2_l = cse(np.dot(ruu, nrml), "fundform2_L")
    ff2_m = cse(np.dot(ruv, nrml), "fundform2_M")
    ff2_n = cse(np.dot(rvv, nrml), "fundform2_N")

    result = np.zeros((2, 2), dtype=object)
    result[0, 0] = ff2_l
    result[0, 1] = result[1, 0] = ff2_m
    result[1, 1] = ff2_n

    return result
Esempio n. 25
0
def surface_normal(ambient_dim, dim=None, dd=None):
    dd = as_dofdesc(dd)
    if dim is None:
        dim = ambient_dim - 1

    # NOTE: Don't be tempted to add a sign here. As it is, it produces
    # exterior normals for positively oriented curves.

    pder = pseudoscalar(ambient_dim, dim, dd=dd) \
            / area_element(ambient_dim, dim, dd=dd)

    # Dorst Section 3.7.2
    return cse(pder << pder.I.inv(),
            "surface_normal",
            cse_scope.DISCRETIZATION)
Esempio n. 26
0
def _equilateral_parametrization_derivative_matrix(ambient_dim, dim=None,
        where=None):
    if dim is None:
        dim = ambient_dim - 1

    pder_mat = parametrization_derivative_matrix(ambient_dim, dim, where)

    # The above procedure works well only when the 'reference' end of the
    # mapping is in equilateral coordinates.
    from modepy.tools import EQUILATERAL_TO_UNIT_MAP
    equi_to_unit = EQUILATERAL_TO_UNIT_MAP[dim].a

    # This is the Jacobian of the (equilateral reference element) -> (global) map.
    return cse(
            np.dot(pder_mat, equi_to_unit),
            "equilateral_pder_mat")
Esempio n. 27
0
def mean_curvature(ambient_dim, dim=None, where=None):
    """(Numerical) mean curvature."""

    if dim is None:
        dim = ambient_dim - 1

    if not (dim == 1 and ambient_dim == 2):
        raise NotImplementedError(
                "only know how to calculate curvature for a curve in 2D")

    xp, yp = parametrization_derivative_matrix(ambient_dim, dim, where)

    xpp, ypp = cse(
            reference_jacobian([xp[0], yp[0]], ambient_dim, dim, where),
            "p2d_matrix", cse_scope.DISCRETIZATION)

    return (xp[0]*ypp[0] - yp[0]*xpp[0]) / (xp[0]**2 + yp[0]**2)**(3/2)
Esempio n. 28
0
def inverse_surface_metric_derivative(rst_axis, xyz_axis,
        ambient_dim, dim=None, dd=None):
    if dim is None:
        dim = ambient_dim - 1

    if ambient_dim == dim:
        return inverse_metric_derivative(rst_axis, xyz_axis,
                ambient_dim, dim=dim, dd=dd)
    else:
        inv_form1 = inverse_first_fundamental_form(ambient_dim, dim=dim, dd=dd)
        imd = sum(
                inv_form1[rst_axis, d]*forward_metric_derivative(xyz_axis, d, dd=dd)
                for d in range(dim))

        return cse(imd,
                prefix="ds%d_dx%d" % (rst_axis, xyz_axis),
                scope=cse_scope.DISCRETIZATION)
Esempio n. 29
0
def _max_curvature(ambient_dim, dim=None, where=None):
    # An attempt at a 'max curvature' criterion.

    if dim is None:
        dim = ambient_dim - 1

    if ambient_dim == 2:
        return abs(mean_curvature(ambient_dim, dim, where=where))
    elif ambient_dim == 3:
        shape_op = shape_operator(ambient_dim, dim, where=where)

        abs_principal_curvatures = [
                abs(x) for x in _small_mat_eigenvalues(shape_op)]
        from pymbolic.primitives import Max
        return cse(Max(tuple(abs_principal_curvatures)))
    else:
        raise NotImplementedError("curvature criterion not implemented in %d "
                "dimensions" % ambient_dim)
Esempio n. 30
0
def int_g_dsource(ambient_dim, dsource, kernel, density,
            qbx_forced_limit, source=None, target=None,
            kernel_arguments=None, **kwargs):
    r"""
    .. math::

        \int_\Gamma \operatorname{dsource} \dot \nabla_y
            \dot g(x-y) \sigma(y) dS_y

    where :math:`\sigma` is *density*, and
    *dsource*, a multivector.
    Note that the first product in the integrand
    is a geometric product.

    .. attribute:: dsource

        A :class:`pymbolic.geometric_algebra.MultiVector`.
    """

    if kernel_arguments is None:
        kernel_arguments = {}

    if isinstance(kernel_arguments, tuple):
        kernel_arguments = dict(kernel_arguments)

    kernel = _insert_source_derivative_into_kernel(kernel)

    from pytools.obj_array import make_obj_array
    nabla = MultiVector(make_obj_array(
        [NablaComponent(axis, None)
            for axis in range(ambient_dim)]))

    def add_dir_vec_to_kernel_args(coeff):
        result = kernel_arguments.copy()
        result[_DIR_VEC_NAME] = _get_dir_vec(coeff, ambient_dim)
        return result

    density = cse(density)
    return (dsource*nabla).map(
            lambda coeff: IntG(
                kernel,
                density, qbx_forced_limit, source, target,
                kernel_arguments=add_dir_vec_to_kernel_args(coeff),
                **kwargs))
Esempio n. 31
0
def inverse_first_fundamental_form(ambient_dim, dim=None, dd=None):
    if dim is None:
        dim = ambient_dim - 1

    if ambient_dim == dim:
        inv_mder = inverse_metric_derivative_mat(ambient_dim, dim=dim, dd=dd)
        inv_form1 = inv_mder.dot(inv_mder.T)
    else:
        form1 = first_fundamental_form(ambient_dim, dim, dd)

        if dim == 1:
            inv_form1 = np.array([[1.0 / form1[0, 0]]], dtype=object)
        elif dim == 2:
            (E, F), (_, G) = form1  # noqa: N806
            inv_form1 = 1.0 / (E * G - F * F) * np.array([[G, -F], [-F, E]],
                                                         dtype=object)
        else:
            raise ValueError("%dD surfaces not supported" % dim)

    return cse(inv_form1, "inv_form1_mat", cse_scope.DISCRETIZATION)
Esempio n. 32
0
def shape_operator(ambient_dim, dim=None, where=None):
    if dim is None:
        dim = ambient_dim - 1

    if ambient_dim != 3 and dim != 2:
        raise NotImplementedError("only available for surfaces in 3D")

    # https://en.wikipedia.org/w/index.php?title=Differential_geometry_of_surfaces&oldid=833587563
    (E, F), (F, G) = first_fundamental_form(ambient_dim, dim, where)
    (e, f), (f, g) = second_fundamental_form(ambient_dim, dim, where)

    result = np.zeros((2, 2), dtype=object)
    result[0, 0] = e*G-f*F
    result[0, 1] = f*G-g*F
    result[1, 0] = f*E-e*F
    result[1, 1] = g*E-f*F

    return cse(
            1/(E*G-F*F)*result,
            "shape_operator")
    def bdry_flux(self, q_bdry, q_null, tag):
        if tag == self.boundaryconditions_tag['stressfree']:
            signP = -1
            signv = 1
        elif tag == self.boundaryconditions_tag['fixed']:
            signP = 1
            signv = -1
        else:
            raise ValueError("Invalid boundary conditions")

        dim = self.dimensions

        P = self.P(q_bdry)
        v = self.v(q_bdry)
        v_null = q_null

        # One entry for each flux direction
        if dim == 1:
            return [cse(join_fields(
                        v_null,
                        signP*P[0], # flux rho_v
                        signv*v[0]  # flux F
                        ), "x_bflux")]
        elif dim == 2:
            return [cse(join_fields(
                        v_null, signP*P[0],signP*P[2],       # flux rho_v
                        signv*v[0],v_null,v_null,signv*v[1]  # flux F
                        ), "x_bflux"),
                    cse(join_fields(
                        v_null, signP*P[2],signP*P[1],       # flux rho_v
                        v_null,signv*v[1],signv*v[0],v_null  # flux F
                        ), "y_bflux")]
        elif dim == 3:
            return [cse(join_fields(
                        v_null, signP*P[0],signP*P[5],signP*P[4],                                  # flux rho_v
                        signv*v[0],v_null,v_null,v_null,v_null,v_null,signv*v[2],v_null,signv*v[1] # flux F
                        ), "x_bflux"),
                    cse(join_fields(
                        v_null, signP*P[5],signP*P[1],signP*P[3],                                  # flux rho_v
                        v_null,signv*v[1],v_null,v_null,signv*v[2],v_null,v_null,signv*v[0],v_null # flux F
                        ), "y_bflux"),
                    cse(join_fields(
                        v_null, signP*P[4],signP*P[3],signP*P[2],                                  # flux rho_v
                        v_null,v_null,signv*v[2],signv*v[1],v_null,signv*v[0],v_null,v_null,v_null # flux F
                        ), "z_bflux")]
        else:
            raise ValueError("Invalid dimension")
Esempio n. 34
0
def mv_normal(dd, ambient_dim, dim=None):
    """Exterior unit normal as a MultiVector."""

    dd = as_dofdesc(dd)

    if not dd.is_trace():
        raise ValueError("may only request normals on boundaries")

    if dim is None:
        dim = ambient_dim - 1

    # Don't be tempted to add a sign here. As it is, it produces
    # exterior normals for positively oriented curves.

    pder = (
        pseudoscalar(ambient_dim, dim, dd=dd) /  # noqa: W504
        area_element(ambient_dim, dim, dd=dd))
    return cse(
        # Dorst Section 3.7.2
        pder << pder.I.inv(),
        "normal",
        cse_scope.DISCRETIZATION)
Esempio n. 35
0
def forward_metric_derivative(xyz_axis, rst_axis, dd=None):
    r"""
    Pointwise metric derivatives representing

    .. math::

        \frac{d x_{\mathtt{xyz\_axis}} }{d r_{\mathtt{rst\_axis}} }
    """
    if dd is None:
        dd = DD_VOLUME

    inner_dd = dd.with_qtag(QTAG_NONE)

    diff_op = _sym().RefDiffOperator(rst_axis, inner_dd)

    result = diff_op(NodeCoordinateComponent(xyz_axis, inner_dd))

    if dd.uses_quadrature():
        result = _sym().interp(inner_dd, dd)(result)

    return cse(result,
               prefix="dx%d_dr%d" % (xyz_axis, rst_axis),
               scope=cse_scope.DISCRETIZATION)
Esempio n. 36
0
def second_fundamental_form(ambient_dim, dim=None, dd=None):
    if dim is None:
        dim = ambient_dim - 1

    normal = surface_normal(ambient_dim, dim=dim, dd=dd).as_vector()
    if dim == 1:
        second_ref_axes = [((0, 2),)]
    elif dim == 2:
        second_ref_axes = [((0, 2),), ((0, 1), (1, 1)), ((1, 2),)]
    else:
        raise ValueError("%dD surfaces not supported" % dim)

    from pytools import flatten
    form2 = np.empty((dim, dim), dtype=object)
    for ref_axes in second_ref_axes:
        i, j = flatten([rst_axis] * n for rst_axis, n in ref_axes)

        ruv = np.array([
            forward_metric_nth_derivative(xyz_axis, ref_axes, dd=dd)
            for xyz_axis in range(ambient_dim)])
        form2[i, j] = form2[j, i] = normal.dot(ruv)

    return cse(form2, "form2_mat", cse_scope.DISCRETIZATION)
    def flux(self, q):
        P = self.P(q)
        v = self.v(q)
        v_null = Field('state_null')
        dim = self.dimensions

        # One entry for each flux direction
        if dim == 1:
            return [cse(join_fields(
                        v_null,
                        P[0],  # flux rho_v
                        v[0]   # flux F
                        ), "x_flux")]

        elif dim == 2:
            return [cse(join_fields(
                        v_null, P[0],P[2],# flux rho_v
                        v[0],v_null,v[1]  # flux F
                        ), "x_flux"),
                    cse(join_fields(
                        v_null, P[2],P[1],# flux rho_v
                        v_null,v[1],v[0]  # flux F
                        ), "y_flux")]
        elif dim == 3:
            return [cse(join_fields(
                        v_null, P[0],P[5],P[4],             # flux rho_v
                        v[0],v_null,v_null,v_null,v[2],v[1] # flux F
                        ), "x_flux"),
                    cse(join_fields(
                        v_null, P[5],P[1],P[3],             # flux rho_v
                        v_null,v[1],v_null,v[2],v_null,v[0] # flux F
                        ), "y_flux"),
                    cse(join_fields(
                        v_null, P[4],P[3],P[2],             # flux rho_v
                        v_null,v_null,v[2],v[1],v[0],v_null # flux F
                        ), "z_flux")]
        else:
            raise ValueError("Invalid dimension")
Esempio n. 38
0
def pseudoscalar(where=None):
    return cse(
            ParametrizationDerivative(where).a.project_max_grade(),
            "pseudoscalar", cse_scope.DISCRETIZATION)
Esempio n. 39
0
def area(where=None):
    return cse(integral(Ones(where), where), "area",
            cse_scope.DISCRETIZATION)
Esempio n. 40
0
def surf_grad_S(kernel, arg, dim):
    """
    :arg dim: The dimension of the ambient space.
    """

    return project_to_tangential(cse(grad_S(kernel, arg, dim)))
Esempio n. 41
0
def area_element(ambient_dim, dim=None, dd=None):
    return cse(sqrt(pseudoscalar(ambient_dim, dim, dd=dd).norm_squared()),
               "area_element", cse_scope.DISCRETIZATION)
Esempio n. 42
0
def project_to_tangential(xyz_vec, which=None):
    return tangential_to_xyz(
            cse(xyz_to_tangential(xyz_vec, which), which))
Esempio n. 43
0
def project_to_tangential(xyz_vec, where=None):
    return tangential_to_xyz(
            cse(xyz_to_tangential(xyz_vec, where), where))
Esempio n. 44
0
def area_element(ambient_dim, dim=None, where=None):
    return cse(
            sqrt(pseudoscalar(ambient_dim, dim, where).norm_squared()),
            "area_element", cse_scope.DISCRETIZATION)
Esempio n. 45
0
def sqrt_jac_q_weight(ambient_dim, dim=None, where=None):
    return cse(
            sqrt(
                area_element(ambient_dim, dim, where)
                * QWeight(where)),
            "sqrt_jac_q_weight", cse_scope.DISCRETIZATION)
Esempio n. 46
0
def first_fundamental_form(ambient_dim, dim=None, dd=None):
    if dim is None:
        dim = ambient_dim - 1

    mder = forward_metric_derivative_mat(ambient_dim, dim=dim, dd=dd)
    return cse(mder.T.dot(mder), "form1_mat", cse_scope.DISCRETIZATION)
Esempio n. 47
0
def area_element(where=None):
    return cse(
            sqrt(pseudoscalar(where).a.norm_squared()),
            "area_element", cse_scope.DISCRETIZATION)
Esempio n. 48
0
def area(ambient_dim, dim, where=None):
    return cse(integral(ambient_dim, dim, Ones(where), where), "area",
            cse_scope.DISCRETIZATION)