예제 #1
0
    def circumradius(self, o):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()
        if not domain.is_piecewise_linear_simplex_domain():
            # Don't lower for non-affine cells, instead leave it to
            # form compiler
            warning(
                "Only know how to compute the circumradius of an affine cell.")
            return o

        cellname = domain.ufl_cell().cellname()
        cellvolume = self.cell_volume(CellVolume(domain))

        if cellname == "interval":
            r = 0.5 * cellvolume

        elif cellname == "triangle":
            J = self.jacobian(Jacobian(domain))
            trev = CellEdgeVectors(domain)
            num_edges = 3
            i, j, k = indices(3)
            elen = [
                sqrt((J[i, j] * trev[edge, j]) * (J[i, k] * trev[edge, k]))
                for edge in range(num_edges)
            ]

            r = (elen[0] * elen[1] * elen[2]) / (4.0 * cellvolume)

        elif cellname == "tetrahedron":
            J = self.jacobian(Jacobian(domain))
            trev = CellEdgeVectors(domain)
            num_edges = 6
            i, j, k = indices(3)
            elen = [
                sqrt((J[i, j] * trev[edge, j]) * (J[i, k] * trev[edge, k]))
                for edge in range(num_edges)
            ]

            # elen[3] = length of edge 3
            # la, lb, lc = lengths of the sides of an intermediate triangle
            la = elen[3] * elen[2]
            lb = elen[4] * elen[1]
            lc = elen[5] * elen[0]
            # p = perimeter
            p = (la + lb + lc)
            # s = semiperimeter
            s = p / 2
            # area of intermediate triangle with Herons formula
            triangle_area = sqrt(s * (s - la) * (s - lb) * (s - lc))
            r = triangle_area / (6.0 * cellvolume)

        else:
            error("Unhandled cell type %s." % cellname)

        return r
예제 #2
0
def apply_known_single_pullback(r, element):
    """Apply pullback with given mapping.

    :arg r: Expression wrapped in ReferenceValue
    :arg element: The element defining the mapping
    """
    # Need to pass in r rather than the physical space thing, because
    # the latter may be a ListTensor or similar, rather than a
    # Coefficient/Argument (in the case of mixed elements, see below
    # in apply_single_function_pullbacks), to which we cannot apply ReferenceValue
    mapping = element.mapping()
    domain = r.ufl_domain()
    if mapping == "physical":
        return r
    elif mapping == "identity":
        return r
    elif mapping == "contravariant Piola":
        J = Jacobian(domain)
        detJ = JacobianDeterminant(J)
        transform = (1.0 / detJ) * J
        # Apply transform "row-wise" to TensorElement(PiolaMapped, ...)
        *k, i, j = indices(len(r.ufl_shape) + 1)
        kj = (*k, j)
        f = as_tensor(transform[i, j] * r[kj], (*k, i))
        return f
    elif mapping == "covariant Piola":
        K = JacobianInverse(domain)
        # Apply transform "row-wise" to TensorElement(PiolaMapped, ...)
        *k, i, j = indices(len(r.ufl_shape) + 1)
        kj = (*k, j)
        f = as_tensor(K[j, i] * r[kj], (*k, i))
        return f
    elif mapping == "L2 Piola":
        detJ = JacobianDeterminant(domain)
        return r / detJ
    elif mapping == "double contravariant Piola":
        J = Jacobian(domain)
        detJ = JacobianDeterminant(J)
        transform = (1.0 / detJ) * J
        # Apply transform "row-wise" to TensorElement(PiolaMapped, ...)
        *k, i, j, m, n = indices(len(r.ufl_shape) + 2)
        kmn = (*k, m, n)
        f = as_tensor((1.0 / detJ)**2 * J[i, m] * r[kmn] * J[j, n], (*k, i, j))
        return f
    elif mapping == "double covariant Piola":
        K = JacobianInverse(domain)
        # Apply transform "row-wise" to TensorElement(PiolaMapped, ...)
        *k, i, j, m, n = indices(len(r.ufl_shape) + 2)
        kmn = (*k, m, n)
        f = as_tensor(K[m, i] * r[kmn] * K[n, j], (*k, i, j))
        return f
    else:
        error("Should never be reached!")
예제 #3
0
파일: fem.py 프로젝트: jmv2009/tsfc
 def jacobian_at(self, point):
     ps = PointSingleton(point)
     expr = Jacobian(self.mt.terminal.ufl_domain())
     assert ps.expression.shape == (
         expr.ufl_domain().topological_dimension(), )
     if self.mt.restriction == '+':
         expr = PositiveRestricted(expr)
     elif self.mt.restriction == '-':
         expr = NegativeRestricted(expr)
     config = {"point_set": PointSingleton(point)}
     config.update(self.config)
     context = PointSetContext(**config)
     expr = self.preprocess(expr, context)
     return map_expr_dag(context.translator, expr)
예제 #4
0
    def cell_normal(self, o):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()
        gdim = domain.geometric_dimension()
        tdim = domain.topological_dimension()

        if tdim == gdim - 1:  # n-manifold embedded in n-1 space
            i = Index()
            J = self.jacobian(Jacobian(domain))

            if tdim == 2:
                # Surface in 3D
                t0 = as_vector(J[i, 0], i)
                t1 = as_vector(J[i, 1], i)
                cell_normal = cross_expr(t0, t1)
            elif tdim == 1:
                # Line in 2D (cell normal is 'up' for a line pointing
                # to the 'right')
                cell_normal = as_vector((-J[1, 0], J[0, 0]))
            else:
                error("Cell normal not implemented for tdim %d, gdim %d" %
                      (tdim, gdim))

            # Return normalized vector, sign corrected by cell
            # orientation
            co = CellOrientation(domain)
            return co * cell_normal / sqrt(cell_normal[i] * cell_normal[i])
        else:
            error("What do you want cell normal in gdim={0}, tdim={1} to be?".
                  format(gdim, tdim))
예제 #5
0
    def max_facet_edge_length(self, o):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()
        if not domain.is_piecewise_linear_simplex_domain():
            # Don't lower for non-affine cells, instead leave it to
            # form compiler
            warning(
                "Only know how to compute the max_facet_edge_length of an affine cell."
            )
            return o

        cellname = domain.ufl_cell().cellname()

        if cellname == "triangle":
            return self.facet_area(FacetArea(domain))
        elif cellname == "tetrahedron":
            J = self.jacobian(Jacobian(domain))
            trev = FacetEdgeVectors(domain)
            num_edges = 3
            i, j, k = indices(3)
            elen = [
                sqrt((J[i, j] * trev[edge, j]) * (J[i, k] * trev[edge, k]))
                for edge in range(num_edges)
            ]
            return max_value(elen[0], max_value(elen[1], elen[2]))
        else:
            error("Unhandled cell type %s." % cellname)
예제 #6
0
    def facet_jacobian(self, o):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()
        J = self.jacobian(Jacobian(domain))
        RFJ = CellFacetJacobian(domain)
        i, j, k = indices(3)
        return as_tensor(J[i, k] * RFJ[k, j], (i, j))
