def _collapse(x): # Works by calling the matmat(), transpose() and add() methods of # diag_op/matrix_op, depending on the input types. The input is a tree # structure of block.block_mul objects, which is collapsed recursively. # This method knows too much about the internal variables of the # block_mul objects... should convert to accessor functions. from block.block_compose import block_mul, block_add, block_sub, block_transpose from block.block_mat import block_mat from block.block_util import isscalar from dolfin import GenericMatrix if isinstance(x, (matrix_op, diag_op)): return x elif isinstance(x, GenericMatrix): return matrix_op(x.down_cast().mat()) elif isinstance(x, block_mat): if x.blocks.shape != (1, 1): raise NotImplementedError( "collapse() for block_mat with shape != (1,1)") return _collapse(x[0, 0]) elif isinstance(x, block_mul): factors = map(_collapse, x) while len(factors) > 1: A = factors.pop(0) B = factors[0] if isscalar(A): C = A * B if isscalar(B) else B.matmat(A) else: C = A.matmat(B) factors[0] = C return factors[0] elif isinstance(x, block_add): A, B = map(_collapse, x) if isscalar(A) and isscalar(B): return A + B else: return B.add(A) if isscalar(A) else A.add(B) elif isinstance(x, block_sub): A, B = map(_collapse, x) if isscalar(A) and isscalar(B): return A - B else: return B.add(A, lscale=-1.0) if isscalar(A) else A.add(B, rscale=-1.0) elif isinstance(x, block_transpose): A = _collapse(x.A) return A if isscalar(A) else A._transpose() elif isscalar(x): return x else: raise NotImplementedError("collapse() for type '%s'" % str(type(x)))
def matmat(self, other): from block.block_util import isscalar try: if isscalar(other): C = self.M.copy() C.scale(other) return matrix_op(C, self.transposed) other = other.down_cast() if hasattr(other, 'mat'): if self.transposed and other.transposed: return other.transpose().matmat( self.transpose()).transpose() if self.transposed: C = self.M.transposeMatMult(other.mat()) elif other.transposed: C = self.M.matTransposeMult(other.mat()) else: C = self.M.matMult(other.mat()) return matrix_op(C) else: C = self.M.copy() if self.transposed: C.diagonalScale(L=other.vec()) else: C.diagonalScale(R=other.vec()) return matrix_op(C, self.transposed) except AttributeError: raise TypeError("can't extract matrix data from type '%s'" % str(type(other)))
def matmat(self, other): from PyTrilinos import Epetra from block.block_util import isscalar try: if isscalar(other): C = type(self.M)(self.M) if other != 1: C.Scale(other) return matrix_op(C, self.transposed) other = other.down_cast() if hasattr(other, 'mat'): from PyTrilinos import EpetraExt # Create result matrix C. This is done in a contorted way, to # ensure all diagonals are present. A = type(self.M)(self.M) A.PutScalar(1.0) B = type(other.mat())(other.mat()) B.PutScalar(1.0) C = Epetra.FECrsMatrix(Epetra.Copy, self.rowmap(), 100) EpetraExt.Multiply(A, self.transposed, B, other.transposed, C) C.OptimizeStorage() # C is now finalised, we can use it to store the real mat-mat product. assert (0 == EpetraExt.Multiply(self.M, self.transposed, other.mat(), other.transposed, C)) result = matrix_op(C) # Sanity check. Cannot trust EpetraExt.Multiply always, it seems. from block import block_vec v = result.create_vec() block_vec([v]).randomize() xv = self * other * v err = (xv - result * v).norm('l2') / (xv).norm('l2') if (err > 1e-3): print('++ a :', self * other) print('++ a\':', result) print( '++ EpetraExt.Multiply computed wrong result; ||(a-a\')x||/||ax|| = %g' % err) return result else: C = type(self.M)(self.M) if self.transposed: C.LeftScale(other.vec()) else: C.RightScale(other.vec()) return matrix_op(C, self.transposed) except AttributeError: raise TypeError("can't extract matrix data from type '%s'" % str(type(other)))
def __init__(self, A, prectype, parameters=None, pdes=1, nullspace=None): from dolfin import info from time import time T = time() Ad = A.down_cast().mat() if nullspace: from block.block_util import isscalar ns = PETSc.NullSpace() if isscalar(nullspace): ns.create(constant=True) else: ns.create(constant=False, vectors=[v.down_cast().vec() for v in nullspace]) try: Ad.setNearNullSpace(ns) except: info( 'failed to set near null space (not supported in petsc4py version)' ) self.A = A self.petsc_prec = PETSc.PC() self.petsc_prec.create() self.petsc_prec.setType(prectype) # self.petsc_prec.setOperators(Ad, Ad, PETSc.Mat.Structure.SAME_PRECONDITIONER) self.petsc_prec.setOperators(Ad, Ad) # Merge parameters into the options database if parameters: origOptions = PETSc.Options().getAll() for key, val in iter(parameters.items()): PETSc.Options().setValue(key, val) # Create preconditioner based on the options database self.petsc_prec.setFromOptions() self.petsc_prec.setUp() # Reset the options database if parameters: for key in iter(parameters.keys()): PETSc.Options().delValue(key) for key, val in iter(origOptions.items()): PETSc.Options().setValue(key, val) info('constructed %s preconditioner in %.2f s' % (self.__class__.__name__, time() - T))
def add(self, other, lscale=1.0, rscale=1.0): from block.block_util import isscalar try: if isscalar(other): x = self.v.copy() x.shift(rscale * other) return diag_op(x) other = other.down_cast() if isinstance(other, matrix_op): return other.add(self, lscale=rscale, rscale=lscale) else: x = self.v.copy() x.axpby(lscale, rscale, other.vec()) return diag_op(x) except AttributeError: raise TypeError("can't extract matrix data from type '%s'" % str(type(other)))
def add(self, other, lscale=1.0, rscale=1.0): from block.block_util import isscalar from PyTrilinos import Epetra try: if isscalar(other): x = Epetra.Vector(self.v.Map()) x.PutScalar(other) other = diag_op(x) other = other.down_cast() if isinstance(other, matrix_op): return other.add(self, lscale=rscale, rscale=lscale) else: x = Epetra.Vector(self.v) x.Update(rscale, other.vec(), lscale) return diag_op(x) except AttributeError: raise TypeError("can't extract matrix data from type '%s'"%str(type(other)))
def matmat(self, other): try: from block.block_util import isscalar if isscalar(other): x = self.v.copy() x.scale(other) return diag_op(x) other = other.down_cast() if hasattr(other, 'mat'): C = other.mat().copy() C.diagonalScale(L=self.v) return matrix_op(C) else: x = self.v.copy() x.pointwiseMult(self.v, other.vec()) return diag_op(x) except AttributeError: raise TypeError("can't extract matrix data from type '%s'" % str(type(other)))
def matmat(self, other): from PyTrilinos import Epetra try: from block.block_util import isscalar if isscalar(other): x = Epetra.Vector(self.v) x.Scale(other) return diag_op(x) other = other.down_cast() if hasattr(other, 'mat'): C = Epetra.FECrsMatrix(other.mat()) C.LeftScale(self.v) return matrix_op(C) else: x = Epetra.Vector(self.v.Map()) x.Multiply(1.0, self.v, other.vec(), 0.0) return diag_op(x) except AttributeError: raise TypeError("can't extract matrix data from type '%s'"%str(type(other)))