示例#1
0
def test_facet_areas(dim):
    """
    Check that the computation of facet areas
    sums to the expected value for a uniform
    (isotropic) triangular or tetrahedral mesh.
    """
    mesh = uniform_mesh(dim, 1)
    fa = get_facet_areas(mesh)
    expected = 5.414214 if dim == 2 else 10.242641
    assert np.isclose(sum(fa.dat.data), expected)
示例#2
0
def test_hessian_metric(dim):
    """
    Check that `hessian_metric` does not have
    an effect on a matrix that is already SPD.
    """
    mesh = uniform_mesh(dim, 4)
    P1_ten = TensorFunctionSpace(mesh, "CG", 1)
    M = interpolate(Identity(dim), P1_ten)
    M -= hessian_metric(M)
    assert np.isclose(norm(M), 0.0)
示例#3
0
def test_is_spd(dim):
    """
    Check `is_spd` correctly identifies an SPD
    and a non-SPD matrix.
    """
    mesh = uniform_mesh(dim, 4)
    P1_ten = TensorFunctionSpace(mesh, "CG", 1)
    M = interpolate(Identity(dim), P1_ten)
    assert is_spd(M)
    M.sub(0).assign(-1)
    assert not is_spd(M)
示例#4
0
def test_consistency_3d(measure):
    """
    Check that the C++ and Python implementations of the
    quality measures are consistent for a non-uniform
    3D tetrahedral mesh.
    """
    np.random.seed(0)
    measure = getattr(mq, measure)
    mesh = uniform_mesh(3, 4)
    mesh.coordinates.dat.data[:] += np.random.rand(
        *mesh.coordinates.dat.data.shape)
    quality_cpp = measure(mesh)
    quality_py = measure(mesh, python=True)
    assert np.allclose(quality_cpp.dat.data, quality_py.dat.data)
示例#5
0
def test_uniform_quality_3d(measure, expected):
    """
    Check that the computation of each quality measure
    gives the expected value for a uniform (isotropic)
    3D tetrahedral mesh.
    """
    measure = getattr(mq, measure)
    mesh = uniform_mesh(3, 4)
    if measure.__name__ == "get_quality_metrics3d":
        P1_ten = TensorFunctionSpace(mesh, "CG", 1)
        M = interpolate(Identity(3), P1_ten)
        q = measure(mesh, M)
    else:
        q = measure(mesh)
    true_vals = np.array([expected for _ in q.dat.data])
    assert np.allclose(true_vals, q.dat.data)
    if measure.__name__ == "get_volumes3d":
        assert np.isclose(sum(q.dat.data), 1)
示例#6
0
def test_intersection(dim):
    """
    Check that metric intersection DTRT when
    applied to two isotropic metrics.
    """
    mesh = uniform_mesh(dim, 3)
    P1_ten = TensorFunctionSpace(mesh, "CG", 1)
    I = Identity(dim)
    M1 = interpolate(2 * I, P1_ten)
    M2 = interpolate(1 * I, P1_ten)
    M = metric_intersection(M1, M2)
    assert np.isclose(norm(Function(M).assign(M - M1)), 0.0)
    M2.interpolate(2 * I)
    M = metric_intersection(M1, M2)
    assert np.isclose(norm(Function(M).assign(M - M1)), 0.0)
    assert np.isclose(norm(Function(M).assign(M - M2)), 0.0)
    M2.interpolate(4 * I)
    M = metric_intersection(M1, M2)
    assert np.isclose(norm(Function(M).assign(M - M2)), 0.0)
示例#7
0
def test_metric_math(dim):
    """
    Check that the metric exponential and
    metric logarithm are indeed inverses.
    """
    mesh = uniform_mesh(dim, 1)
    P0_ten = firedrake.TensorFunctionSpace(mesh, "DG", 0)
    I = ufl.Identity(dim)
    M = firedrake.interpolate(2 * I, P0_ten)
    logM = metric_logarithm(M)
    expected = firedrake.interpolate(np.log(2) * I, P0_ten)
    assert np.allclose(logM.dat.data, expected.dat.data)
    M_ = metric_exponential(logM)
    assert np.allclose(M.dat.data, M_.dat.data)
    expM = metric_exponential(M)
    expected = firedrake.interpolate(np.exp(2) * I, P0_ten)
    assert np.allclose(expM.dat.data, expected.dat.data)
    M_ = metric_logarithm(expM)
    assert np.allclose(M.dat.data, M_.dat.data)
示例#8
0
def test_eigendecomposition(dim, reorder):
    """
    Check decomposition of a metric into its eigenvectors
    and eigenvalues.

      * The eigenvectors should be orthonormal.
      * Applying `compute_eigendecomposition` followed by
        `set_eigendecomposition` should get back the metric.
    """
    mesh = uniform_mesh(dim, 20)

    # Recover Hessian metric for some arbitrary sensor
    f = np.prod([sin(pi * xi) for xi in SpatialCoordinate(mesh)])
    metric = hessian_metric(recover_hessian(f, mesh=mesh))
    P1_ten = metric.function_space()

    # Extract the eigendecomposition
    evectors, evalues = compute_eigendecomposition(metric, reorder=reorder)

    # Check eigenvectors are orthonormal
    err = Function(P1_ten)
    err.interpolate(dot(evectors, transpose(evectors)) - Identity(dim))
    if not np.isclose(norm(err), 0.0):
        raise ValueError("Eigenvectors are not orthonormal")

    # Check eigenvalues are in descending order
    if reorder:
        P1 = FunctionSpace(mesh, "CG", 1)
        for i in range(dim - 1):
            f = interpolate(evalues[i], P1)
            f -= interpolate(evalues[i + 1], P1)
            if f.vector().gather().min() < 0.0:
                raise ValueError("Eigenvalues are not in descending order")

    # Reassemble it and check the two match
    metric -= assemble_eigendecomposition(evectors, evalues)
    if not np.isclose(norm(metric), 0.0):
        raise ValueError("Reassembled metric does not match")
示例#9
0
def test_density_quotients_decomposition(dim, reorder):
    """
    Check decomposition of a metric into its density
    and anisotropy quotients.

    Reassembling should get back the metric.
    """
    mesh = uniform_mesh(dim, 20)

    # Recover Hessian metric for some arbitrary sensor
    f = np.prod([sin(pi * xi) for xi in SpatialCoordinate(mesh)])
    metric = hessian_metric(recover_hessian(f, mesh=mesh))

    # Extract the density, anisotropy quotients and eigenvectors
    density, quotients, evectors = density_and_quotients(metric,
                                                         reorder=reorder)
    quotients.interpolate(
        as_vector([pow(density / Q, 2 / dim) for Q in quotients]))

    # Reassemble the matrix and check the two match
    metric -= assemble_eigendecomposition(evectors, quotients)
    if not np.isclose(norm(metric), 0.0):
        raise ValueError("Reassembled metric does not match")
示例#10
0
def mesh_for_sensors(dim, n):
    mesh = uniform_mesh(dim, n, l=2)
    coords = firedrake.Function(mesh.coordinates)
    coords -= 1.0
    return firedrake.Mesh(coords)