예제 #7
0
    def jacobian_inverse(self, o):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()
        J = self.jacobian(Jacobian(domain))
        # TODO: This could in principle use
        # preserve_types[JacobianDeterminant] with minor refactoring:
        K = inverse_expr(J)
        return K
예제 #8
0
파일: fem.py 프로젝트: knut0815/tsfc
    def jacobian_at(self, point):
        expr = Jacobian(self.mt.terminal.ufl_domain())
        if self.mt.restriction == '+':
            expr = PositiveRestricted(expr)
        elif self.mt.restriction == '-':
            expr = NegativeRestricted(expr)
        expr = preprocess_expression(expr)

        config = {"point_set": PointSingleton(point)}
        config.update(self.config)
        context = PointSetContext(**config)
        return map_expr_dag(context.translator, expr)
예제 #9
0
    def jacobian_determinant(self, o):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()
        J = self.jacobian(Jacobian(domain))
        detJ = determinant_expr(J)

        # TODO: Is "signing" the determinant for manifolds the
        #       cleanest approach?  The alternative is to have a
        #       specific type for the unsigned pseudo-determinant.
        if domain.topological_dimension() < domain.geometric_dimension():
            co = CellOrientation(domain)
            detJ = co * detJ

        return detJ
예제 #10
0
    def facet_normal(self, o):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()
        tdim = domain.topological_dimension()

        if tdim == 1:
            # Special-case 1D (possibly immersed), for which we say
            # that n is just in the direction of J.
            J = self.jacobian(Jacobian(domain))  # dx/dX
            ndir = J[:, 0]

            gdim = domain.geometric_dimension()
            if gdim == 1:
                nlen = abs(ndir[0])
            else:
                i = Index()
                nlen = sqrt(ndir[i] * ndir[i])

            rn = ReferenceNormal(domain)  # +/- 1.0 here
            n = rn[0] * ndir / nlen
            r = n
        else:
            # Recall that the covariant Piola transform u -> J^(-T)*u
            # preserves tangential components. The normal vector is
            # characterised by having zero tangential component in
            # reference and physical space.
            Jinv = self.jacobian_inverse(JacobianInverse(domain))
            i, j = indices(2)

            rn = ReferenceNormal(domain)
            # compute signed, unnormalised normal; note transpose
            ndir = as_vector(Jinv[j, i] * rn[j], i)

            # normalise
            i = Index()
            n = ndir / sqrt(ndir[i] * ndir[i])
            r = n

        if r.ufl_shape != o.ufl_shape:
            error("Inconsistent dimensions (in=%d, out=%d)." %
                  (o.ufl_shape[0], r.ufl_shape[0]))
        return r
예제 #11
0
def apply_single_function_pullbacks(g):
    element = g.ufl_element()
    mapping = element.mapping()

    r = ReferenceValue(g)
    gsh = g.ufl_shape
    rsh = r.ufl_shape

    if mapping == "physical":
        # TODO: Is this right for immersed things
        assert gsh == rsh
        return r

    # Shortcut the "identity" case which includes Expression and
    # Constant from dolfin that may be ill-formed without a domain
    # (until we get that fixed)
    if mapping == "identity":
        assert rsh == gsh
        return r

    gsize = product(gsh)
    rsize = product(rsh)

    # Create some geometric objects for reuse
    domain = g.ufl_domain()
    J = Jacobian(domain)
    detJ = JacobianDeterminant(domain)
    Jinv = JacobianInverse(domain)

    # Create contravariant transform for reuse (note that detJ is the
    # _signed_ (pseudo-)determinant)
    transform_hdiv = (1.0 / detJ) * J

    # Shortcut simple cases for a more efficient representation,
    # including directly Piola-mapped elements and mixed elements of
    # any combination of affinely mapped elements without symmetries
    if mapping == "symmetries":
        fcm = element.flattened_sub_element_mapping()
        assert gsize >= rsize
        assert len(fcm) == gsize
        assert sorted(set(fcm)) == sorted(range(rsize))
        g_components = [r[fcm[i]] for i in range(gsize)]
        g_components = reshape_to_nested_list(g_components, gsh)
        f = as_tensor(g_components)
        assert f.ufl_shape == g.ufl_shape
        return f
    elif mapping == "contravariant Piola":
        assert transform_hdiv.ufl_shape == (gsize, rsize)
        i, j = indices(2)
        f = as_vector(transform_hdiv[i, j] * r[j], i)
        # f = as_tensor(transform_hdiv[i, j]*r[k,j], (k,i)) # FIXME: Handle Vector(Piola) here?
        assert f.ufl_shape == g.ufl_shape
        return f
    elif mapping == "covariant Piola":
        assert Jinv.ufl_shape == (rsize, gsize)
        i, j = indices(2)
        f = as_vector(Jinv[j, i] * r[j], i)
        # f = as_tensor(Jinv[j, i]*r[k,j], (k,i)) # FIXME: Handle Vector(Piola) here?
        assert f.ufl_shape == g.ufl_shape
        return f
    elif mapping == "double covariant Piola":
        i, j, m, n = indices(4)
        f = as_tensor(Jinv[m, i] * r[m, n] * Jinv[n, j], (i, j))
        assert f.ufl_shape == g.ufl_shape
        return f
    elif mapping == "double contravariant Piola":
        i, j, m, n = indices(4)
        f = as_tensor(
            (1.0 / detJ) * (1.0 / detJ) * J[i, m] * r[m, n] * J[j, n], (i, j))
        assert f.ufl_shape == g.ufl_shape
        return f
    elif mapping == "L2 Piola":
        assert rsh == gsh
        return r / detJ

    # By placing components in a list and using as_vector at the end,
    # we're assuming below that both global function g and its
    # reference value r have vector shape, which is the case for most
    # elements with the exceptions:
    # - TensorElements
    #   - All cases with scalar subelements and without symmetries
    #     are covered by the shortcut above
    #     (ONLY IF REFERENCE VALUE SHAPE PRESERVES TENSOR RANK)
    #   - All cases with scalar subelements and without symmetries are
    #     covered by the shortcut above

    g_components = [None] * gsize
    gpos = 0
    rpos = 0

    r = as_vector([r[idx] for idx in numpy.ndindex(r.ufl_shape)])
    for subelm in sub_elements_with_mappings(element):
        gm = product(subelm.value_shape())
        rm = product(subelm.reference_value_shape())

        mp = subelm.mapping()
        if mp == "identity":
            assert gm == rm
            for i in range(gm):
                g_components[gpos + i] = r[rpos + i]

        elif mp == "symmetries":
            """
            tensor_element.value_shape() == (2,2)
            tensor_element.reference_value_shape() == (3,)
            tensor_element.symmetry() == { (1,0): (0,1) }
            tensor_element.component_mapping() == { (0,0): 0, (0,1): 1, (1,0): 1, (1,1): 2 }
            tensor_element.flattened_component_mapping() == { 0: 0, 1: 1, 2: 1, 3: 2 }
            """
            fcm = subelm.flattened_sub_element_mapping()
            assert gm >= rm
            assert len(fcm) == gm
            assert sorted(set(fcm)) == sorted(range(rm))
            for i in range(gm):
                g_components[gpos + i] = r[rpos + fcm[i]]

        elif mp == "contravariant Piola":
            assert transform_hdiv.ufl_shape == (gm, rm)
            # Get reference value vector corresponding to this subelement:
            rv = as_vector([r[rpos + k] for k in range(rm)])
            # Apply transform with IndexSum over j for each row
            j = Index()
            for i in range(gm):
                g_components[gpos + i] = transform_hdiv[i, j] * rv[j]

        elif mp == "covariant Piola":
            assert Jinv.ufl_shape == (rm, gm)
            # Get reference value vector corresponding to this subelement:
            rv = as_vector([r[rpos + k] for k in range(rm)])
            # Apply transform with IndexSum over j for each row
            j = Index()
            for i in range(gm):
                g_components[gpos + i] = Jinv[j, i] * rv[j]

        elif mp == "double covariant Piola":
            # components are flatten, map accordingly
            rv = as_vector([r[rpos + k] for k in range(rm)])
            (gdim, _) = subelm.value_shape()
            (rdim, _) = subelm.reference_value_shape()
            for i in range(gdim):
                for j in range(gdim):
                    gv = 0
                    # int times Index is not allowed. so sum by hand
                    for m in range(rdim):
                        for n in range(rdim):
                            gv += Jinv[m, i] * rv[m * rdim + n] * Jinv[n, j]
                    g_components[gpos + i * gdim + j] = gv

        elif mp == "double contravariant Piola":
            # components are flatten, map accordingly
            rv = as_vector([r[rpos + k] for k in range(rm)])
            (gdim, _) = subelm.value_shape()
            (rdim, _) = subelm.reference_value_shape()
            for i in range(gdim):
                for j in range(gdim):
                    gv = 0
                    # int times Index is not allowed. so sum by hand
                    for m in range(rdim):
                        for n in range(rdim):
                            gv += ((1.0 / detJ) * (1.0 / detJ) * J[i, m] *
                                   rv[m * rdim + n] * J[j, n])
                    g_components[gpos + i * gdim + j] = gv

        elif mp == "L2 Piola":
            assert gm == rm
            for i in range(gm):
                g_components[gpos + i] = r[rpos + i] / detJ

        else:
            error("Unknown subelement mapping type %s for element %s." %
                  (mp, str(subelm)))

        gpos += gm
        rpos += rm

    # Wrap up components in a vector, must return same shape as input
    # function g
    f = as_tensor(numpy.asarray(g_components).reshape(gsh))
    assert f.ufl_shape == g.ufl_shape
    return f
