def cm1dvelocity(img: np.array, vel: np.array, alpha0: float, alpha1: float) -> np.array: """Computes the source for a L2-H1 mass conserving flow for a 1D image sequence and a given velocity. Takes a one-dimensional image sequence and a velocity, and returns a minimiser of the L2-H1 mass conservation functional with spatio-temporal regularisation. Args: img (np.array): A 1D image sequence of shape (m, n), where m is the number of time steps and n is the number of pixels. vel (np.array): A 1D image sequence of shape (m, n). alpha0 (float): The spatial regularisation parameter. alpha1 (float): The temporal regularisation parameter. Returns: k: A source array of shape (m, n). """ # Create mesh. [m, n] = img.shape mesh = UnitSquareMesh(m - 1, n - 1) # Define function space and functions. V = FunctionSpace(mesh, 'CG', 1) k = TrialFunction(V) w = TestFunction(V) # Convert image to function. f = Function(V) f.vector()[:] = dh.img2funvec(img) # Convert velocity to function. v = Function(V) v.vector()[:] = dh.img2funvec(vel) # Define derivatives of data. ft = Function(V) ftv = np.diff(img, axis=0) * (m - 1) ftv = np.concatenate((ftv, ftv[-1, :].reshape(1, n)), axis=0) ft.vector()[:] = dh.img2funvec(ftv) fx = Function(V) fxv = np.gradient(img, 1 / (n - 1), axis=1) fx.vector()[:] = dh.img2funvec(fxv) # Define weak formulation. A = k * w * dx + alpha0 * k.dx(1) * w.dx(1) * dx + alpha1 * k.dx(0) * w.dx( 0) * dx b = (ft + v.dx(1) * f + v * fx) * w * dx # Compute solution. k = Function(V) solve(A == b, k) # Convert back to array. k = dh.funvec2img(k.vector().get_local(), m, n) return k
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
def cm1dsource(img: np.array, k: np.array, alpha0: float, alpha1: float) -> np.array: """Computes the L2-H1 mass conserving flow for a 1D image sequence and a given source. Takes a one-dimensional image sequence and a source, and returns a minimiser of the L2-H1 mass conservation functional with spatio-temporal regularisation. Args: img (np.array): A 1D image sequence of shape (m, n), where m is the number of time steps and n is the number of pixels. k (np.array): A 1D image sequence of shape (m, n). alpha0 (float): The spatial regularisation parameter. alpha1 (float): The temporal regularisation parameter. Returns: v: A velocity array of shape (m, n). """ # Create mesh. [m, n] = img.shape mesh = UnitSquareMesh(m - 1, n - 1) # Define function space and functions. V = FunctionSpace(mesh, 'CG', 1) v = TrialFunction(V) w = TestFunction(V) # Convert image to function. f = Function(V) f.vector()[:] = dh.img2funvec(img) # Convert source to function. g = Function(V) g.vector()[:] = dh.img2funvec(k) # Define derivatives of data. ft = Function(V) ftv = np.diff(img, axis=0) * (m - 1) ftv = np.concatenate((ftv, ftv[-1, :].reshape(1, n)), axis=0) ft.vector()[:] = dh.img2funvec(ftv) fx = Function(V) fxv = np.gradient(img, 1 / (n - 1), axis=1) fx.vector()[:] = dh.img2funvec(fxv) ft = f.dx(0) fx = f.dx(1) # Define weak formulation. A = - (fx*v + f*v.dx(1)) * (fx*w + f*w.dx(1))*dx \ - alpha0*v.dx(1)*w.dx(1)*dx - alpha1*v.dx(0)*w.dx(0)*dx b = ft * (fx * w + f * w.dx(1)) * dx - g * (fx * w + f * w.dx(1)) * dx # Compute solution. v = Function(V) solve(A == b, v) # Convert back to array. vel = dh.funvec2img(v.vector().get_local(), m, n) return vel
def test_img2funvec(self): m, n = 3, 4 img = np.ones((m, n)) v = dh.img2funvec(img) np.testing.assert_allclose(v, np.ones(m * n))
def test_img2fun_fun2img(self): m, n = 7, 13 img = np.random.rand(m, n) v = dh.img2funvec(img) np.testing.assert_allclose(dh.funvec2img(v, m, n), img)
def cmscr1dnewton(img: np.array, alpha0: float, alpha1: float, alpha2: float, alpha3: float, beta: float) \ -> (np.array, np.array, float, float, bool): """Same as cmscr1d but doesn't use FEniCS Newton method. Args: img (np.array): A 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 parameter for convective regularisation. Returns: v: A velocity array of shape (m, n). k: A source array of shape (m, n). res (float): The residual. fun (float): The function value. converged (bool): True if Newton's method converged. """ # Create mesh. [m, n] = img.shape mesh = UnitSquareMesh(m - 1, n - 1) # Define function space and functions. W = VectorFunctionSpace(mesh, 'CG', 1, dim=2) w = Function(W) v, k = split(w) w1, w2 = TestFunctions(W) # Convert image to function. V = FunctionSpace(mesh, 'CG', 1) f = Function(V) f.vector()[:] = dh.img2funvec(img) # Define derivatives of data. ft = Function(V) ftv = np.diff(img, axis=0) * (m - 1) ftv = np.concatenate((ftv, ftv[-1, :].reshape(1, n)), axis=0) ft.vector()[:] = dh.img2funvec(ftv) fx = Function(V) fxv = np.gradient(img, 1 / (n - 1), axis=1) fx.vector()[:] = dh.img2funvec(fxv) # Define weak formulation. F = - (ft + fx*v + f*v.dx(1) - k) * (fx*w1 + f*w1.dx(1))*dx \ - alpha0*v.dx(1)*w1.dx(1)*dx - alpha1*v.dx(0)*w1.dx(0)*dx \ - beta*k.dx(1)*(k.dx(0) + k.dx(1)*v)*w1*dx \ + (ft + fx*v + f*v.dx(1) - k)*w2*dx \ - alpha2*k.dx(1)*w2.dx(1)*dx - alpha3*k.dx(0)*w2.dx(0)*dx \ - beta*(k.dx(0)*w2.dx(0) + k.dx(1)*v*w2.dx(0) + k.dx(0)*v*w2.dx(1) + k.dx(1)*v*v*w2.dx(1))*dx # Define derivative. dv, dk = TrialFunctions(W) J = - (fx*dv + f*dv.dx(1) - dk)*(fx*w1 + f*w1.dx(1))*dx \ - alpha0*dv.dx(1)*w1.dx(1)*dx - alpha1*dv.dx(0)*w1.dx(0)*dx \ - beta*(k.dx(1)*(dk.dx(0) + dk.dx(1)*v) + k.dx(1)**2*dv + dk.dx(1)*(k.dx(0) + k.dx(1)*v))*w1*dx \ + (fx*dv + f*dv.dx(1) - dk)*w2*dx \ - alpha2*dk.dx(1)*w2.dx(1)*dx - alpha3*dk.dx(0)*w2.dx(0)*dx \ - beta*(dv*k.dx(1) + dk.dx(0) + v*dk.dx(1))*w2.dx(0)*dx \ - beta*(dv*k.dx(0) + 2*v*dv*k.dx(1) + v*dk.dx(0) + v*v*dk.dx(1))*w2.dx(1)*dx # Alternatively, use automatic differentiation. # J = derivative(F, w) # Define algorithm parameters. tol = 1e-10 maxiter = 100 # Define increment. dw = Function(W) # Run Newton iteration. niter = 0 eps = 1 res = 1 while res > tol and niter < maxiter: niter += 1 # Solve for increment. solve(J == -F, dw) # Update solution. w.assign(w + dw) # Compute norm of increment. eps = dw.vector().norm('l2') # Compute norm of residual. a = assemble(F) res = a.norm('l2') # Print statistics. print("Iteration {0} eps={1} res={2}".format(niter, eps, res)) # Split solution. v, k = w.split(deepcopy=True) # Evaluate and print residual and functional value. res = abs(ft + fx * v + f * v.dx(1) - k) fun = 0.5 * (res**2 + alpha0 * v.dx(1)**2 + alpha1 * v.dx(0)**2 + alpha2 * k.dx(1)**2 + alpha3 * k.dx(0)**2 + beta * (k.dx(0) + k.dx(1) * v)**2) print('Res={0}, Func={1}'.format(assemble(res * dx), assemble(fun * dx))) # Convert back to array. vel = dh.funvec2img(v.vector().get_local(), m, n) k = dh.funvec2img(k.vector().get_local(), m, n) return vel, k, assemble(res * dx), assemble(fun * dx), res > tol
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