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 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 estimate_diagonal_inv2(Asolver, k, d): """ An unbiased stochastic estimator for the diagonal of A^-1. d = [ \sum_{j=1}^k vj .* A^{-1} vj ] ./ [ \sum_{j=1}^k vj .* vj ] where - vj are i.i.d. ~ N(0, I) - .* and ./ represent the element-wise multiplication and division of vectors, respectively. REFERENCE: Costas Bekas, Effrosyni Kokiopoulou, and Yousef Saad, An estimator for the diagonal of a matrix, Applied Numerical Mathematics, 57 (2007), pp. 1214-1229. """ x, b = Vector(), Vector() if hasattr(Asolver, "init_vector"): Asolver.init_vector(x, 1) Asolver.init_vector(b, 0) else: Asolver.get_operator().init_vector(x, 1) Asolver.get_operator().init_vector(b, 0) num = np.zeros(b.array().shape, dtype=b.array().dtype) den = np.zeros(num.shape, dtype=num.dtype) for i in range(k): x.zero() b.set_local(np.random.randn(num.shape[0])) Asolver.solve(x, b) num = num + (x.array() * b.array()) den = den + (b.array() * b.array()) d.set_local(num / den)
def get_diagonal(A, d, solve_mode=True): """ Compute the diagonal of the square operator A or its inverse A^{-1} (if solve_mode=True). """ ej, xj = Vector(), Vector() if hasattr(A, "init_vector"): A.init_vector(ej, 1) A.init_vector(xj, 0) else: A.get_operator().init_vector(ej, 1) A.get_operator().init_vector(xj, 0) ncol = ej.size() da = np.zeros(ncol, dtype=ej.array().dtype) for j in range(ncol): ej[j] = 1. if solve_mode: A.solve(xj, ej) else: A.mult(ej, xj) da[j] = xj[j] ej[j] = 0. d.set_local(da)
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 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')
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