예제 #12
0
def test_triangle_symbolic_geometry(uflacs_representation_only):
    mesh = UnitSquareMesh(1, 1)
    assert mesh.num_cells() == 2
    gdim = mesh.geometry().dim()
    tdim = mesh.topology().dim()

    area = 1.0  # known volume of mesh
    A = area / 2.0  # volume of single cell
    Aref = 1.0 / 2.0  # the volume of the UFC reference triangle

    mf = MeshFunction("size_t", mesh, mesh.topology().dim())
    dx = Measure("dx", domain=mesh, subdomain_data=mf)

    U0 = FunctionSpace(mesh, "DG", 0)
    V0 = VectorFunctionSpace(mesh, "DG", 0)
    U1 = FunctionSpace(mesh, "DG", 1)
    V1 = VectorFunctionSpace(mesh, "DG", 1)

    # Check coordinates with various consistency checking
    x = SpatialCoordinate(mesh)
    X = CellCoordinate(mesh)
    J = Jacobian(mesh)
    detJ = JacobianDeterminant(mesh)
    K = JacobianInverse(mesh)
    vol = CellVolume(mesh)
    h = CellDiameter(mesh)
    R = Circumradius(mesh)

    coordinates = mesh.coordinates()
    cells = mesh.cells()

    for k in range(mesh.num_cells()):
        # Mark current cell
        mf.set_all(0)
        mf[k] = 1

        # Check integration area vs detJ
        # This is not currently implemented in uflacs:
        # x0 = CellOrigin(mesh)
        # But we can extract it from the mesh for a given cell k
        x0 = as_vector(coordinates[cells[k][0]][:])
        # Validate known cell volume A
        assert round(assemble(1.0 * dx(1)) - A, 7) == 0.0
        assert round(assemble(1.0 / A * dx(1)) - 1.0, 7) == 0.0
        assert round(assemble(A * dx(1)) - A**2, 7) == 0.0
        # Compare abs(detJ) to A
        A2 = Aref * abs(detJ)
        assert round(assemble((A - A2)**2 * dx(1)) - 0.0, 7) == 0.0
        assert round(assemble(1.0 / A2 * dx(1)) - 1.0, 7) == 0.0
        assert round(assemble(A2 * dx(1)) - A**2, 7) == 0.0
        # Compare vol to A
        A4 = vol
        assert round(assemble((A - A4)**2 * dx(1)) - 0.0, 7) == 0.0
        assert round(assemble((A2 - A4)**2 * dx(1)) - 0.0, 7) == 0.0
        assert round(assemble(1.0 / A4 * dx(1)) - 1.0, 7) == 0.0
        assert round(assemble(A4 * dx(1)) - A**2, 7) == 0.0

        # Check integral of reference coordinate components over reference
        # triangle:
        Xmp = (1.0 / 6.0, 1.0 / 6.0)
        for j in range(tdim):
            # Scale by detJ^-1 to get reference cell integral
            assert round(assemble(X[j] / abs(detJ) * dx(1)) - Xmp[j], 7) == 0.0

        # Check average of physical coordinate components over each cell:
        # Compute average of vertex coordinates extracted from mesh
        verts = [coordinates[i][:] for i in cells[k]]
        vavg = sum(verts[1:], verts[0]) / len(verts)
        for i in range(gdim):
            # Scale by A^-1 to get average of x, not integral
            assert round(assemble(x[i] / A * dx(1)) - vavg[i], 7) == 0.0

        # Check affine coordinate relations x=x0+J*X, X=K*(x-x0), K*J=I
        x0 = as_vector(coordinates[cells[k][0]][:])
        assert round(assemble((x - (x0 + J * X))**2 * dx(1)), 7) == 0.0
        assert round(assemble((X - K * (x - x0))**2 * dx(1)), 7) == 0.0
        assert round(assemble((K * J - Identity(tdim))**2 / A * dx(1)),
                     7) == 0.0

        # Check cell diameter and circumradius
        assert round(assemble(h / vol * dx(1)) - Cell(mesh, k).h(), 7) == 0.0
        assert round(
            assemble(R / vol * dx(1)) - Cell(mesh, k).circumradius(), 7) == 0.0
