Exemple #1
0
    def test_cm1d_weak_solution_default(self):
        # Define temporal and spatial sample points.
        m, n = 10, 20

        # Define mesh and function space.
        mesh = UnitSquareMesh(m - 1, n - 1)
        V = dh.create_function_space(mesh, 'default')

        # Create zero function.
        f = Function(V)

        # Compute velocity.
        v, res, fun = cm1d_weak_solution(V, f, f.dx(0), f.dx(1), 1.0, 1.0)
        v = v.vector().get_local()

        np.testing.assert_allclose(v.shape, m * n)
        np.testing.assert_allclose(v, np.zeros_like(v))
        np.testing.assert_allclose(res, 0)
        np.testing.assert_allclose(fun, 0)

        V = dh.create_function_space(mesh, 'periodic')

        # Create zero function.
        f = Function(V)

        # Compute velocity.
        v, res, fun = cm1d_weak_solution(V, f, f.dx(0), f.dx(1), 1.0, 1.0)
        v = v.vector().get_local()

        np.testing.assert_allclose(v.shape, m * (n - 1))
        np.testing.assert_allclose(v, np.zeros_like(v))
        np.testing.assert_allclose(res, 0)
        np.testing.assert_allclose(fun, 0)
Exemple #2
0
    def test_cms1d_weak_solution_given_velocity_default(self):
        # Define temporal and spatial sample points.
        m, n = 10, 20

        # Define mesh and function space.
        mesh = UnitSquareMesh(m - 1, n - 1)
        V = dh.create_function_space(mesh, 'default')

        # Create zero function.
        f = Function(V)

        # Create velocity.
        v = Function(V)
        vx = Function(V)

        # Compute source.
        k, res, fun = cms1d_weak_solution_given_velocity(
            V, f, f.dx(0), f.dx(1), v, vx, 1.0, 1.0)
        k = k.vector().get_local()

        np.testing.assert_allclose(k.shape, m * n)
        np.testing.assert_allclose(k, np.zeros_like(k))

        V = dh.create_function_space(mesh, 'periodic')

        # Create zero function.
        f = Function(V)

        # Compute source.
        k, res, fun = cms1d_weak_solution_given_velocity(
            V, f, f.dx(0), f.dx(1), v, vx, 1.0, 1.0)
        k = k.vector().get_local()

        np.testing.assert_allclose(k.shape, m * (n - 1))
        np.testing.assert_allclose(k, np.zeros_like(k))
Exemple #3
0
def cm1d_exp(m: int, n: int, f: Expression, ft: Expression, fx: Expression,
             alpha0: float, alpha1: float) -> np.array:
    """Computes the L2-H1 mass conserving flow for a 1D image sequence.

    Args:
        m (int): Number of temporal sampling points.
        n (int): Number of spatial sampling points.
        f (Expression): 1D image sequence.
        ft (Expression): Partial derivative of f wrt. time.
        fx (Expression): Partial derivative of f wrt. space.
        alpha0 (float): Spatial regularisation parameter.
        alpha1 (float): Temporal regularisation parameter.

    Returns:
        v (np.array): A velocity array of shape (m, n).
        res (float): The residual.
        func (float): The value of the functional.

    """
    # Define mesh and function space.
    mesh = UnitSquareMesh(m - 1, n - 1)
    V = dh.create_function_space(mesh, 'default')

    # Compute velocity.
    v, res, fun = cm1d_weak_solution(V, f, ft, fx, alpha0, alpha1)

    # Convert to array and return.
    return dh.funvec2img(v.vector().get_local(), m, n), res, fun
Exemple #4
0
def of1d_exp(m: int, n: int,
             f: Expression, ft: Expression, fx: Expression,
             alpha0: float, alpha1) -> (np.array, float, float):
    """Computes the L2-H1 optical flow for a 1D image sequence.

    Takes a one-dimensional image sequence and partial derivatives, and returns
    a minimiser of the Horn-Schunck functional with spatio-temporal
    regularisation.

    Args:
        m (int): Number of temporal sampling points.
        n (int): Number of spatial sampling points.
        f (Expression): 1D image sequence.
        ft (Expression): Partial derivative of f wrt. time.
        fx (Expression): Partial derivative of f wrt. space.
        alpha0 (float): Spatial regularisation parameter.
        alpha1 (float): Temporal regularisation parameter.

    Returns:
        v (np.array): A velocity array of shape (m, n).
        res (float): The residual.
        fun (float): The function value.
    """
    # Define mesh and function space.
    mesh = UnitSquareMesh(m - 1, n - 1)
    V = dh.create_function_space(mesh, 'default')

    # Compute velocity.
    v, res, fun = of1d_weak_solution(V, f, ft, fx, alpha0, alpha1)

    # Convert to array and return.
    return dh.funvec2img(v.vector().get_local(), m, n), res, fun
