def AorthogonalityCheck(A, U, d): """ Test the frobenious norm of D^{-1}(U^TAU) - I_k """ V = np.zeros(U.shape) AV = np.zeros(U.shape) Av = Vector() v = Vector() A.init_vector(Av,0) A.init_vector(v,1) nvec = U.shape[1] for i in range(0,nvec): v.set_local(U[:,i]) v *= 1./math.sqrt(d[i]) A.mult(v,Av) AV[:,i] = Av.array() V[:,i] = v.array() VtAV = np.dot(V.T, AV) err = VtAV - np.eye(nvec, dtype=VtAV.dtype) # plt.imshow(np.abs(err)) # plt.colorbar() # plt.show() print "i, ||Vt(i,:)AV(:,i) - I_i||_F, V[:,i] = 1/sqrt(lambda_i) U[:,i]" for i in range(1,nvec+1): print i, np.linalg.norm(err[0:i,0:i], 'fro')
def to_dense(A): """ Convert a sparse matrix A to dense. For debugging only. """ if hasattr(A, "getrow"): n = A.size(0) m = A.size(1) B = np.zeros((n, m), dtype=np.float64) for i in range(0, n): [j, val] = A.getrow(i) B[i, j] = val return B else: x = Vector() Ax = Vector() A.init_vector(x, 1) A.init_vector(Ax, 0) n = Ax.array().shape[0] m = x.array().shape[0] B = np.zeros((n, m), dtype=np.float64) for i in range(0, m): i_ind = np.array([i], dtype=np.intc) x.set_local(np.ones(i_ind.shape), i_ind) A.mult(x, Ax) B[:, i] = Ax.array() x.set_local(np.zeros(i_ind.shape), i_ind) return B
def trace2(self, W=None): """ Compute the trace of A*A (Note this is the square of Frob norm, since A is symmetic). If the weight W is provided, it will compute the trace of (AW)^2. This is equivalent to tr_W(A) = \sum_i lambda_i^2, where lambda_i are the generalized eigenvalues of A x = lambda W^-1 x. Note if U is a W-orthogonal matrix then tr_W(A) = \sum_i D(i,i)^2. """ if W is None: UtU = np.dot(self.U.T, self.U) dUtU = self.d[:, None] * UtU #diag(d)*UtU. tr2 = np.sum(dUtU * dUtU) else: WU = np.zeros(self.U.shape, dtype=self.U.dtype) u, wu = Vector(), Vector() W.init_vector(u, 1) W.init_vector(wu, 0) for i in range(self.U.shape[1]): u.set_local(self.U[:, i]) W.mult(u, wu) WU[:, i] = wu.array() UtWU = np.dot(self.U.T, WU) dUtWU = self.d[:, None] * UtWU #diag(d)*UtU. tr2 = np.power(np.linalg.norm(dUtWU), 2) return tr2
def trace(self, W=None): """ Compute the trace of A. If the weight W is given compute the trace of W^1/2AW^1/2. This is equivalent to tr_W(A) = \sum_i lambda_i, where lambda_i are the generalized eigenvalues of A x = lambda W^-1 x. Note if U is a W-orthogonal matrix then tr_W(A) = \sum_i D(i,i). """ if W is None: diagUtU = np.sum(self.U * self.U, 0) tr = np.sum(self.d * diagUtU) else: WU = np.zeros(self.U.shape, dtype=self.U.dtype) u, wu = Vector(), Vector() W.init_vector(u, 1) W.init_vector(wu, 0) for i in range(self.U.shape[1]): u.set_local(self.U[:, i]) W.mult(u, wu) WU[:, i] = wu.array() diagWUtU = np.sum(WU * self.U, 0) tr = np.sum(self.d * diagWUtU) return tr
def Rn_basis(vec): if not isinstance(vec, block_vec): vec = Vector(mpi_comm_world(), vec.local_size()) values = np.zeros(vec.local_size()) for i in range(len(values)): values[i] = 1. vec.set_local(values) yield vec values[:] *= 0. else: assert False
def doublePass(A,Omega,k): """ The double pass algorithm for the HEP as presented in [1]. Inputs: - A: the operator for which we need to estimate the dominant eigenpairs. - Omega: a random gassian matrix with m >= k columns. - k: the number of eigenpairs to extract. Outputs: - d: the estimate of the k dominant eigenvalues of A - U: the estimate of the k dominant eigenvectors of A. U^T U = I_k """ w = Vector() y = Vector() A.init_vector(w,1) A.init_vector(y,0) nvec = Omega.shape[1] assert(nvec >= k ) Y = np.zeros(Omega.shape) for ivect in range(0,nvec): w.set_local(Omega[:,ivect]) A.mult(w,y) Y[:,ivect] = y.array() Q,_ = np.linalg.qr(Y) AQ = np.zeros(Omega.shape) for ivect in range(0,nvec): w.set_local(Q[:,ivect]) A.mult(w,y) AQ[:,ivect] = y.array() T = np.dot(Q.T, AQ) d, V = np.linalg.eigh(T) sort_perm = d.argsort() sort_perm = sort_perm[::-1] d = d[sort_perm[0:k]] V = V[:, sort_perm[0:k]] U = np.dot(Q,V) return d, U
def BorthogonalityTest(B, U): """ Test the frobenious norm of U^TBU - I_k """ BU = np.zeros(U.shape) Bu = Vector() u = Vector() B.init_vector(Bu,0) B.init_vector(u,1) nvec = U.shape[1] for i in range(0,nvec): u.set_local(U[:,i]) B.mult(u,Bu) BU[:,i] = Bu.array() UtBU = np.dot(U.T, BU) err = UtBU - np.eye(nvec, dtype=UtBU.dtype) print "|| UtBU - I ||_F = ", np.linalg.norm(err, 'fro')
def to_dense(A): """ Convert a sparse matrix A to dense. For debugging only. """ v = Vector() A.init_vector(v) mpi_comm = v.mpi_comm() nprocs = MPI.size(mpi_comm) if nprocs > 1: raise Exception("to_dense is only serial") if hasattr(A, "getrow"): n = A.size(0) m = A.size(1) B = np.zeros((n, m), dtype=np.float64) for i in range(0, n): [j, val] = A.getrow(i) B[i, j] = val return B else: x = Vector() Ax = Vector() A.init_vector(x, 1) A.init_vector(Ax, 0) n = Ax.array().shape[0] m = x.array().shape[0] B = np.zeros((n, m), dtype=np.float64) for i in range(0, m): i_ind = np.array([i], dtype=np.intc) x.set_local(np.ones(i_ind.shape), i_ind) x.apply("sum_values") A.mult(x, Ax) B[:, i] = Ax.array() x.set_local(np.zeros(i_ind.shape), i_ind) x.apply("sum_values") return B
def to_dense(A, mpi_comm = mpi_comm_world() ): """ Convert a sparse matrix A to dense. For debugging only. """ v = Vector(mpi_comm) A.init_vector(v) nprocs = MPI.size(mpi_comm) if nprocs > 1: raise Exception("to_dense is only serial") if hasattr(A, "getrow"): n = A.size(0) m = A.size(1) B = np.zeros( (n,m), dtype=np.float64) for i in range(0,n): [j, val] = A.getrow(i) B[i,j] = val return B else: x = Vector(mpi_comm) Ax = Vector(mpi_comm) A.init_vector(x,1) A.init_vector(Ax,0) n = Ax.get_local().shape[0] m = x.get_local().shape[0] B = np.zeros( (n,m), dtype=np.float64) for i in range(0,m): i_ind = np.array([i], dtype=np.intc) x.set_local(np.ones(i_ind.shape), i_ind) x.apply("sum_values") A.mult(x,Ax) B[:,i] = Ax.get_local() x.set_local(np.zeros(i_ind.shape), i_ind) x.apply("sum_values") return B
def gradab(self, m1, m2): self.gradm1 = project(nabla_grad(m1), self.Vd) self.gradm2 = project(nabla_grad(m2), self.Vd) uwv00, uwv00ind = [], [] uwv10, uwv10ind = [], [] uwv01, uwv01ind = [], [] uwv11, uwv11ind = [], [] for ii, x in enumerate(self.x): G = np.array([self.gradm1(x), self.gradm2(x)]).T u, s, v = np.linalg.svd(G) sqrts2eps = np.sqrt(s**2 + self.eps) W = np.diag(s / sqrts2eps) uwv = u.dot(W.dot(v)) uwv00.append(uwv[0][0]) uwv00ind.append(ii) uwv10.append(uwv[1][0]) uwv10ind.append(ii) uwv01.append(uwv[0][1]) uwv01ind.append(ii) uwv11.append(uwv[1][1]) uwv11ind.append(ii) Grad = Function(self.VV) grad = Grad.vector() rhsG = Vector() self.Gx1test.init_vector(rhsG, 1) rhsG.set_local(np.array(uwv00), np.array(uwv00ind, dtype=np.intc)) rhsG.apply('insert') grad.axpy(1.0, self.Gx1test * rhsG) rhsG.zero() rhsG.set_local(np.array(uwv10), np.array(uwv10ind, dtype=np.intc)) rhsG.apply('insert') grad.axpy(1.0, self.Gy1test * rhsG) rhsG.set_local(np.array(uwv01), np.array(uwv01ind, dtype=np.intc)) rhsG.apply('insert') grad.axpy(1.0, self.Gx2test * rhsG) rhsG.set_local(np.array(uwv11), np.array(uwv11ind, dtype=np.intc)) rhsG.apply('insert') grad.axpy(1.0, self.Gy2test * rhsG) return grad * self.k
def residualCheck(A, B, U, d): """ Test the l2 norm of the residual: r[:,i] = d[i] B U[:,i] - A U[:,i] """ u = Vector() Au = Vector() Bu = Vector() Binv_r = Vector() A.init_vector(u, 1) A.init_vector(Au, 0) B.init_vector(Bu, 0) B.init_vector(Binv_r, 0) nvec = d.shape[0] print("lambda", "||Au - lambdaBu||") for i in range(0, nvec): u.set_local(U[:, i]) A.mult(u, Au) B.mult(u, Bu) Au.axpy(-d[i], Bu) print(d[i], Au.norm("l2"))
class CGSampler: """ This class implements the CG sampler algorithm to generate samples from N(0, A^-1). REFERENCE: Albert Parker and Colin Fox Sampling Gaussian Distributions in Krylov Spaces with Conjugate Gradient SIAM J SCI COMPUT, Vol 34, No. 3 pp. B312-B334 """ def __init__(self): """ Construct the solver with default parameters tolerance = 1e-4 print_level = 0 verbose = 0 """ self.parameters = {} self.parameters["tolerance"] = 1e-4 self.parameters["print_level"] = 0 self.parameters["verbose"] = 0 self.A = None self.converged = False self.iter = 0 self.b = Vector() self.r = Vector() self.p = Vector() self.Ap = Vector() def set_operator(self, A): """ Set the operator A, such that x ~ N(0, A^-1). Note A is any object that provides the methods init_vector and mult. """ self.A = A self.A.init_vector(self.r, 0) self.A.init_vector(self.p, 0) self.A.init_vector(self.Ap, 0) self.A.init_vector(self.b, 0) self.b.set_local(np.random.randn(self.b.array().shape[0])) def sample(self, noise, s): """ Generate a sample s ~ N(0, A^-1). noise is a numpy.array of i.i.d. normal variables used as input. For a fixed realization of noise the algorithm is fully deterministic. The size of noise determine the maximum number of CG iterations. """ s.zero() self.iter = 0 self.converged = False # r0 = b self.r.zero() self.r.axpy(1., self.b) #p0 = r0 self.p.zero() self.p.axpy(1., self.r) self.A.mult(self.p, self.Ap) d = self.p.inner(self.Ap) tol2 = self.parameters["tolerance"] * self.parameters["tolerance"] rnorm2_old = self.r.inner(self.r) if self.parameters["verbose"] > 0: print "initial residual = ", math.sqrt(rnorm2_old) while (not self.converged) and (self.iter < noise.shape[0]): gamma = rnorm2_old / d s.axpy(noise[self.iter] / math.sqrt(d), self.p) self.r.axpy(-gamma, self.Ap) rnorm2 = self.r.inner(self.r) beta = rnorm2 / rnorm2_old # p_new = r + beta p self.p *= beta self.p.axpy(1., self.r) self.A.mult(self.p, self.Ap) d = self.p.inner(self.Ap) rnorm2_old = rnorm2 if rnorm2 < tol2: self.converged = True else: rnorm2_old = rnorm2 self.iter = self.iter + 1 if self.parameters["verbose"] > 0: print "Final residual {0} after {1} iterations".format( math.sqrt(rnorm2_old), self.iter)
def doublePassG(A, B, Binv, Omega,k, check_Bortho = False, check_Aortho=False, check_residual = False): """ The double pass algorithm for the GHEP as presented in [2]. B-orthogonalization is achieved using the PreCholQR algorithm. Inputs: - A: the operator for which we need to estimate the dominant generalized eigenpairs. - B: the rhs operator - Omega: a random gassian matrix with m >= k columns. - k: the number of eigenpairs to extract. Outputs: - d: the estimate of the k dominant eigenvalues of A - U: the estimate of the k dominant eigenvectors of A. U^T B U = I_k """ w = Vector() ybar = Vector() y = Vector() A.init_vector(w,1) A.init_vector(y,0) A.init_vector(ybar,0) nvec = Omega.shape[1] assert(nvec >= k ) Ybar = np.zeros(Omega.shape) Y = np.zeros(Omega.shape) for ivect in range(0,nvec): w.set_local(Omega[:,ivect]) A.mult(w,ybar) Binv.solve(y, ybar) Ybar[:,ivect] = ybar.array() Y[:,ivect] = y.array() Z,_ = np.linalg.qr(Y) BZ = np.zeros(Omega.shape) for ivect in range(0,nvec): w.set_local(Z[:,ivect]) B.mult(w,y) BZ[:, ivect] = y.array() R = np.linalg.cholesky( np.dot(Z.T,BZ )) Q = np.linalg.solve(R, Z.T).T AQ = np.zeros(Q.shape, dtype=Q.dtype) for ivect in range(0,nvec): w.set_local(Q[:,ivect]) A.mult(w,ybar) AQ[:,ivect] = ybar.array() T = np.dot(Q.T, AQ) d, V = np.linalg.eigh(T) sort_perm = d.argsort() sort_perm = sort_perm[::-1] d = d[sort_perm[0:k]] V = V[:, sort_perm[0:k]] U = np.dot(Q,V) if check_Bortho: BorthogonalityTest(B, U) if check_Aortho: AorthogonalityCheck(A, U, d) if check_residual: residualCheck(A,B, U, d) return d, U
class TraceEstimator: """ An unbiased stochastic estimator for the trace of A. d = \sum_{j=1}^k (vj, A vj) where - vj are i.i.d. Rademacher or Gaussian random vectors - (.,.) represents the inner product The number of samples k is estimated at run time based on the variance of the estimator. REFERENCE: Haim Avron and Sivan Toledo, Randomized algorithms for estimating the trace of an implicit symmetric positive semi-definite matrix, Journal of the ACM (JACM), 58 (2011), p. 17. """ def __init__(self, A, solve_mode=False, accurancy=1e-1, init_vector=None, random_engine=rademacher_engine): """ Constructor: - A: an operator - solve_mode: if True we estimate the trace of A^{-1}, otherwise of A. - accurancy: we stop when the standard deviation of the estimator is less then accurancy*tr(A). - init_vector: use a custom function to initialize a vector compatible with the range/domain of A - random_engine: which type of i.i.d. random variables to use (Rademacher or Gaussian) """ self.A = A self.accurancy = accurancy self.random_engine = random_engine self.iter = 0 self.z = Vector() self.Az = Vector() if solve_mode: self._apply = self._apply_solve else: self._apply = self._apply_mult if init_vector is None: A.init_vector(self.z, 0) A.init_vector(self.Az, 0) else: init_vector(self.z, 0) init_vector(self.Az, 0) def _apply_mult(self, z, Az): self.A.mult(z, Az) def _apply_solve(self, z, Az): self.A.solve(Az, z) def __call__(self, min_iter=5, max_iter=100): """ Estimate the trace of A (or A^-1) using at least min_iter and at most max_iter samples. """ sum_tr = 0 sum_tr2 = 0 self.iter = 0 size = len(self.z.array()) while self.iter < min_iter: self.iter += 1 self.z.set_local(self.random_engine(size)) self._apply(self.z, self.Az) tr = self.z.inner(self.Az) sum_tr += tr sum_tr2 += tr * tr exp_tr = sum_tr / float(self.iter) exp_tr2 = sum_tr2 / float(self.iter) var_tr = exp_tr2 - exp_tr * exp_tr # print(exp_tr, math.sqrt( var_tr ), self.accurancy*exp_tr) self.converged = True while (math.sqrt(var_tr) > self.accurancy * exp_tr): self.iter += 1 self.z.set_local(self.random_engine(size)) self._apply(self.z, self.Az) tr = self.z.inner(self.Az) sum_tr += tr sum_tr2 += tr * tr exp_tr = sum_tr / float(self.iter) exp_tr2 = sum_tr2 / float(self.iter) var_tr = exp_tr2 - exp_tr * exp_tr # print(exp_tr, math.sqrt( var_tr ), self.accurancy*exp_tr) if (self.iter > max_iter): self.converged = False break return exp_tr, var_tr