예제 #13
0
def test_manifold_piola_mapped_functions(square3d, any_representation):
    mesh = square3d
    area = sqrt(3.0)  # known area of mesh
    A = area / 2.0

    mf = MeshFunction("size_t", mesh, mesh.topology().dim())
    mf[0] = 0
    mf[1] = 1
    dx = Measure("dx", domain=mesh, subdomain_data=mf)

    x = SpatialCoordinate(mesh)

    J = Jacobian(mesh)
    detJ = JacobianDeterminant(mesh)  # pseudo-determinant
    K = JacobianInverse(mesh)  # pseudo-inverse

    Q1 = VectorFunctionSpace(mesh, "CG", 1)
    U1 = VectorFunctionSpace(mesh, "DG", 1)
    V1 = FunctionSpace(mesh, "N1div", 1)
    W1 = FunctionSpace(mesh, "N1curl", 1)

    dq = TestFunction(Q1)
    du = TestFunction(U1)
    dv = TestFunction(V1)
    dw = TestFunction(W1)

    assert U1.ufl_element().mapping() == "identity"
    assert V1.ufl_element().mapping() == "contravariant Piola"
    assert W1.ufl_element().mapping() == "covariant Piola"

    if any_representation != "uflacs":
        return

    # Check that projection test fails if it should fail:
    vec = Constant((0.0, 0.0, 0.0))
    q1 = project(vec, Q1)
    u1 = project(vec, U1)
    v1 = project(vec, V1)
    w1 = project(vec, W1)
    # Projection of zero gets us zero for all spaces
    assert assemble(q1**2 * dx) == 0.0
    assert assemble(u1**2 * dx) == 0.0
    assert assemble(v1**2 * dx) == 0.0
    assert assemble(w1**2 * dx) == 0.0
    # Changing vec to nonzero, check that dM/df != 0 at f=0
    vec = Constant((2.0, 2.0, 2.0))
    assert round(
        assemble(derivative((q1 - vec)**2 * dx, q1)).norm('l2') -
        assemble(-4.0 * sum(dq) * dx).norm('l2'), 7) == 0.0
    assert round(
        assemble(derivative((u1 - vec)**2 * dx, u1)).norm('l2') -
        assemble(-4.0 * sum(du) * dx).norm('l2'), 7) == 0.0
    assert round(
        assemble(derivative((v1 - vec)**2 * dx, v1)).norm('l2') -
        assemble(-4.0 * sum(dv) * dx).norm('l2'), 7) == 0.0
    assert round(
        assemble(derivative((w1 - vec)**2 * dx, w1)).norm('l2') -
        assemble(-4.0 * sum(dw) * dx).norm('l2'), 7) == 0.0

    # Project piecewise linears to scalar and vector CG1 spaces on
    # manifold
    vec = Constant((1.0, 1.0, 1.0))
    q1 = project(vec, Q1)
    u1 = project(vec, U1)
    v1 = project(vec, V1)
    w1 = project(vec, W1)

    # If vec can be represented exactly in space this should be zero:
    assert round(assemble((q1 - vec)**2 * dx), 7) == 0.0
    assert round(assemble((u1 - vec)**2 * dx), 7) == 0.0
    assert round(assemble(
        (v1 - vec)**2 * dx), 7) > 0.0  # Exact representation not possible?
    assert round(assemble(
        (w1 - vec)**2 * dx), 7) > 0.0  # Exact representation not possible?

    # In the l2norm projection is correct these should be zero:
    assert round(assemble(derivative((q1 - vec)**2 * dx, v1)).norm('l2'),
                 7) == 0.0
    assert round(assemble(derivative((u1 - vec)**2 * dx, w1)).norm('l2'),
                 7) == 0.0
    assert round(assemble(derivative((v1 - vec)**2 * dx, v1)).norm('l2'),
                 7) == 0.0
    assert round(assemble(derivative((w1 - vec)**2 * dx, w1)).norm('l2'),
                 7) == 0.0

    # Hdiv mapping of a local constant vector should be representable
    # in hdiv conforming space
    vec = (1.0 / detJ) * J * as_vector((3.0, 5.0))
    q1 = project(vec, Q1)
    u1 = project(vec, U1)
    v1 = project(vec, V1)
    w1 = project(vec, W1)
    assert round(assemble(
        (q1 - vec)**2 * dx), 7) > 0.0  # Exact representation not possible?
    assert round(assemble((u1 - vec)**2 * dx), 7) == 0.0
    assert round(assemble((v1 - vec)**2 * dx), 7) == 0.0
    assert round(assemble(
        (w1 - vec)**2 * dx), 7) > 0.0  # Exact representation not possible?

    # Hcurl mapping of a local constant vector should be representable
    # in hcurl conforming space
    vec = K.T * as_vector((5.0, 2.0))
    q1 = project(vec, Q1)
    u1 = project(vec, U1)
    v1 = project(vec, V1)
    w1 = project(vec, W1)
    assert round(assemble(
        (q1 - vec)**2 * dx), 7) > 0.0  # Exact representation not possible?
    assert round(assemble((u1 - vec)**2 * dx), 7) == 0.0
    assert round(assemble(
        (v1 - vec)**2 * dx), 7) > 0.0  # Exact representation not possible?
    assert round(assemble((w1 - vec)**2 * dx), 7) == 0.0
