def block_collapse(self): """Create a block_mat of block_adds from a block_add of block_mats. See block_transform.block_collapse.""" from block_mat import block_mat from block_util import isscalar from block_transform import block_collapse, block_simplify A = block_collapse(self.A) B = block_collapse(self.B) # The self.__class__(A,B) used below works for both block_sub and # block_add, and any scalar terms are combined by the final call to # block_simplify(). if isinstance(A, block_mat) and isinstance(B, block_mat): m, n = A.blocks.shape C = block_mat(m, n) for row in range(m): for col in range(n): C[row, col] = self.__class__(A[row, col], B[row, col]) elif isinstance(A, block_mat) and isscalar(B): m, n = A.blocks.shape C = block_mat(m, n) for row in range(m): for col in range(n): C[row, col] = self.__class__( A[row, col], B) if row == col else A[row, col] elif isinstance(B, block_mat) and isscalar(A): m, n = B.blocks.shape C = block_mat(m, n) for row in range(m): for col in range(n): C[row, col] = self.__class__( A, B[row, col]) if row == col else B[row, col] else: C = self.__class__(A, B) return block_simplify(C)
def block_simplify(self): """Try to combine scalar terms and remove additive identities, recursively. A fuller explanation is found in block_transform.block_simplify. """ from block_util import isscalar from block_transform import block_simplify A = block_simplify(self.A) B = block_simplify(self.B) if isscalar(A) and A == 0: return B if isscalar(B) and B == 0: return A return A + B
def block_simplify(self): """Try to convert identities to scalars, recursively. A fuller explanation is found in block_transform.block_simplify. """ from block_util import isscalar from block_transform import block_simplify m, n = self.blocks.shape res = block_mat(m, n) # Recursive call for i in range(m): for j in range(n): res[i, j] = block_simplify(self[i, j]) # Check if result after recursive conversion is the (scaled) identity v0 = res.blocks[0, 0] if m != n: return res for i in range(m): for j in range(n): block = res.blocks[i, j] if not isscalar(block): return res if i == j: if block != v0: return res else: if block != 0: return res return v0
def transpmult(self, x): from block_util import isscalar for op in self.chain: if isscalar(op): if op != 1: x = op * x else: x = op.transpmult(x) return x
def block_collapse(self): """Create a block_mat of block_muls from a block_mul of block_mats. See block_transform.block_collapse.""" from block_mat import block_mat from block_util import isscalar from block_transform import block_collapse, block_simplify # Reduce all composed objects ops = map(block_collapse, self.chain) # Do the matrix multiply, blockwise. Note that we use # block_mul(A,B) rather than A*B to avoid any implicit calculations # (e.g., scalar*matrix->matrix) -- the result will be transformed by # block_simplify() in the end to take care of any stray scalars. while len(ops) > 1: B = ops.pop() A = ops.pop() if isinstance(A, block_mat) and isinstance(B, block_mat): m, n = A.blocks.shape p, q = B.blocks.shape C = block_mat(m, q) for row in range(m): for col in range(q): for i in range(n): C[row, col] += block_mul(A[row, i], B[i, col]) elif isinstance(A, block_mat) and isscalar(B): m, n = A.blocks.shape C = block_mat(m, n) for row in range(m): for col in range(n): C[row, col] = block_mul(A[row, col], B) elif isinstance(B, block_mat) and isscalar(A): m, n = B.blocks.shape C = block_mat(m, n) for row in range(m): for col in range(n): C[row, col] = block_mul(A, B[row, col]) else: C = block_mul(A, B) ops.append(C) return block_simplify(ops[0])
def matvec(self, x): from dolfin import GenericVector, GenericMatrix from block_util import isscalar import numpy m, n = self.blocks.shape y = block_vec(m) for i in range(m): for j in range(n): if isinstance(self[i, j], (numpy.ndarray, numpy.matrix)): z = numpy.matrix(self[i, j]) * numpy.matrix( x[j].array()).transpose() z = numpy.array(z).flatten() if y[i] is None: y[i] = x[j].copy() y[i].resize(len(z)) y[i][:] += z[:] continue if self[i,j] is None or self[i,j]==0 \ or x[j] is None or (isscalar(x[j]) and x[j]==0): # Skip multiply if zero continue if self[i, j] == 1: # Skip multiply if identity z = x[j] else: # Do the block multiply if isinstance(self[i, j], GenericMatrix): z = self[i, j].create_vec(dim=0) self[i, j].mult(x[j], z) else: z = self[i, j] * x[j] if z == NotImplemented: return NotImplemented if not isinstance(z, (GenericVector, block_vec)): # Typically, this happens when for example a # block_vec([0,0]) is used without calling allocate() or # setting BCs. The result is a Matrix*scalar=Matrix. One # way to fix this issue is to convert all scalars to a # proxy class in block_vec.__init__, and let this proxy # class have a __rmul__ that converts to vector on demand. # (must also stop conversion anything->blockcompose for # this case) raise RuntimeError( 'unexpected result in matvec, %s\n-- possibly because a block_vec contains scalars ' \ 'instead of vectors, use create_vec() or allocate()' % type(z)) if y[i] is None: y[i] = z else: if len(y[i]) != len(z): raise RuntimeError( 'incompatible dimensions in block (%d,%d) -- %d, was %d' % (i, j, len(z), len(y[i]))) y[i] += z return y
def block_simplify(self): """Try to simplify the transpose, recursively. A fuller explanation is found in block_transform.block_simplify. """ from block_util import isscalar from block_transform import block_simplify A = block_simplify(self.A) if isscalar(A): return A if isinstance(A, block_transpose): return A.A return block_transpose(A)
def block_kronecker(A, B): """Create the Kronecker (tensor) product of two matrices. The result is returned as a product of two block matrices, (A x Ib) * (Ia x B), where Ia and Ib are the appropriate identities (A=A*Ia, B=Ib*B). This will often limit the number of repeated applications of the inner operators. Note: If the A operator is a DOLFIN matrix and B is not a block_mat, A will be expanded into a block_mat (one block per matrix entry). This will not work in parallel. Apart from that, no blocks are expanded, i.e., only the block matrices are expanded in the Kronecker product. To form the Kronecker sum, you can extract (A x Ib) and (Ia x B) like this: C,D = block_kronecker(A,B); sum=C+D Similarly, it may be wise to do the inverse separately: C,D = block_kronecker(A,B); inverse = some_invert(D)*ConjGrad(C) """ from block_util import isscalar import dolfin if isinstance(A, dolfin.GenericMatrix) and not isinstance(B, block_mat): A = block_mat(A.array()) assert isinstance(A, block_mat) or isinstance(B, block_mat) if isinstance(B, block_mat): p,q = B.blocks.shape C = block_mat.diag(A, n=p) else: C = A if isinstance(A, block_mat): m,n = A.blocks.shape if isinstance(B, block_mat): print "Note: block_kronecker with two block matrices is probably buggy" D = block_mat(n,n) for i in range(n): for j in range(n): # A scalar can represent the scaled identity of any # dimension, so no need to diagonal-expand it here. We # don't do this check on the outer diagonal expansions, # because it is clearer to always return two block matrices # of equal dimension rather than sometimes a scalar. b = B[i,j] D[i,j] = b if isscalar(b) else block_mat.diag(b, n=m) else: D = block_mat.diag(B, n=m) else: D = B return block_mul(C,D)
def block_simplify(self): """Try to combine scalar terms and remove multiplicative identities, recursively. A fuller explanation is found in block_transform.block_simplify. """ from block_util import isscalar from block_transform import block_simplify operators = [] scalar = 1.0 for op in self.chain: op = block_simplify(op) if isscalar(op): scalar *= op else: operators.append(op) if scalar == 0: return 0 if scalar != 1 or len(operators) == 0: operators.insert(0, scalar) if len(operators) == 1: return operators[0] ret = block_mul(None, None) ret.chain = operators return ret
def block_assemble(lhs, rhs=None, bcs=None, symmetric=False, signs=None, symmetric_mod=None): """ Assembles block matrices, block vectors or block systems. Input can be arrays of variational forms or block matrices/vectors. Arguments: symmetric : Boundary conditions are applied so that symmetry of the system is preserved. If only the left hand side of the system is given, then a matrix represententing the rhs corrections is returned along with a symmetric matrix. symmetric_mod : Matrix describing symmetric corrections for assembly of the of the rhs of a variational system. signs : An array to specify the signs of diagonal blocks. The sign of the blocks are computed if the argument is not provided. """ error_msg = {'incompatibility' : 'A and b do not have compatible dimensions.', 'symm_mod error' : 'symmetric_mod argument only accepted when assembling a vector', 'not square' : 'A must be square for symmetric assembling', 'invalid bcs' : 'Expecting a list or list of lists of DirichletBC.', 'invalid signs' : 'signs should be a list of length n containing only 1 or -1', 'mpi and symm' : 'Symmetric application of BC not yet implemented in parallel'} # Check arguments from numpy import ndarray has_rhs = True if isinstance(rhs, ndarray) else rhs != None has_lhs = True if isinstance(rhs, ndarray) else rhs != None if symmetric: from dolfin import MPI, mpi_comm_world if MPI.size(mpi_comm_world()) > 1: raise NotImplementedError(error_msg['mpi and symm']) if has_lhs and has_rhs: A, b = map(block_tensor,[lhs,rhs]) n, m = A.blocks.shape if not ( isinstance(b,block_vec) and len(b.blocks) is m): raise TypeError(error_msg['incompatibility']) else: A, b = block_tensor(lhs), None if isinstance(A,block_vec): A, b = None, A n, m = 0, len(b.blocks) else: n,m = A.blocks.shape if A and symmetric and (m is not n): raise RuntimeError(error_msg['not square']) if symmetric_mod and ( A or not b ): raise RuntimeError(error_msg['symmetric_mod error']) # First assemble everything needing assembling. from dolfin import assemble assemble_if_form = lambda x: assemble(x, keep_diagonal=True) if _is_form(x) else x if A: A.blocks.flat[:] = map(assemble_if_form,A.blocks.flat) if b: #b.blocks.flat[:] = map(assemble_if_form, b.blocks.flat) b = block_vec(map(assemble_if_form, b.blocks.flat)) # If there are no boundary conditions then we are done. if bcs is None: if A: return [A,b] if b else A else: return b # check if arguments are forms, in which case bcs have to be split from ufl import Form if isinstance(lhs, Form): from splitting import split_bcs bcs = split_bcs(bcs, m) # Otherwise check that boundary conditions are valid. if not hasattr(bcs,'__iter__'): raise TypeError(error_msg['invalid bcs']) if len(bcs) is not m: raise TypeError(error_msg['invalid bcs']) from dolfin import DirichletBC for bc in bcs: if isinstance(bc,DirichletBC) or bc is None: pass else: if not hasattr(bc,'__iter__'): raise TypeError(error_msg['invalid bcs']) else: for bc_i in bc: if isinstance(bc_i,DirichletBC): pass else: raise TypeError(error_msg['invalid bcs']) bcs = [bc if hasattr(bc,'__iter__') else [bc] if bc else bc for bc in bcs] # Apply BCs if we are only assembling the righ hand side if not A: if symmetric_mod: b.allocate(symmetric_mod) for i in xrange(m): if bcs[i]: if isscalar(b[i]): b[i], val = create_vec_from(bcs[i][0]), b[i] b[i][:] = val for bc in bcs[i]: bc.apply(b[i]) if symmetric_mod: b.allocate(symmetric_mod) b -= symmetric_mod*b return b # If a signs argument is passed, check if it is valid. # Otherwise guess. if signs and symmetric: if ( hasattr(signs,'__iter__') and len(signs)==m ): for sign in signs: if sign not in (-1,1): raise TypeError(error_msg['invalid signs']) else: raise TypeError(error_msg['invalid signs']) elif symmetric: from numpy.random import random signs = [0]*m for i in xrange(m): if isscalar(A[i,i]): signs[i] = -1 if A[i,i] < 0 else 1 else: x = A[i,i].create_vec(dim=1) x.set_local(random(x.local_size())) signs[i] = -1 if x.inner(A[i,i]*x) < 0 else 1 # Now apply boundary conditions. if b: b.allocate(A) elif symmetric: # If we are preserving symmetry but don't have the rhs b, # then we need to store the symmetric corretions to b # as a matrix which we call A_mod b, A_mod = A.create_vec(), A.copy() for i in xrange(n): if bcs[i]: for bc in bcs[i]: # Apply BCs to the diagonal block. if isscalar(A[i,i]): A[i,i] = _new_square_matrix(bc,A[i,i]) if symmetric: A_mod[i,i] = A[i,i].copy() if symmetric: bc.zero_columns(A[i,i],b[i],signs[i]) bc.apply(A_mod[i,i]) elif b: bc.apply(A[i,i],b[i]) else: bc.apply(A[i,i]) # Zero out the rows corresponding to BC dofs. for j in range(i) + range(i+1,n): if A[i,j] is 0: continue assert not isscalar(A[i,j]) bc.zero(A[i,j]) # If we are not preserving symmetry then we are done at this point. # Otherwise, we need to zero out the columns as well if symmetric: for j in range(i) + range(i+1,n): if A[j,i] is 0: continue assert not isscalar(A[j,i]) bc.zero_columns(A[j,i],b[j]) bc.zero(A_mod[i,j]) result = [A] if symmetric: for i in range(n): for j in range(n): A_mod[i,j] -= A[i,j] result += [A_mod] if b: result += [b] return result[0] if len(result)==1 else result