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))
    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
示例#3
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
示例#4
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