예제 #14
0
def test_manifold_symbolic_geometry(square3d, uflacs_representation_only):
    mesh = square3d
    assert mesh.num_cells() == 2
    gdim = mesh.geometry().dim()
    tdim = mesh.topology().dim()

    area = sqrt(3.0)  # known area of mesh
    A = area / 2.0  # area of single cell
    Aref = 0.5  # 0.5 is the area of the UFC reference triangle

    mf = MeshFunction("size_t", mesh, mesh.topology().dim())
    mf[0] = 0
    mf[1] = 1
    dx = Measure("dx", domain=mesh, subdomain_data=mf)

    U0 = FunctionSpace(mesh, "DG", 0)
    V0 = VectorFunctionSpace(mesh, "DG", 0)

    # 0 means up=+1.0, 1 means down=-1.0
    orientations = mesh.cell_orientations()
    assert orientations[0] == 1  # down
    assert orientations[1] == 0  # up

    # Check cell orientation, should be -1.0 (down) and +1.0 (up) on
    # the two cells respectively by construction
    co = CellOrientation(mesh)
    co0 = assemble(co / A * dx(0))
    co1 = assemble(co / A * dx(1))
    assert round(abs(co0) - 1.0, 7) == 0.0  # should be +1 or -1
    assert round(abs(co1) - 1.0, 7) == 0.0  # should be +1 or -1
    assert round(co1 + co0, 7) == 0.0  # should cancel out
    assert round(co0 - -1.0, 7) == 0.0  # down
    assert round(co1 - +1.0, 7) == 0.0  # up

    # Check cell normal directions component for component
    cn = CellNormal(mesh)
    assert assemble(cn[0] / A * dx(0)) > 0.0
    assert assemble(cn[0] / A * dx(1)) < 0.0
    assert assemble(cn[1] / A * dx(0)) < 0.0
    assert assemble(cn[1] / A * dx(1)) > 0.0
    assert assemble(cn[2] / A * dx(0)) > 0.0
    assert assemble(cn[2] / A * dx(1)) > 0.0
    # Check cell normal normalization
    assert round(assemble(cn**2 / A * dx(0)) - 1.0, 7) == 0.0
    assert round(assemble(cn**2 / A * dx(1)) - 1.0, 7) == 0.0

    # Check coordinates with various consistency checking
    x = SpatialCoordinate(mesh)
    X = CellCoordinate(mesh)
    J = Jacobian(mesh)
    detJ = JacobianDeterminant(mesh)  # pseudo-determinant
    K = JacobianInverse(mesh)  # pseudo-inverse
    vol = CellVolume(mesh)
    h = CellDiameter(mesh)
    R = Circumradius(mesh)

    # This is not currently implemented in uflacs:
    # x0 = CellOrigin(mesh)
    # But by happy accident, x0 is the same vertex for both our triangles:
    x0 = as_vector((0.0, 0.0, 1.0))

    # Checks on each cell separately
    for k in range(mesh.num_cells()):
        # Mark current cell
        mf.set_all(0)
        mf[k] = 1

        # Check integration area vs detJ
        # Validate known cell area A
        assert round(assemble(1.0 * dx(1)) - A, 7) == 0.0
        assert round(assemble(1.0 / A * dx(1)) - 1.0, 7) == 0.0
        assert round(assemble(A * dx(1)) - A**2, 7) == 0.0
        # Compare abs(detJ) to A
        A2 = Aref * abs(detJ)
        assert round(assemble((A - A2)**2 * dx(1)) - 0.0, 7) == 0.0
        assert round(assemble(1.0 / A2 * dx(1)) - 1.0, 7) == 0.0
        assert round(assemble(A2 * dx(1)) - A**2, 7) == 0.0
        # Validate cell orientation
        assert round(assemble(co * dx(1)) - A * (1 if k == 1 else -1),
                     7) == 0.0
        # Compare co*detJ to A (detJ is pseudo-determinant with sign
        # restored, *co again is equivalent to abs())
        A3 = Aref * co * detJ
        assert round(assemble((A - A3)**2 * dx(1)) - 0.0, 7) == 0.0
        assert round(assemble((A2 - A3)**2 * dx(1)) - 0.0, 7) == 0.0
        assert round(assemble(1.0 / A3 * dx(1)) - 1.0, 7) == 0.0
        assert round(assemble(A3 * dx(1)) - A**2, 7) == 0.0
        # Compare vol to A
        A4 = vol
        assert round(assemble((A - A4)**2 * dx(1)) - 0.0, 7) == 0.0
        assert round(assemble((A2 - A4)**2 * dx(1)) - 0.0, 7) == 0.0
        assert round(assemble((A3 - A4)**2 * dx(1)) - 0.0, 7) == 0.0
        assert round(assemble(1.0 / A4 * dx(1)) - 1.0, 7) == 0.0
        assert round(assemble(A4 * dx(1)) - A**2, 7) == 0.0

        # Check cell diameter and circumradius
        assert round(assemble(h / vol * dx(1)) - Cell(mesh, k).h(), 7) == 0.0
        assert round(
            assemble(R / vol * dx(1)) - Cell(mesh, k).circumradius(), 7) == 0.0

        # Check integral of reference coordinate components over reference
        # triangle: \int_0^1 \int_0^{1-x} x dy dx = 1/6
        Xmp = (1.0 / 6.0, 1.0 / 6.0)
        for j in range(tdim):
            # Scale by detJ^-1 to get reference cell integral
            assert round(assemble(X[j] / abs(detJ) * dx(1)) - Xmp[j], 7) == 0.0

        # Check average of physical coordinate components over each cell:
        xmp = [
            (2.0 / 3.0, 1.0 / 3.0, 2.0 / 3.0),  # midpoint of cell 0
            (1.0 / 3.0, 2.0 / 3.0, 2.0 / 3.0),  # midpoint of cell 1
        ]
        for i in range(gdim):
            # Scale by A^-1 to get average of x, not integral
            assert round(assemble(x[i] / A * dx(1)) - xmp[k][i], 7) == 0.0

    # Check affine coordinate relations x=x0+J*X, X=K*(x-x0), K*J=I
    assert round(assemble((x - (x0 + J * X))**2 * dx), 7) == 0.0
    assert round(assemble((X - K * (x - x0))**2 * dx), 7) == 0.0
    assert round(assemble((K * J - Identity(2))**2 / A * dx), 7) == 0.0
예제 #15
0
def test_manifold_line_geometry(mesh, uflacs_representation_only):
    assert uflacs_representation_only == "uflacs"
    assert parameters["form_compiler"]["representation"] == "uflacs"

    gdim = mesh.geometry().dim()
    tdim = mesh.topology().dim()

    # Create cell markers and integration measure
    mf = MeshFunction("size_t", mesh, mesh.topology().dim())
    dx = Measure("dx", domain=mesh, subdomain_data=mf)

    # Create symbolic geometry for current mesh
    x = SpatialCoordinate(mesh)
    X = CellCoordinate(mesh)
    co = CellOrientation(mesh)
    cn = CellNormal(mesh)
    J = Jacobian(mesh)
    detJ = JacobianDeterminant(mesh)
    K = JacobianInverse(mesh)
    vol = CellVolume(mesh)
    h = CellDiameter(mesh)
    R = Circumradius(mesh)

    # Check that length computed via integral doesn't change with
    # refinement
    length = assemble(1.0 * dx)
    mesh2 = refine(mesh)
    assert mesh2.num_cells() == 2 * mesh.num_cells()
    dx2 = Measure("dx")
    length2 = assemble(1.0 * dx2(mesh2))
    assert round(length - length2, 7) == 0.0

    # Check that number of cells can be computed correctly by scaling
    # integral by |detJ|^-1
    num_cells = assemble(1.0 / abs(detJ) * dx)
    assert round(num_cells - mesh.num_cells(), 7) == 0.0

    # Check that norm of Jacobian column matches detJ and volume
    assert round(length - assemble(sqrt(J[:, 0]**2) / abs(detJ) * dx),
                 7) == 0.0
    assert round(assemble((vol - abs(detJ)) * dx), 7) == 0.0
    assert round(length - assemble(vol / abs(detJ) * dx), 7) == 0.0

    coords = mesh.coordinates()
    cells = mesh.cells()

    # Checks on each cell separately
    for k in range(mesh.num_cells()):
        # Mark current cell
        mf.set_all(0)
        mf[k] = 1

        x0 = Constant(tuple(coords[cells[k][0], :]))

        # Integrate x components over a cell and compare with midpoint
        # computed from coords
        for j in range(gdim):
            xm = 0.5 * (coords[cells[k][0], j] + coords[cells[k][1], j])
            assert round(assemble(x[j] / abs(detJ) * dx(1)) - xm, 7) == 0.0

        # Jacobian column is pointing away from x0
        assert assemble(dot(J[:, 0], x - x0) * dx(1)) > 0.0

        # Check affine coordinate relations x=x0+J*X, X=K*(x-x0), K*J=I
        assert round(assemble((x - (x0 + J * X))**2 * dx(1)), 7) == 0.0
        assert round(assemble((X - K * (x - x0))**2 * dx(1)), 7) == 0.0
        assert round(assemble((K * J - Identity(tdim))**2 * dx(1)), 7) == 0.0

        # Check cell diameter and circumradius
        assert round(assemble(h / vol * dx(1)) - Cell(mesh, k).h(), 7) == 0.0
        assert round(
            assemble(R / vol * dx(1)) - Cell(mesh, k).circumradius(), 7) == 0.0

        # Jacobian column is orthogonal to cell normal
        if gdim == 2:
            assert round(assemble(dot(J[:, 0], cn) * dx(1)), 7) == 0.0

            # Create 3d tangent and cell normal vectors
            tangent = as_vector((J[0, 0], J[1, 0], 0.0))
            tangent = co * tangent / sqrt(tangent**2)
            normal = as_vector((cn[0], cn[1], 0.0))
            up = cross(tangent, normal)

            # Check that t,n,up are orthogonal
            assert round(assemble(dot(tangent, normal) * dx(1)), 7) == 0.0
            assert round(assemble(dot(tangent, up) * dx(1)), 7) == 0.0
            assert round(assemble(dot(normal, up) * dx(1)), 7) == 0.0
            assert round(assemble((cross(up, tangent) - normal)**2 * dx(1)),
                         7) == 0.0

            assert round(assemble(up**2 * dx(1)), 7) > 0.0
            assert round(assemble((up[0]**2 + up[1]**2) * dx(1)), 7) == 0.0
            assert round(assemble(up[2] * dx(1)), 7) > 0.0