Exemple #5
0
    def test_cmscr1d_weak_solution_zero_Dirichlet(self):
        print("Running test 'test_cmscr1d_weak_solution_zero_Dirichlet'")
        # Define temporal and spatial sample points.
        m, n = 10, 20

        # Define mesh and function space.
        mesh = UnitSquareMesh(m - 1, n - 1)
        V = dh.create_function_space(mesh, 'default')
        W = dh.create_vector_function_space(mesh, 'default')

        # Define boundary conditions for velocity.
        bc = DirichletBC(W.sub(0), Constant(0), dh.DirichletBoundary())

        # Create zero function.
        f = Function(V)

        # Compute velocity.
        v, k, res, fun, converged = cmscr1d_weak_solution(W, f,
                                                          f.dx(0), f.dx(1),
                                                          1.0, 1.0,
                                                          1.0, 1.0, 1.0,
                                                          bcs=bc)
        v = v.vector().get_local()
        k = k.vector().get_local()

        np.testing.assert_allclose(v.shape, m * n)
        np.testing.assert_allclose(v, np.zeros_like(v))
        np.testing.assert_allclose(k.shape, m * n)
        np.testing.assert_allclose(k, np.zeros_like(k))
Exemple #6
0
def cmscr1d_img_pb(img: np.array,
                   alpha0: float, alpha1: float,
                   alpha2: float, alpha3: float,
                   beta: float, deriv='mesh') \
                   -> (np.array, np.array, float, float, bool):
    """Computes the L2-H1 mass conserving flow with source for a 1D image
    sequence with spatio-temporal and convective regularisation with periodic
    spatial boundary.

    Args:
        img (np.array): 1D image sequence of shape (m, n), where m is the
                        number of time steps and n is the number of pixels.
        alpha0 (float): The spatial regularisation parameter for v.
        alpha1 (float): The temporal regularisation parameter for v.
        alpha2 (float): The spatial regularisation parameter for k.
        alpha3 (float): The temporal regularisation parameter for k.
        beta (float): The convective regularisation parameter.
        deriv (str): Specifies how to approximate pertial derivatives.
                     When set to 'mesh' it uses FEniCS built in function.

    Returns:
        v (np.array): A velocity array of shape (m, n).
        k (np.array): A source array of shape (m, n).
        res (float): The residual.
        fun (float): The function value.
        converged (bool): True if Newton's method converged.

    """
    # Check for valid arguments.
    valid = {'mesh'}
    if deriv not in valid:
        raise ValueError("Argument 'deriv' must be one of %r." % valid)

    # Create mesh.
    m, n = img.shape
    mesh = UnitSquareMesh(m - 1, n - 1)

    # Define function space.
    V = dh.create_function_space(mesh, 'periodic')
    W = dh.create_vector_function_space(mesh, 'periodic')

    # Convert array to function.
    f = Function(V)
    f.vector()[:] = dh.img2funvec_pb(img)

    # Compute partial derivatives.
    ft, fx = f.dx(0), f.dx(1)

    # Compute velocity.
    v, k, res, fun, converged = cmscr1d_weak_solution(W, f, ft, fx, alpha0,
                                                      alpha1, alpha2, alpha3,
                                                      beta)

    # Convert back to array and return.
    v = dh.funvec2img(v.vector().get_local(), m, n)
    k = dh.funvec2img(k.vector().get_local(), m, n)
    return v, k, res, fun, converged
Exemple #7
0
    def test_cmscr1d_weak_solution_default(self):
        print("Running test 'test_cmscr1d_weak_solution_default'")
        # Define temporal and spatial sample points.
        m, n = 10, 20

        # Define mesh and function space.
        mesh = UnitSquareMesh(m - 1, n - 1)
        V = dh.create_function_space(mesh, 'default')
        W = dh.create_vector_function_space(mesh, 'default')

        # Create zero function.
        f = Function(V)

        # Compute velocity.
        v, k, res, fun, converged = cmscr1d_weak_solution(W, f,
                                                          f.dx(0), f.dx(1),
                                                          1.0, 1.0,
                                                          1.0, 1.0, 1.0)
        v = v.vector().get_local()
        k = k.vector().get_local()

        np.testing.assert_allclose(v.shape, m * n)
        np.testing.assert_allclose(v, np.zeros_like(v))
        np.testing.assert_allclose(k.shape, m * n)
        np.testing.assert_allclose(k, np.zeros_like(k))

        V = dh.create_function_space(mesh, 'periodic')

        # Create zero function.
        f = Function(V)

        # Compute velocity.
        v, k, res, fun, converged = cmscr1d_weak_solution(W, f,
                                                          f.dx(0), f.dx(1),
                                                          1.0, 1.0, 1.0,
                                                          1.0, 1.0)
        v = v.vector().get_local()
        k = k.vector().get_local()

        np.testing.assert_allclose(v.shape, m * n)
        np.testing.assert_allclose(v, np.zeros_like(v))
        np.testing.assert_allclose(k.shape, m * n)
        np.testing.assert_allclose(k, np.zeros_like(k))
