def test_DoGIP_vs_FEniCS(self): print( '\n== testing DoGIP vs. FEniCS for problem of weighted projection ====' ) for dim, pol_order in itertools.product([2, 3], [1, 2]): print('dim={}; pol_order={}'.format(dim, pol_order)) N = 2 # creating MESH, defining MATERIAL and SOURCE if dim == 2: mesh = UnitSquareMesh(N, N) m = Expression("1+10*16*x[0]*(1-x[0])*x[1]*(1-x[1])", degree=3) # material coefficients f = Expression("x[0]*x[0]*x[1]", degree=2) elif dim == 3: mesh = UnitCubeMesh(N, N, N) m = Expression("1+100*x[0]*(1-x[0])*x[1]*x[2]", degree=2) # material coefficients f = Expression("(1-x[0])*x[1]*x[2]", degree=2) mesh.coordinates()[:] += 0.1 * np.random.random( mesh.coordinates().shape) # mesh perturbation ## standard approach with FEniCS ############################################# V = FunctionSpace(mesh, "CG", pol_order) # original FEM space u, v = TrialFunction(V), TestFunction(V) u_fenics = Function(V) solve(m * u * v * dx == m * f * v * dx, u_fenics) ## DoGIP - double-grid integration with interpolation-projection ############# W = FunctionSpace(mesh, "CG", 2 * pol_order) # double-grid space w = TestFunction(W) A_dogip = assemble( m * w * dx).get_local() # diagonal matrix of material coefficients b = assemble(m * f * v * dx) # vector of right-hand side # assembling interpolation-projection matrix B B = get_B(V, W, problem=0) # # linear solver on double grid, standard Afun = lambda x: B.T.dot(A_dogip * B.dot(x)) Alinoper = linalg.LinearOperator((V.dim(), V.dim()), matvec=Afun, dtype=np.float) x, info = linalg.cg(Alinoper, b.get_local(), x0=np.zeros(V.dim()), tol=1e-10, maxiter=1e3, callback=None) # testing the difference between DoGIP and FEniCS self.assertAlmostEqual( 0, np.linalg.norm(u_fenics.vector().get_local() - x)) print('...ok')
def get_Bhat(dim=2, pol_order=1, problem=1, W_space=True): """ Calculate interpolation matrix between V and W on a reference element Parameters ---------- dim : int topological dimension of the problem por_order : int polynomial order of finite element space problem : {int, string} parameter defining problem: 0 - weighted projection, 1 - elliptic problem Returns ------- B : ndarray interpolation matrix between original and double-grid finite element basis """ mesh = get_reference_element(dim) cell = Cell(mesh, 0) V = FunctionSpace(mesh, "CG", pol_order) # original FEM space if W_space: if problem in [1, 'elliptic']: W = FunctionSpace(mesh, "DG", 2 * (pol_order - 1)) # double-grid space B = np.zeros([W.dim(), V.dim(), dim]) elif problem in [0, 'projection']: W = FunctionSpace(mesh, "CG", 2 * pol_order) # double-grid space B = np.zeros([W.dim(), V.dim()]) dofcoors = W.element().tabulate_dof_coordinates(cell) else: if problem in [1, 'elliptic']: dofcoors = get_dof_coordinates(dim, 2 * (pol_order - 1)) B = np.zeros([dofcoors.shape[0], V.dim(), dim]) elif problem in [0, 'projection']: dofcoors = get_dof_coordinates(dim, 2 * pol_order) B = np.zeros([dofcoors.shape[0], V.dim()]) elementV = V.element() if problem in [1, 'elliptic']: eval_basis = lambda *args: elementV.evaluate_basis_derivatives_all( *args) der = (1, ) postprocess = lambda B: np.einsum('jkr->rjk', B) val_shape = (V.dim(), dim) elif problem in [0, 'projection']: eval_basis = lambda *args: elementV.evaluate_basis_all(*args) der = () postprocess = lambda B: B val_shape = (V.dim()) for jj, dofcoor in enumerate(dofcoors): args = der + (dofcoor, cell.get_vertex_coordinates(), cell.orientation()) B[jj] = eval_basis(*args).reshape(val_shape) return postprocess(B)
def test_multiplication(self): print( '\n== testing multiplication of system matrix for problem of weighted projection ====' ) for dim, pol_order in itertools.product([2, 3], [1, 2]): N = 2 # no. of elements print('dim={0}, pol_order={1}, N={2}'.format(dim, pol_order, N)) # creating MESH and defining MATERIAL if dim == 2: mesh = UnitSquareMesh(N, N) m = Expression("1+10*16*x[0]*(1-x[0])*x[1]*(1-x[1])", degree=2) # material coefficients elif dim == 3: mesh = UnitCubeMesh(N, N, N) m = Expression("1+10*16*x[0]*(1-x[0])*(1-x[1])*x[2]", degree=2) # material coefficients mesh.coordinates()[:] += 0.1 * np.random.random( mesh.coordinates().shape) # mesh perturbation V = FunctionSpace(mesh, "CG", pol_order) # original FEM space W = FunctionSpace(mesh, "CG", 2 * pol_order) # double-grid space print('assembling local matrices for DoGIP...') Bhat = get_Bhat( dim, pol_order, problem=0) # projection between V on W on a reference element AT_dogip = get_A_T(m, V, W, problem=0) dofmapV = V.dofmap() def system_multiplication_DoGIP(AT_dogip, Bhat, u_vec): # mutliplication with DoGIP decomposition Au = np.zeros_like(u_vec) for ii, cell in enumerate(cells(mesh)): ind = dofmapV.cell_dofs(ii) # local to global map Au[ind] += Bhat.T.dot(AT_dogip[ii] * Bhat.dot(u_vec[ind])) return Au print('assembling FEM sparse matrix') u, v = TrialFunction(V), TestFunction(V) Asp = assemble(m * u * v * dx, tensor=EigenMatrix()) # Asp = Asp.sparray() print('multiplication...') ur = Function(V) # creating random vector ur_vec = 5 * np.random.random(V.dim()) ur.vector().set_local(ur_vec) Au_DoGIP = system_multiplication_DoGIP( AT_dogip, Bhat, ur_vec) # DoGIP multiplication Auex = Asp.dot(ur_vec) # FEM multiplication with sparse matrix # testing the difference between DoGIP and FEniCS self.assertAlmostEqual(0, np.linalg.norm(Auex - Au_DoGIP)) print('...ok')
class Solution(): def __init__(self, prefix, noperiodic=False): """Access a KSDG solution. Required argument: prefix: This is the prefix for the names of the files in which the solution is stored. (Typically, it would be the value of the --save option to ksdgsolver<d>.py. optional argument: noperiodic=False: Treat solution as periodic, even if it wasn't. """ self.prefix = os.path.expanduser(prefix) self.prefix = os.path.expandvars(self.prefix) self.noperiodic = noperiodic self.meshfile = self.prefix + '_omesh.xml.gz' if os.path.isfile(self.meshfile): self.mesh = self.omesh = Mesh(self.meshfile) else: self.meshfile = prefix + '_mesh.xml.gz' self.mesh = Mesh(self.meshfile) self.optionsfile = self.prefix + '_options.txt' with open(self.optionsfile, 'r') as optsf: self.options = optsf.readlines() self.fsinfofile = fsinfo_filename(self.prefix, 0) with h5py.File(self.fsinfofile, 'r') as fsf: self.dim, self.degree = ( int(fsf['dim'][()]), int(fsf['degree'][()]), ) try: # self.periodic = fsf['periodic'].value #deprecated self.periodic = fsf['periodic'][()] except KeyError: self.periodic = False try: self.ligands = pickle.loads(fsf['ligands'][()]) except KeyError: self.ligands = False try: self.param_names = pickle.loads(fsf['param_names'][()]) self.variable = True except KeyError: self.param_names = [] self.variable = False try: self.params0 = pickle.loads(fsf['params0'][()]) except KeyError: self.params0 = {} try: pfuncs = fsf['param_funcs'][()] self.param_funcs = dill.loads(pfuncs.tobytes()) except (KeyError, ValueError): self.param_funcs = {} def identity(t, params={}): return t self.param_funcs['t'] = identity if self.params0 and self.ligands: capp = [s for s in self.options if '--cappotential' in s] if 'tophat' in capp[0]: self.cappotential = 'tophat' elif 'witch' in capp[0]: self.cappotential = 'witch' self.Vgroups = copy.deepcopy(self.ligands) self.Vparams = ParameterList(default_parameters) self.Vparams.add(self.Vgroups.params()) def Vfunc(Us, params={}): self.Vparams.update(params) # copy params into ligands return self.Vgroups.V(Us) # compute V def Vtophat(rho, params={}): tanh = ufl.tanh((rho - params['rhomax']) / params['cushion']) return params['maxscale'] * params['sigma']**2 / 2 * (tanh + 1) def Vwitch(rho, params={}): tanh = ufl.tanh((rho - params['rhomax']) / params['cushion']) return (params['maxscale'] * params['sigma']**2 / 2 * (tanh + 1) * (rho / params['rhomax'])) Vcap = Vwitch if self.cappotential == 'witch' else Vtophat def V2(Us, rho, params={}): return Vfunc(Us, params=params) + Vcap(rho, params=params) self.V = V2 self.tsfile = prefix + '_ts.h5' self.ts = self.timeSeries = KSDGTimeSeries(self.tsfile, 'r') self.tstimes = self.ts.sorted_times() self.tmin, self.tmax = self.tstimes[0], self.tstimes[-1] kwargs = dict(mesh=self.mesh, dim=self.dim, degree=self.degree, t0=self.tmin, ligands=self.ligands, parameters=self.params0, V=V2, periodic=self.periodic) if self.variable: kwargs['param_funcs'] = self.param_funcs self.ksdg = makeKSDGSolver(**kwargs) # self.V = self.ksdg.V if self.periodic and self.noperiodic: self.mesh = self.ksdg.mesh self.meshstats = mesh_stats(self.ksdg.mesh) else: try: self.meshstats = mesh_stats(self.ksdg.omesh) except AttributeError: self.meshstats = mesh_stats(self.ksdg.mesh) if self.periodic and not self.noperiodic: self.ex_mesh = ExpandedMesh(self.ksdg.sol) self.function = self.ex_mesh.expanded_function else: self.ex_mesh = None self.function = self.ksdg.sol self.fs = self.function.function_space() self.transform = integerify_transform( np.reshape(self.fs.tabulate_dof_coordinates(), (-1, self.dim))) def params(self, t): pd = collections.OrderedDict([ (name, self.param_funcs[name](t, params=self.params0)) for name in self.param_names ]) return pd def load(self, t): """Load solution for time t.""" vec = self.ts.retrieve_by_time(t) if self.periodic and not self.noperiodic: self.ex_mesh.sub_function.vector()[:] = vec self.ex_mesh.sub_function.vector().apply('insert') self.ex_mesh.expand() else: self.ksdg.sol.vector()[:] = vec self.ksdg.sol.vector().apply('insert') def project(self, t=None): """Project solution onto CG subspace. soln.project(t) project projects the solution at time t onto a CG FunctionSpace. The projected function is left in self.CGfunction. If argument t is not supplied, the currently loaded solution will be used. (In this case, you must have called load(t) before calling project(). project returns self.CGfunction as its value. """ solver_type = 'petsc' if t is not None: self.load(t) if not hasattr(self, 'CS'): # # create CG FunctionSpace # self.CE = CGelement(self.fs) self.CS = FunctionSpace(self.fs.mesh(), self.CE) self.CGfunction = fe.project(self.function, self.CS, solver_type=solver_type) return self.CGfunction def images(self, t=None): self.project(t) if not hasattr(self, 'ims'): self.CGdofs = local_dofs(self.CS) order = np.argsort(self.CGdofs[:, -1]) self.CGdofs = self.CGdofs[order] self.CGintdofs = np.empty_like(self.CGdofs, dtype=int) self.CGintdofs[:, :2] = np.rint(self.CGdofs[:, :2], casting='unsafe') self.CGintdofs[:, -1] = np.rint(self.CGdofs[:, -1], casting='unsafe') self.CGintdofs[:, 2:-1] = integerify(self.CGdofs[:, 2:-1], transform=self.transform) self.sss = np.unique(self.CGintdofs[:, 1]) self.nss = len(self.sss) assert np.alltrue(self.sss == np.arange(self.nss)) self.dims = np.zeros(shape=(self.dim, ), dtype=int) for i in range(self.dim): ixs = np.unique(self.CGintdofs[:, 2 + i]) self.dims[i] = len(ixs) assert np.alltrue(ixs == np.arange(self.dims[i])) self.imshape = (self.nss, ) + tuple(self.dims) assert self.CS.dim() == np.prod(self.imshape) checknums = 1 + np.arange(self.CS.dim()) self.iims = np.zeros(shape=self.imshape, dtype=int) self.CGixs = ((self.CGintdofs[:, 1], ) + tuple(self.CGintdofs[:, 2:-1].transpose())) self.iims[self.CGixs] = checknums assert np.alltrue(self.iims > 0) del (self.iims) self.ims = np.empty(shape=self.imshape) self.ims[self.CGixs] = self.CGfunction.vector()[:] return self.ims def Vufl(self, t): "Make a UFL object for plotting V" self.project(t) split = fe.split(self.function) irho = split[0] iUs = split[1:] self.ksdg.set_time(t) V = (self.V(iUs, irho, params=self.ksdg.iparams) / (self.ksdg.iparams['sigma']**2 / 2)) fs = self.function.function_space() cell = fs.ufl_cell() self.SE = fe.FiniteElement('CG', cell, self.degree) self.SS = fe.FunctionSpace(fs.mesh(), self.SE) pV = fe.project(V, self.SS, solver_type='petsc') return pV def Vimage(self, t=None): pV = self.Vufl(t) self.Vimshape = tuple(self.dims) self.VCGcoords = np.reshape(self.SS.tabulate_dof_coordinates(), (-1, self.dim)) self.VCGintcoords = integerify(self.VCGcoords, transform=self.transform) self.Vixs = (tuple(self.VCGintcoords.transpose())) self.Vimg = np.zeros(self.Vimshape, dtype=float) self.Vimg[self.Vixs] = pV.vector()[:] return self.Vimg