예제 #16
0
def test_apply_single_function_pullbacks_triangle():
    cell = triangle
    domain = as_domain(cell)

    Ul2 = FiniteElement("DG L2", cell, 1)
    U = FiniteElement("CG", cell, 1)
    V = VectorElement("CG", cell, 1)
    Vd = FiniteElement("RT", cell, 1)
    Vc = FiniteElement("N1curl", cell, 1)
    T = TensorElement("CG", cell, 1)
    S = TensorElement("CG", cell, 1, symmetry=True)

    Uml2 = Ul2 * Ul2
    Um = U * U
    Vm = U * V
    Vdm = V * Vd
    Vcm = Vd * Vc
    Tm = Vc * T
    Sm = T * S

    W = S * T * Vc * Vd * V * U

    ul2 = Coefficient(Ul2)
    u = Coefficient(U)
    v = Coefficient(V)
    vd = Coefficient(Vd)
    vc = Coefficient(Vc)
    t = Coefficient(T)
    s = Coefficient(S)

    uml2 = Coefficient(Uml2)
    um = Coefficient(Um)
    vm = Coefficient(Vm)
    vdm = Coefficient(Vdm)
    vcm = Coefficient(Vcm)
    tm = Coefficient(Tm)
    sm = Coefficient(Sm)

    w = Coefficient(W)

    rul2 = ReferenceValue(ul2)
    ru = ReferenceValue(u)
    rv = ReferenceValue(v)
    rvd = ReferenceValue(vd)
    rvc = ReferenceValue(vc)
    rt = ReferenceValue(t)
    rs = ReferenceValue(s)

    ruml2 = ReferenceValue(uml2)
    rum = ReferenceValue(um)
    rvm = ReferenceValue(vm)
    rvdm = ReferenceValue(vdm)
    rvcm = ReferenceValue(vcm)
    rtm = ReferenceValue(tm)
    rsm = ReferenceValue(sm)

    rw = ReferenceValue(w)

    assert len(w) == 4 + 4 + 2 + 2 + 2 + 1
    assert len(rw) == 3 + 4 + 2 + 2 + 2 + 1
    assert len(w) == 15
    assert len(rw) == 14

    # Geometric quantities we need:
    J = Jacobian(domain)
    detJ = JacobianDeterminant(domain)
    Jinv = JacobianInverse(domain)
    i, j, k, l = indices(4)

    # Contravariant H(div) Piola mapping:
    M_hdiv = (1.0 / detJ) * J
    # Covariant H(curl) Piola mapping: Jinv.T

    mappings = {
        # Simple elements should get a simple representation
        ul2:
        rul2 / detJ,
        u:
        ru,
        v:
        rv,
        vd:
        as_vector(M_hdiv[i, j] * rvd[j], i),
        vc:
        as_vector(Jinv[j, i] * rvc[j], i),
        t:
        rt,
        s:
        as_tensor([[rs[0], rs[1]], [rs[1], rs[2]]]),
        # Mixed elements become a bit more complicated
        uml2:
        as_vector([ruml2[0] / detJ, ruml2[1] / detJ]),
        um:
        rum,
        vm:
        rvm,
        vdm:
        as_vector([
            # V
            rvdm[0],
            rvdm[1],
            # Vd
            M_hdiv[0, j] * as_vector([rvdm[2], rvdm[3]])[j],
            M_hdiv[1, j] * as_vector([rvdm[2], rvdm[3]])[j],
        ]),
        vcm:
        as_vector([
            # Vd
            M_hdiv[0, j] * as_vector([rvcm[0], rvcm[1]])[j],
            M_hdiv[1, j] * as_vector([rvcm[0], rvcm[1]])[j],
            # Vc
            Jinv[i, 0] * as_vector([rvcm[2], rvcm[3]])[i],
            Jinv[i, 1] * as_vector([rvcm[2], rvcm[3]])[i],
        ]),
        tm:
        as_vector([
            # Vc
            Jinv[i, 0] * as_vector([rtm[0], rtm[1]])[i],
            Jinv[i, 1] * as_vector([rtm[0], rtm[1]])[i],
            # T
            rtm[2],
            rtm[3],
            rtm[4],
            rtm[5],
        ]),
        sm:
        as_vector([
            # T
            rsm[0],
            rsm[1],
            rsm[2],
            rsm[3],
            # S
            rsm[4],
            rsm[5],
            rsm[5],
            rsm[6],
        ]),
        # This combines it all:
        w:
        as_vector([
            # S
            rw[0],
            rw[1],
            rw[1],
            rw[2],
            # T
            rw[3],
            rw[4],
            rw[5],
            rw[6],
            # Vc
            Jinv[i, 0] * as_vector([rw[7], rw[8]])[i],
            Jinv[i, 1] * as_vector([rw[7], rw[8]])[i],
            # Vd
            M_hdiv[0, j] * as_vector([rw[9], rw[10]])[j],
            M_hdiv[1, j] * as_vector([rw[9], rw[10]])[j],
            # V
            rw[11],
            rw[12],
            # U
            rw[13],
        ]),
    }

    # Check functions of various elements outside a mixed context
    check_single_function_pullback(ul2, mappings)
    check_single_function_pullback(u, mappings)
    check_single_function_pullback(v, mappings)
    check_single_function_pullback(vd, mappings)
    check_single_function_pullback(vc, mappings)
    check_single_function_pullback(t, mappings)
    check_single_function_pullback(s, mappings)

    # Check functions of various elements inside a mixed context
    check_single_function_pullback(uml2, mappings)
    check_single_function_pullback(um, mappings)
    check_single_function_pullback(vm, mappings)
    check_single_function_pullback(vdm, mappings)
    check_single_function_pullback(vcm, mappings)
    check_single_function_pullback(tm, mappings)
    check_single_function_pullback(sm, mappings)

    # Check the ridiculous mixed element W combining it all
    check_single_function_pullback(w, mappings)