Exemple #8
0
def of1d_img(img: np.array, alpha0: float, alpha1: float, deriv) \
             -> (np.array, float, float):
    """Computes the L2-H1 optical flow for a 1D image sequence.

    Takes a one-dimensional image sequence and returns a minimiser of the
    Horn-Schunck functional with spatio-temporal regularisation.

    Allows to specify how to approximate partial derivatives of f numerically.

    Args:
        img (np.array): 1D image sequence of shape (m, n), where m is the
                        number of time steps and n is the number of pixels.
        alpha0 (float): Spatial regularisation parameter.
        alpha1 (float): Temporal regularisation parameter.
        deriv (str): Specifies how to approximate pertial derivatives.
                     When set to 'mesh' it uses FEniCS built in function.
                     When set to 'fd' it uses finite differences.

    Returns:
        v (np.array): A velocity array of shape (m, n).
        res (float): The residual.
        fun (float): The function value.
    """
    # Check for valid arguments.
    valid = {'mesh', 'fd'}
    if deriv not in valid:
        raise ValueError("Argument 'deriv' must be one of %r." % valid)

    # Create mesh.
    m, n = img.shape
    mesh = UnitSquareMesh(m - 1, n - 1)

    # Define function space.
    V = dh.create_function_space(mesh, 'default')

    # Convert array to function.
    f = Function(V)
    f.vector()[:] = dh.img2funvec(img)

    # Compute partial derivatives.
    if deriv is 'mesh':
        ft, fx = f.dx(0), f.dx(1)
    if deriv is 'fd':
        imgt, imgx = nh.partial_derivatives(img)
        ft, fx = Function(V), Function(V)
        ft.vector()[:] = dh.img2funvec(imgt)
        fx.vector()[:] = dh.img2funvec(imgx)

    # Compute velocity.
    v, res, fun = of1d_weak_solution(V, f, ft, fx, alpha0, alpha1)

    # Convert to array and return.
    return dh.funvec2img(v.vector().get_local(), m, n), res, fun
Exemple #9
0
def mpi1d_exp_pb(m: int, n: int, f: Expression, ft: Expression,
                 fx: Expression) -> (np.array, np.array):

    # Define mesh and function space.
    mesh = UnitSquareMesh(m - 1, n - 1)
    V = dh.create_function_space(mesh, 'periodic')

    # Compute velocity.
    v, k = mpi1d_weak_solution(V, f, ft, fx)

    # Convert back to array and return.
    v = dh.funvec2img_pb(v.vector().get_local(), m, n)
    k = dh.funvec2img_pb(k.vector().get_local(), m, n)
    return v, k
Exemple #10
0
def cm1d_img_pb(img: np.array,
                alpha0: float,
                alpha1: float,
                deriv='mesh') -> np.array:
    """Computes the L2-H1 mass conserving flow for a 1D image sequence with
    periodic spatial boundary.

    Allows to specify how to approximate partial derivatives of f numerically.

    Note that the last column of img is ignored.

    Args:
        img (np.array): 1D image sequence of shape (m, n), where m is the
                        number of time steps and n is the number of pixels.
        alpha0 (float): Spatial regularisation parameter.
        alpha1 (float): Temporal regularisation parameter.
        deriv (str): Specifies how to approximate pertial derivatives.
                     When set to 'mesh' it uses FEniCS built in function.

    Returns:
        v (np.array): A velocity array of shape (m, n).
        res (float): The residual.
        func (float): The value of the functional.

    """
    # Check for valid arguments.
    valid = {'mesh'}
    if deriv not in valid:
        raise ValueError("Argument 'deriv' must be one of %r." % valid)

    # Create mesh.
    m, n = img.shape
    mesh = UnitSquareMesh(m - 1, n - 1)

    # Define function space.
    V = dh.create_function_space(mesh, 'periodic')

    # Convert array to function.
    f = Function(V)
    f.vector()[:] = dh.img2funvec_pb(img)

    # Compute partial derivatives.
    ft, fx = f.dx(0), f.dx(1)

    # Compute velocity.
    v, res, fun = cm1d_weak_solution(V, f, ft, fx, alpha0, alpha1)

    # Convert to array and return.
    return dh.funvec2img_pb(v.vector().get_local(), m, n), res, fun
