def check_Q_transfer_minimal(collclass): """ A simple test program to check the order of the Q interpolation/restriction for only 2 coarse nodes """ Mcoarse = 2 coll_coarse = collclass(Mcoarse, 0, 1) for M in range(3, 9): Mfine = M coll_fine = collclass(Mfine, 0, 1) assert coll_fine.left_is_node == coll_coarse.left_is_node, 'ERROR: should be using the same class for coarse and fine Q' if not coll_fine.left_is_node: fine_grid = np.concatenate(([0], coll_fine.nodes)) coarse_grid = np.concatenate(([0], coll_coarse.nodes)) else: fine_grid = coll_fine.nodes coarse_grid = coll_coarse.nodes Pcoll = th.interpolation_matrix_1d(fine_grid, coarse_grid, k=2, pad=0) Rcoll = th.restriction_matrix_1d(fine_grid, coarse_grid, k=2, pad=0) for polyorder in range(1, 3): coeff = np.random.rand(polyorder) ufine = polyval(fine_grid, coeff) ucoarse = polyval(coarse_grid, coeff) uinter = Pcoll.dot(ucoarse) urestr = Rcoll.dot(ufine) err_inter = np.linalg.norm(uinter - ufine, np.inf) err_restr = np.linalg.norm(urestr - ucoarse, np.inf) if polyorder <= 2: assert err_inter < 2E-15, "ERROR: Q-interpolation order is not reached, got %s" % err_inter assert err_restr < 2E-15, "ERROR: Q-restriction order is not reached, got %s" % err_restr else: assert err_inter > 2E-15, "ERROR: Q-interpolation order is higher than expected, got %s" % polyorder
def check_Q_transfer(collclass): """ A simple test program to check the order of the Q interpolation/restriction """ for M in range(3, 9): Mfine = M Mcoarse = int((Mfine+1)/2.0) coll_fine = collclass(Mfine, 0, 1) coll_coarse = collclass(Mcoarse, 0, 1) assert coll_fine.left_is_node == coll_coarse.left_is_node, 'ERROR: should be using the same class for coarse and fine Q' fine_grid = coll_fine.nodes coarse_grid = coll_coarse.nodes for order in range(2,coll_coarse.num_nodes+1): Pcoll = th.interpolation_matrix_1d(fine_grid, coarse_grid, k=order, pad=0, equidist_nested=False) Rcoll = th.restriction_matrix_1d(fine_grid, coarse_grid, k=order, pad=0) for polyorder in range(1,order+2): coeff = np.random.rand(polyorder) ufine = polyval(fine_grid,coeff) ucoarse = polyval(coarse_grid,coeff) uinter = Pcoll.dot(ucoarse) urestr = Rcoll.dot(ufine) err_inter = np.linalg.norm(uinter-ufine, np.inf) err_restr = np.linalg.norm(urestr-ucoarse, np.inf) if polyorder <= order: assert err_inter < 2E-15, "ERROR: Q-interpolation order is not reached, got %s" %err_inter assert err_restr < 2E-15, "ERROR: Q-restriction order is not reached, got %s" % err_restr else: assert err_inter > 2E-15, "ERROR: Q-interpolation order is higher than expected, got %s" % polyorder
def Jacobi(): N = 127 dx = 1.0 / (N + 1) nu = 1.0 K = 20 stencil = [-1, 2, -1] A = sp.sparse.diags(stencil, [-1, 0, 1], shape=(N, N), format='csc') A *= nu / (dx ** 2) D = sp.sparse.diags(2.0 * A.diagonal(), 0, shape=(N, N), format='csc') Dinv = sp.sparse.diags(0.5 * 1.0/A.diagonal(), 0, shape=(N, N), format='csc') f = np.ones(N) f = np.zeros(N) u = np.sin([int(3.0 * N / 4.0) * np.pi * (i+1) * dx for i in range(N)]) res = f - A.dot(u) for k in range(1, K+1): u += Dinv.dot(res) res = f - A.dot(u) print(k, np.linalg.norm(res, np.inf)) print() # print(u) tol = np.linalg.norm(res, np.inf) Nc = int((N + 1) / 2 - 1) dxc = 1.0 / (Nc + 1) Ac = sp.sparse.diags(stencil, [-1, 0, 1], shape=(Nc, Nc), format='csc') Ac *= nu / (dxc ** 2) Dc = sp.sparse.diags(2.0 * Ac.diagonal(), 0, shape=(Nc, Nc), format='csc') Dcinv = sp.sparse.diags(0.5 * 1.0 / Ac.diagonal(), 0, shape=(Nc, Nc), format='csc') fine_grid = np.array([(i + 1) * dx for i in range(N)]) coarse_grid = np.array([(i + 1) * dxc for i in range(Nc)]) I = sp.sparse.csc_matrix(interpolation_matrix_1d(fine_grid, coarse_grid, k=6, periodic=False, equidist_nested=True)) R = sp.sparse.csc_matrix(I.T) T = sp.sparse.csc_matrix(sp.sparse.eye(N) - Dinv.dot(A)) Tc = sp.sparse.csc_matrix(I.dot(sp.sparse.eye(Nc) - Dcinv.dot(Ac)).dot(R)) fvec = np.kron(np.ones(K), Dinv.dot(f)) u = np.zeros(N * K) u = np.kron(np.ones(K), np.sin([int(3.0 * N / 4.0) * np.pi * (i + 1) * dx for i in range(N)])) # u[0: N] = np.sin([int(3.0 * N / 4.0) * np.pi * (i + 1) * dx for i in range(N)]) res = f - A.dot(u[0:N]) uold = u.copy() l = 0 while np.linalg.norm(res, np.inf) > tol and l < K: for k in range(1, K): u[k * N: (k+1) * N] = T.dot(uold[(k-1) * N: k * N]) + Tc.dot(u[(k-1) * N: k * N] - uold[(k-1) * N: k * N]) + fvec[(k-1) * N: k * N] res = f - A.dot(u[-N:]) l += 1 uold = u.copy() print(l, np.linalg.norm(res, np.inf)) print() # print(u[-N:]) E = np.zeros((K, K)) np.fill_diagonal(E[1:, :], 1) E = sp.sparse.csc_matrix(E) Rfull = sp.sparse.kron(sp.sparse.eye(K), R) Ifull = sp.sparse.kron(sp.sparse.eye(K), I) S = sp.sparse.kron(sp.sparse.eye(K), D) - sp.sparse.kron(E, D - A) Sc = Rfull.dot(S).dot(Ifull) # Sc = sp.sparse.kron(sp.sparse.eye(K), Dcinv) - sp.sparse.kron(E, Dc - Ac) Scinv = np.linalg.inv(Sc.todense()) # Sinv = np.linalg.inv(S.todense()) Sdiaginv = sp.sparse.kron(sp.sparse.eye(K), Dinv) Scdiaginv = sp.sparse.kron(sp.sparse.eye(K), Dcinv) u = np.zeros(N * K) # u[0: N] = np.sin([int(3.0 * N / 4.0) * np.pi * (i + 1) * dx for i in range(N)]) u = np.kron(np.ones(K), np.sin([int(3.0 * N / 4.0) * np.pi * (i + 1) * dx for i in range(N)])) l = 0 fvec = np.kron(np.ones(K), f) res = f - A.dot(u[0: N]) while np.linalg.norm(res, np.inf) > tol and l < K: # u += Sainv.dot(u0vec - S.dot(u)) u += Sdiaginv.dot(fvec - S.dot(u)) uH = Rfull.dot(u) uHold = uH.copy() rhsH = Rfull.dot(fvec) + Sc.dot(uH) - Rfull.dot(S.dot(u)) uH += np.ravel(Scinv.dot(rhsH)) # uH += Scdiaginv.dot(rhsH - Sc.dot(uH)) # uH2 = R2full.dot(uH) # uH2old = uH2.copy() # rhsH2 = R2full.dot(rhsH) + Sc2.dot(uH2) - R2full.dot(Sc.dot(uH)) # uH2 = Sc2inv.dot(rhsH2) # uH += I2full.dot(uH2 - uH2old) # uH += Scdiaginv.dot(rhsH - Sc.dot(uH)) u += Ifull.dot(uH - uHold) u += Sdiaginv.dot(fvec - S.dot(u)) res = f - A.dot(u[-N:]) l += 1 print(l, np.linalg.norm(res, np.inf)) # print(u[-N:]) print()
def __init__(self, fine_prob, coarse_prob, params): """ Initialization routine Args: fine_prob: fine problem coarse_prob: coarse problem params: parameters for the transfer operators """ if 'iorder' not in params: raise TransferError('Need iorder parameter for spatial transfer') if 'rorder' not in params: raise TransferError('Need rorder parameter for spatial transfer') # invoke super initialization super(mesh_to_mesh, self).__init__(fine_prob, coarse_prob, params) if type(self.fine_prob.params.nvars) is tuple: if type(self.coarse_prob.params.nvars) is not tuple: raise TransferError( 'nvars parameter of coarse problem needs to be a tuple') if not len(self.fine_prob.params.nvars) == len( self.coarse_prob.params.nvars): raise TransferError( 'nvars parameter of fine and coarse level needs to have the same length' ) elif type(self.fine_prob.params.nvars) is int: if type(self.coarse_prob.params.nvars) is not int: raise TransferError( 'nvars parameter of coarse problem needs to be an int') else: raise TransferError("unknow type of nvars for transfer, got %s" % self.fine_prob.params.nvars) # we have a 1d problem if type(self.fine_prob.params.nvars) is int: # if number of variables is the same on both levels, Rspace and Pspace are identity if self.coarse_prob.params.nvars == self.fine_prob.params.nvars: self.Rspace = sp.eye(self.coarse_prob.params.nvars) self.Pspace = sp.eye(self.fine_prob.params.nvars) # assemble restriction as transpose of interpolation else: if not self.params.periodic: fine_grid = np.array([ (i + 1) * self.fine_prob.dx for i in range(self.fine_prob.params.nvars) ]) coarse_grid = np.array([ (i + 1) * self.coarse_prob.dx for i in range(self.coarse_prob.params.nvars) ]) else: fine_grid = np.array([ i * self.fine_prob.dx for i in range(self.fine_prob.params.nvars) ]) coarse_grid = np.array([ i * self.coarse_prob.dx for i in range(self.coarse_prob.params.nvars) ]) if self.params.rorder <= 1: self.Rspace = th.restriction_matrix_1d( fine_grid, coarse_grid, k=self.params.rorder, periodic=self.params.periodic) else: self.Rspace = 0.5 * th.interpolation_matrix_1d( fine_grid, coarse_grid, k=self.params.rorder, periodic=self.params.periodic).T self.Pspace = th.interpolation_matrix_1d( fine_grid, coarse_grid, k=self.params.iorder, periodic=self.params.periodic) # we have an n-d problem else: Rspace = [] Pspace = [] for i in range(len(self.fine_prob.params.nvars)): # if number of variables is the same on both levels, Rspace and Pspace are identity if self.coarse_prob.params.nvars == self.fine_prob.params.nvars: Rspace.append(sp.eye(self.coarse_prob.params.nvars[i])) Pspace.append(sp.eye(self.fine_prob.params.nvars[i])) # assemble restriction as transpose of interpolation else: if not self.params.periodic: fine_grid = np.array([ (j + 1) * self.fine_prob.dx for j in range(self.fine_prob.params.nvars[i]) ]) coarse_grid = np.array([ (j + 1) * self.coarse_prob.dx for j in range(self.coarse_prob.params.nvars[i]) ]) else: fine_grid = np.array([ j * self.fine_prob.dx for j in range(self.fine_prob.params.nvars[i]) ]) coarse_grid = np.array([ j * self.coarse_prob.dx for j in range(self.coarse_prob.params.nvars[i]) ]) if self.params.iorder <= 1: Rspace.append( th.restriction_matrix_1d( fine_grid, coarse_grid, k=self.params.iorder, periodic=self.params.periodic).T) else: Rspace.append(0.5 * th.interpolation_matrix_1d( fine_grid, coarse_grid, k=self.params.iorder, periodic=self.params.periodic).T) Pspace.append( th.interpolation_matrix_1d( fine_grid, coarse_grid, k=self.params.iorder, periodic=self.params.periodic)) # kronecker 1-d operators for n-d self.Pspace = Pspace[0] for i in range(1, len(Pspace)): self.Pspace = sp.kron(self.Pspace, Pspace[i], format='csc') self.Rspace = Rspace[0] for i in range(1, len(Rspace)): self.Rspace = sp.kron(self.Rspace, Rspace[i], format='csc')
def __init__(self, fine_level, coarse_level, base_transfer_params, space_transfer_class, space_transfer_params): """ Initialization routine Args: fine_level (pySDC.Level.level): fine level connected with the base_transfer operations coarse_level (pySDC.Level.level): coarse level connected with the base_transfer operations base_transfer_params (dict): parameters for the base_transfer operations space_transfer_class: class to perform spatial transfer space_transfer_params (dict): parameters for the space_transfer operations """ # short helper class to add params as attributes class __Pars(FrozenClass): def __init__(self, pars): self.finter = False self.coll_iorder = 2 self.coll_rorder = 1 for k, v in pars.items(): setattr(self, k, v) self._freeze() self.params = __Pars(base_transfer_params) # set up logger self.logger = logging.getLogger('transfer') # just copy by object self.fine = fine_level self.coarse = coarse_level # for Q-based transfer, check if we need to add 0 to the list of nodes if not self.fine.sweep.coll.left_is_node: fine_grid = np.concatenate(([0], self.fine.sweep.coll.nodes)) coarse_grid = np.concatenate(([0], self.coarse.sweep.coll.nodes)) else: fine_grid = self.fine.sweep.coll.nodes coarse_grid = self.coarse.sweep.coll.nodes if self.params.coll_iorder > len(coarse_grid): self.logger.warning( 'requested order of Q-interpolation is not valid, resetting to %s' % len(coarse_grid)) self.params.coll_iorder = len(coarse_grid) if self.params.coll_rorder != 1: self.logger.warning( 'requested order of Q-restriction is != 1, can lead to weird behavior!' ) # set up preliminary transfer matrices for Q-based coarsening Pcoll = th.interpolation_matrix_1d(fine_grid, coarse_grid, k=self.params.coll_iorder, pad=0).toarray() Rcoll = th.restriction_matrix_1d(fine_grid, coarse_grid, k=self.params.coll_rorder, pad=0).toarray() # pad transfer matrices if necessary if self.fine.sweep.coll.left_is_node: self.Pcoll = np.zeros((self.fine.sweep.coll.num_nodes + 1, self.coarse.sweep.coll.num_nodes + 1)) self.Rcoll = np.zeros((self.coarse.sweep.coll.num_nodes + 1, self.fine.sweep.coll.num_nodes + 1)) self.Pcoll[1:, 1:] = Pcoll self.Rcoll[1:, 1:] = Rcoll else: self.Pcoll = Pcoll self.Rcoll = Rcoll # set up spatial transfer self.space_transfer = space_transfer_class( fine_prob=self.fine.prob, coarse_prob=self.coarse.prob, params=space_transfer_params)