예제 #17
0
def test_apply_single_function_pullbacks_triangle3d():
    triangle3d = Cell("triangle", geometric_dimension=3)
    cell = triangle3d
    domain = as_domain(cell)

    UL2 = FiniteElement("DG L2", cell, 1)
    U0 = FiniteElement("DG", cell, 0)
    U = FiniteElement("CG", cell, 1)
    V = VectorElement("CG", cell, 1)
    Vd = FiniteElement("RT", cell, 1)
    Vc = FiniteElement("N1curl", cell, 1)
    T = TensorElement("CG", cell, 1)
    S = TensorElement("CG", cell, 1, symmetry=True)
    COV2T = FiniteElement("Regge", cell, 0)  # (0, 2)-symmetric tensors
    CONTRA2T = FiniteElement("HHJ", cell, 0)  # (2, 0)-symmetric tensors

    Uml2 = UL2 * UL2
    Um = U * U
    Vm = U * V
    Vdm = V * Vd
    Vcm = Vd * Vc
    Tm = Vc * T
    Sm = T * S

    Vd0 = Vd * U0  # case from failing ffc demo

    W = S * T * Vc * Vd * V * U

    ul2 = Coefficient(UL2)
    u = Coefficient(U)
    v = Coefficient(V)
    vd = Coefficient(Vd)
    vc = Coefficient(Vc)
    t = Coefficient(T)
    s = Coefficient(S)
    cov2t = Coefficient(COV2T)
    contra2t = Coefficient(CONTRA2T)

    uml2 = Coefficient(Uml2)
    um = Coefficient(Um)
    vm = Coefficient(Vm)
    vdm = Coefficient(Vdm)
    vcm = Coefficient(Vcm)
    tm = Coefficient(Tm)
    sm = Coefficient(Sm)

    vd0m = Coefficient(Vd0)  # case from failing ffc demo

    w = Coefficient(W)

    rul2 = ReferenceValue(ul2)
    ru = ReferenceValue(u)
    rv = ReferenceValue(v)
    rvd = ReferenceValue(vd)
    rvc = ReferenceValue(vc)
    rt = ReferenceValue(t)
    rs = ReferenceValue(s)
    rcov2t = ReferenceValue(cov2t)
    rcontra2t = ReferenceValue(contra2t)

    ruml2 = ReferenceValue(uml2)
    rum = ReferenceValue(um)
    rvm = ReferenceValue(vm)
    rvdm = ReferenceValue(vdm)
    rvcm = ReferenceValue(vcm)
    rtm = ReferenceValue(tm)
    rsm = ReferenceValue(sm)

    rvd0m = ReferenceValue(vd0m)

    rw = ReferenceValue(w)
    assert len(w) == 9 + 9 + 3 + 3 + 3 + 1
    assert len(rw) == 6 + 9 + 2 + 2 + 3 + 1
    assert len(w) == 28
    assert len(rw) == 23

    assert len(vd0m) == 4
    assert len(rvd0m) == 3

    # Geometric quantities we need:
    J = Jacobian(domain)
    detJ = JacobianDeterminant(domain)
    Jinv = JacobianInverse(domain)
    # o = CellOrientation(domain)
    i, j, k, l = indices(4)

    # Contravariant H(div) Piola mapping:
    M_hdiv = ((1.0 / detJ) * J)  # Not applying cell orientation here
    # Covariant H(curl) Piola mapping: Jinv.T

    mappings = {
        # Simple elements should get a simple representation
        ul2:
        rul2 / detJ,
        u:
        ru,
        v:
        rv,
        vd:
        as_vector(M_hdiv[i, j] * rvd[j], i),
        vc:
        as_vector(Jinv[j, i] * rvc[j], i),
        t:
        rt,
        s:
        as_tensor([[rs[0], rs[1], rs[2]], [rs[1], rs[3], rs[4]],
                   [rs[2], rs[4], rs[5]]]),
        cov2t:
        as_tensor(Jinv[k, i] * rcov2t[k, l] * Jinv[l, j], (i, j)),
        contra2t:
        as_tensor(
            (1.0 / detJ) * (1.0 / detJ) * J[i, k] * rcontra2t[k, l] * J[j, l],
            (i, j)),
        # Mixed elements become a bit more complicated
        uml2:
        as_vector([ruml2[0] / detJ, ruml2[1] / detJ]),
        um:
        rum,
        vm:
        rvm,
        vdm:
        as_vector([
            # V
            rvdm[0],
            rvdm[1],
            rvdm[2],
            # Vd
            M_hdiv[0, j] * as_vector([rvdm[3], rvdm[4]])[j],
            M_hdiv[1, j] * as_vector([rvdm[3], rvdm[4]])[j],
            M_hdiv[2, j] * as_vector([rvdm[3], rvdm[4]])[j],
        ]),
        vcm:
        as_vector([
            # Vd
            M_hdiv[0, j] * as_vector([rvcm[0], rvcm[1]])[j],
            M_hdiv[1, j] * as_vector([rvcm[0], rvcm[1]])[j],
            M_hdiv[2, j] * as_vector([rvcm[0], rvcm[1]])[j],
            # Vc
            Jinv[i, 0] * as_vector([rvcm[2], rvcm[3]])[i],
            Jinv[i, 1] * as_vector([rvcm[2], rvcm[3]])[i],
            Jinv[i, 2] * as_vector([rvcm[2], rvcm[3]])[i],
        ]),
        tm:
        as_vector([
            # Vc
            Jinv[i, 0] * as_vector([rtm[0], rtm[1]])[i],
            Jinv[i, 1] * as_vector([rtm[0], rtm[1]])[i],
            Jinv[i, 2] * as_vector([rtm[0], rtm[1]])[i],
            # T
            rtm[2],
            rtm[3],
            rtm[4],
            rtm[5],
            rtm[6],
            rtm[7],
            rtm[8],
            rtm[9],
            rtm[10],
        ]),
        sm:
        as_vector([
            # T
            rsm[0],
            rsm[1],
            rsm[2],
            rsm[3],
            rsm[4],
            rsm[5],
            rsm[6],
            rsm[7],
            rsm[8],
            # S
            rsm[9],
            rsm[10],
            rsm[11],
            rsm[10],
            rsm[12],
            rsm[13],
            rsm[11],
            rsm[13],
            rsm[14],
        ]),
        # Case from failing ffc demo:
        vd0m:
        as_vector([
            M_hdiv[0, j] * as_vector([rvd0m[0], rvd0m[1]])[j],
            M_hdiv[1, j] * as_vector([rvd0m[0], rvd0m[1]])[j],
            M_hdiv[2, j] * as_vector([rvd0m[0], rvd0m[1]])[j], rvd0m[2]
        ]),
        # This combines it all:
        w:
        as_vector([
            # S
            rw[0],
            rw[1],
            rw[2],
            rw[1],
            rw[3],
            rw[4],
            rw[2],
            rw[4],
            rw[5],
            # T
            rw[6],
            rw[7],
            rw[8],
            rw[9],
            rw[10],
            rw[11],
            rw[12],
            rw[13],
            rw[14],
            # Vc
            Jinv[i, 0] * as_vector([rw[15], rw[16]])[i],
            Jinv[i, 1] * as_vector([rw[15], rw[16]])[i],
            Jinv[i, 2] * as_vector([rw[15], rw[16]])[i],
            # Vd
            M_hdiv[0, j] * as_vector([rw[17], rw[18]])[j],
            M_hdiv[1, j] * as_vector([rw[17], rw[18]])[j],
            M_hdiv[2, j] * as_vector([rw[17], rw[18]])[j],
            # V
            rw[19],
            rw[20],
            rw[21],
            # U
            rw[22],
        ]),
    }

    # Check functions of various elements outside a mixed context
    check_single_function_pullback(ul2, mappings)
    check_single_function_pullback(u, mappings)
    check_single_function_pullback(v, mappings)
    check_single_function_pullback(vd, mappings)
    check_single_function_pullback(vc, mappings)
    check_single_function_pullback(t, mappings)
    check_single_function_pullback(s, mappings)
    check_single_function_pullback(cov2t, mappings)
    check_single_function_pullback(contra2t, mappings)

    # Check functions of various elements inside a mixed context
    check_single_function_pullback(uml2, mappings)
    check_single_function_pullback(um, mappings)
    check_single_function_pullback(vm, mappings)
    check_single_function_pullback(vdm, mappings)
    check_single_function_pullback(vcm, mappings)
    check_single_function_pullback(tm, mappings)
    check_single_function_pullback(sm, mappings)

    # Check the ridiculous mixed element W combining it all
    check_single_function_pullback(w, mappings)