Exemple #11
0
alpha1 = 1e-3  # v_t
alpha2 = 1e-4  # k_x
alpha3 = 1e-4  # k_t
beta = 1e-3  # D_v k
gamma = 1e-4  # k

# Set parameters of data.
w = 0.1
lambdap = 1 / (4 * np.pi)
tau = 1.0
c0 = 0.0

# Create mesh and function spaces.
m, n = 40, 100
mesh = UnitSquareMesh(m - 1, n - 1)
V = dh.create_function_space(mesh, 'default')
W = dh.create_function_space(mesh, 'periodic')

# Run experiments with decaying data.
f = f_decay(degree=2)
ft = f_decay_t(degree=1)
fx = f_decay_x(degree=1)
datastr = DecayingData().string()

# Interpolate function.
fa = interpolate(f, V)
fa = dh.funvec2img(fa.vector().get_local(), m, n)

fa_pb = interpolate(f, W)
fa_pb = dh.funvec2img_pb(fa_pb.vector().get_local(), m, n)
Exemple #12
0
def cmscr1d_img(img: np.array,
                alpha0: float, alpha1: float,
                alpha2: float, alpha3: float,
                beta: float, deriv, mesh=None, bc='natural') \
                -> (np.array, np.array, float, float, bool):
    """Computes the L2-H1 mass conserving flow with source for a 1D image
    sequence with spatio-temporal and convective regularisation.

    Allows to specify how to approximate partial derivatives of f numerically.

    Args:
        img (np.array): 1D image sequence of shape (m, n), where m is the
                        number of time steps and n is the number of pixels.
        alpha0 (float): The spatial regularisation parameter for v.
        alpha1 (float): The temporal regularisation parameter for v.
        alpha2 (float): The spatial regularisation parameter for k.
        alpha3 (float): The temporal regularisation parameter for k.
        beta (float): The convective regularisation parameter.
        deriv (str): Specifies how to approximate pertial derivatives.
                     When set to 'mesh' it uses FEniCS built in function.
                     When set to 'fd' it uses finite differences.
        mesh: A custom mesh (optional). Must have (m - 1, n - 1) cells.
        bc (str): One of {'natural', 'zero', 'zerospace'} for boundary
                    conditions for the velocity v (optional).

    Returns:
        v (np.array): A velocity array of shape (m, n).
        k (np.array): A source array of shape (m, n).
        res (float): The residual.
        fun (float): The function value.
        converged (bool): True if Newton's method converged.

    """
    # Check for valid arguments.
    valid = {'mesh', 'fd'}
    if deriv not in valid:
        raise ValueError("Argument 'deriv' must be one of %r." % valid)

    # Create mesh.
    m, n = img.shape
    if mesh is None:
        mesh = UnitSquareMesh(m - 1, n - 1)

    # Define function spaces.
    V = dh.create_function_space(mesh, 'default')
    W = dh.create_vector_function_space(mesh, 'default')

    # Convert array to function.
    f = Function(V)
    f.vector()[:] = dh.img2funvec(img)

    # Compute partial derivatives.
    if deriv is 'mesh':
        ft, fx = f.dx(0), f.dx(1)
    if deriv is 'fd':
        imgt, imgx = nh.partial_derivatives(img)
        ft, fx = Function(V), Function(V)
        ft.vector()[:] = dh.img2funvec(imgt)
        fx.vector()[:] = dh.img2funvec(imgx)

    # Check for valid arguments.
    valid = {'natural', 'zero', 'zerospace'}
    if bc not in valid:
        raise ValueError("Argument 'bc' must be one of %r." % valid)

    # Define boundary conditions for velocity.
    if bc is 'natural':
        bc = []
    if bc is 'zero':
        bc = DirichletBC(W.sub(0), Constant(0), dh.DirichletBoundary())
    if bc is 'zerospace':
        bc = DirichletBC(W.sub(0), Constant(0), dh.DirichletBoundarySpace())

    # Compute velocity.
    v, k, res, fun, converged = cmscr1d_weak_solution(W,
                                                      f,
                                                      ft,
                                                      fx,
                                                      alpha0,
                                                      alpha1,
                                                      alpha2,
                                                      alpha3,
                                                      beta,
                                                      bcs=bc)

    # Convert back to array and return.
    v = dh.funvec2img(v.vector().get_local(), m, n)
    k = dh.funvec2img(k.vector().get_local(), m, n)
    return v, k, res, fun, converged