예제 #18
0
    def _mapped(self, t):
        # Check that we have a valid input object
        if not isinstance(t, Terminal):
            error("Expecting a Terminal.")

        # Get modifiers accumulated by previous handler calls
        ngrads = self._ngrads
        restricted = self._restricted
        avg = self._avg
        if avg != "":
            error("Averaging not implemented.")  # FIXME

        # These are the global (g) and reference (r) values
        if isinstance(t, FormArgument):
            g = t
            r = ReferenceValue(g)
        elif isinstance(t, GeometricQuantity):
            g = t
            r = g
        else:
            error("Unexpected type {0}.".format(type(t).__name__))

        # Some geometry mapping objects we may need multiple times below
        domain = t.ufl_domain()
        J = Jacobian(domain)
        detJ = JacobianDeterminant(domain)
        K = JacobianInverse(domain)

        # Restrict geometry objects if applicable
        if restricted:
            J = J(restricted)
            detJ = detJ(restricted)
            K = K(restricted)

        # Create Hdiv mapping from possibly restricted geometry objects
        Mdiv = (1.0 / detJ) * J

        # Get component indices of global and reference terminal objects
        gtsh = g.ufl_shape
        # rtsh = r.ufl_shape
        gtcomponents = compute_indices(gtsh)
        # rtcomponents = compute_indices(rtsh)

        # Create core modified terminal, with eventual
        # layers of grad applied directly to the terminal,
        # then eventual restriction applied last
        for i in range(ngrads):
            g = Grad(g)
            r = ReferenceGrad(r)
        if restricted:
            g = g(restricted)
            r = r(restricted)

        # Get component indices of global and reference objects with
        # grads applied
        gsh = g.ufl_shape
        # rsh = r.ufl_shape
        # gcomponents = compute_indices(gsh)
        # rcomponents = compute_indices(rsh)

        # Get derivative component indices
        dsh = gsh[len(gtsh):]
        dcomponents = compute_indices(dsh)

        # Create nested array to hold expressions for global
        # components mapped from reference values
        def ndarray(shape):
            if len(shape) == 0:
                return [None]
            elif len(shape) == 1:
                return [None] * shape[-1]
            else:
                return [ndarray(shape[1:]) for i in range(shape[0])]
        global_components = ndarray(gsh)

        # Compute mapping from reference values for each global component
        for gtc in gtcomponents:

            if isinstance(t, FormArgument):

                # Find basic subelement and element-local component
                # ec, element, eoffset = t.ufl_element().extract_component2(gtc) # FIXME: Translate this correctly
                eoffset = 0
                ec, element = t.ufl_element().extract_reference_component(gtc)

                # Select mapping M from element, pick row emapping =
                # M[ec,:], or emapping = [] if no mapping
                if isinstance(element, MixedElement):
                    error("Expecting a basic element here.")
                mapping = element.mapping()
                if mapping == "contravariant Piola":  # S == HDiv:
                    # Handle HDiv elements with contravariant piola
                    # mapping contravariant_hdiv_mapping = (1/det J) *
                    # J * PullbackOf(o)
                    ec, = ec
                    emapping = Mdiv[ec, :]
                elif mapping == "covariant Piola":  # S == HCurl:
                    # Handle HCurl elements with covariant piola mapping
                    # covariant_hcurl_mapping = JinvT * PullbackOf(o)
                    ec, = ec
                    emapping = K[:, ec]  # Column of K is row of K.T
                elif mapping == "identity":
                    emapping = None
                else:
                    error("Unknown mapping {0}".format(mapping))

            elif isinstance(t, GeometricQuantity):
                eoffset = 0
                emapping = None

            else:
                error("Unexpected type {0}.".format(type(t).__name__))

            # Create indices
            # if rtsh:
            #     i = Index()
            if len(dsh) != ngrads:
                error("Mismatch between derivative shape and ngrads.")
            if ngrads:
                ii = indices(ngrads)
            else:
                ii = ()

            # Apply mapping row to reference object
            if emapping:  # Mapped, always nonscalar terminal Not
                # using IndexSum for the mapping row dot product to
                # keep it simple, because we don't have a slice type
                emapped_ops = [emapping[s] * Indexed(r, MultiIndex((FixedIndex(eoffset + s),) + ii))
                               for s in range(len(emapping))]
                emapped = sum(emapped_ops[1:], emapped_ops[0])
            elif gtc:  # Nonscalar terminal, unmapped
                emapped = Indexed(r, MultiIndex((FixedIndex(eoffset),) + ii))
            elif ngrads:  # Scalar terminal, unmapped, with derivatives
                emapped = Indexed(r, MultiIndex(ii))
            else:  # Scalar terminal, unmapped, no derivatives
                emapped = r

            for di in dcomponents:
                # Multiply derivative mapping rows, parameterized by
                # free column indices
                dmapping = as_ufl(1)
                for j in range(ngrads):
                    dmapping *= K[ii[j], di[j]]  # Row of K is column of JinvT

                # Compute mapping from reference values for this
                # particular global component
                global_value = dmapping * emapped

                # Apply index sums
                # if rtsh:
                #     global_value = IndexSum(global_value, MultiIndex((i,)))
                # for j in range(ngrads): # Applied implicitly in the dmapping * emapped above
                #     global_value = IndexSum(global_value, MultiIndex((ii[j],)))

                # This is the component index into the full object
                # with grads applied
                gc = gtc + di

                # Insert in nested list
                comp = global_components
                for i in gc[:-1]:
                    comp = comp[i]
                comp[0 if gc == () else gc[-1]] = global_value

        # Wrap nested list in as_tensor unless we have a scalar
        # expression
        if gsh:
            tensor = as_tensor(global_components)
        else:
            tensor, = global_components
        return tensor