def step(self, u1, u0, t, dt, tol=1.0e-10, maxiter=100, verbose=True ): v = TestFunction(self.problem.V) L, b = self.problem.get_system(t) Lu0 = Function(self.problem.V) L.mult(u0.vector(), Lu0.vector()) rhs = assemble(u0 * v * self.problem.dx_multiplier * self.problem.dx) rhs.axpy(dt, -Lu0.vector() + b) A = self.M # Apply boundary conditions. for bc in self.problem.get_bcs(t + dt): bc.apply(A, rhs) # Both Jacobi-and AMG-preconditioners are order-optimal for the mass # matrix, Jacobi is a little more lightweight. solver = KrylovSolver('cg', 'jacobi') solver.parameters['relative_tolerance'] = tol solver.parameters['absolute_tolerance'] = 0.0 solver.parameters['maximum_iterations'] = maxiter solver.parameters['monitor_convergence'] = verbose solver.set_operator(A) solver.solve(u1.vector(), rhs) return
def step(self, u1, u0, t, dt, tol=1.0e-10, verbose=True, maxiter=1000, krylov='gmres', preconditioner='ilu' ): v = TestFunction(self.problem.V) L0, b0 = self.problem.get_system(t) L1, b1 = self.problem.get_system(t + dt) Lu0 = Function(self.problem.V) L0.mult(u0.vector(), Lu0.vector()) rhs = assemble(u0 * v * self.problem.dx_multiplier * self.problem.dx) rhs.axpy(-dt * 0.5, Lu0.vector()) rhs.axpy(+dt * 0.5, b0 + b1) A = self.M + 0.5 * dt * L1 # Apply boundary conditions. for bc in self.problem.get_bcs(t + dt): bc.apply(A, rhs) solver = KrylovSolver(krylov, preconditioner) solver.parameters['relative_tolerance'] = tol solver.parameters['absolute_tolerance'] = 0.0 solver.parameters['maximum_iterations'] = maxiter solver.parameters['monitor_convergence'] = verbose solver.set_operator(A) solver.solve(u1.vector(), rhs) return
def step(self, u1, u0, t, dt, tol=1.0e-10, maxiter=1000, verbose=True, krylov='gmres', preconditioner='ilu' ): L, b = self.problem.get_system(t + dt) A = self.M + dt * L v = TestFunction(self.problem.V) rhs = assemble(u0 * v * self.problem.dx_multiplier * self.problem.dx) rhs.axpy(dt, b) # Apply boundary conditions. for bc in self.problem.get_bcs(t + dt): bc.apply(A, rhs) solver = KrylovSolver(krylov, preconditioner) solver.parameters['relative_tolerance'] = tol solver.parameters['absolute_tolerance'] = 0.0 solver.parameters['maximum_iterations'] = maxiter solver.parameters['monitor_convergence'] = verbose P = self.problem.get_preconditioner(t + dt) if P: solver.set_operators(A, self.M + dt * P) else: solver.set_operator(A) solver.solve(u1.vector(), rhs) return
class Test(): def __init__(self): mesh = UnitSquareMesh(200, 200) self.V = FunctionSpace(mesh, 'Lagrange', 2) test, trial = TestFunction(self.V), TrialFunction(self.V) M = assemble(inner(test, trial) * dx) #self.M = assemble(inner(test, trial)*dx) self.solverM = KrylovSolver("cg", "amg") self.solverM.set_operator(M)
def solve_alpha_M_beta_F(self, alpha, beta, b, t): # Solve alpha * M * u + beta * F(u, t) = b for u. A = alpha * self.M + beta * self.A rhs = b - beta * self.b self.bcs.apply(A, rhs) solver = KrylovSolver('gmres', 'ilu') solver.parameters['relative_tolerance'] = 1.0e-13 solver.parameters['absolute_tolerance'] = 0.0 solver.parameters['maximum_iterations'] = 100 solver.parameters['monitor_convergence'] = True solver.set_operator(A) u = Function(self.V) solver.solve(u.vector(), rhs) return u
def solve_alpha_M_beta_F(self, alpha, beta, b, t): # Solve alpha * M * u + beta * F(u, t) = b for u. A = alpha * self.M + beta * self.A f.t = t v = TestFunction(self.V) b2 = assemble(f * v * dx) rhs = b.vector() - beta * b2 self.bcs.apply(A, rhs) solver = \ KrylovSolver('gmres', 'ilu') if alpha < 0.0 or beta > 0.0 \ else KrylovSolver('cg', 'amg') solver.parameters['relative_tolerance'] = 1.0e-13 solver.parameters['absolute_tolerance'] = 0.0 solver.parameters['maximum_iterations'] = 100 solver.parameters['monitor_convergence'] = False solver.set_operator(A) u = Function(self.V) solver.solve(u.vector(), rhs) return u
def stokes_solve( up_out, mu, u_bcs, p_bcs, f, dx=dx, verbose=True, tol=1.0e-10, maxiter=1000 ): # Some initial sanity checks. assert mu > 0.0 WP = up_out.function_space() # Translate the boundary conditions into the product space. new_bcs = [] for k, bcs in enumerate([u_bcs, p_bcs]): for bc in bcs: space = bc.function_space() C = space.component() if len(C) == 0: new_bcs.append(DirichletBC(WP.sub(k), bc.value(), bc.domain_args[0])) elif len(C) == 1: new_bcs.append(DirichletBC(WP.sub(k).sub(int(C[0])), bc.value(), bc.domain_args[0])) else: raise RuntimeError('Illegal number of subspace components.') # TODO define p*=-1 and reverse sign in the end to get symmetric system? # Define variational problem (u, p) = TrialFunctions(WP) (v, q) = TestFunctions(WP) r = Expression('x[0]', degree=1, domain=WP.mesh()) print("mu = %e" % mu) # build system a = mu * inner(r * grad(u), grad(v)) * 2 * pi * dx \ - ((r * v[0]).dx(0) + (r * v[1]).dx(1)) * p * 2 * pi * dx \ + ((r * u[0]).dx(0) + (r * u[1]).dx(1)) * q * 2 * pi * dx #- div(r*v)*p* 2*pi*dx \ #+ q*div(r*u)* 2*pi*dx L = inner(f, v) * 2 * pi * r * dx A, b = assemble_system(a, L, new_bcs) mode = 'lu' if mode == 'lu': solve(A, up_out.vector(), b, 'lu') elif mode == 'gmres': # For preconditioners for the Stokes system, see # # Fast iterative solvers for discrete Stokes equations; # J. Peters, V. Reichelt, A. Reusken. # prec = mu * inner(r * grad(u), grad(v)) * 2 * pi * dx \ - p * q * 2 * pi * r * dx P, btmp = assemble_system(prec, L, new_bcs) solver = KrylovSolver('tfqmr', 'amg') #solver = KrylovSolver('gmres', 'amg') solver.set_operators(A, P) solver.parameters['monitor_convergence'] = verbose solver.parameters['report'] = verbose solver.parameters['absolute_tolerance'] = 0.0 solver.parameters['relative_tolerance'] = tol solver.parameters['maximum_iterations'] = maxiter # Solve solver.solve(up_out.vector(), b) elif mode == 'fieldsplit': raise NotImplementedError('Fieldsplit solver not yet implemented.') # For an assortment of preconditioners, see # # Performance and analysis of saddle point preconditioners # for the discrete steady-state Navier-Stokes equations; # H.C. Elman, D.J. Silvester, A.J. Wathen; # Numer. Math. (2002) 90: 665-688; # <http://citeseerx.ist.psu.edu/viewdoc/summary?doi=10.1.1.145.3554>. # # Set up field split. W = SubSpace(WP, 0) P = SubSpace(WP, 1) u_dofs = W.dofmap().dofs() p_dofs = P.dofmap().dofs() prec = PETScPreconditioner() prec.set_fieldsplit([u_dofs, p_dofs], ['u', 'p']) PETScOptions.set('pc_type', 'fieldsplit') PETScOptions.set('pc_fieldsplit_type', 'additive') PETScOptions.set('fieldsplit_u_pc_type', 'lu') PETScOptions.set('fieldsplit_p_pc_type', 'jacobi') # Create Krylov solver with custom preconditioner. solver = PETScKrylovSolver('gmres', prec) solver.set_operator(A) return
class Magnitude(MetaField): """ Compute the magnitude of a Function-evaluated Field. Supports function spaces where all subspaces are equal. """ def before_first_compute(self, get): u = get(self.valuename) if isinstance(u, Function): if LooseVersion(dolfin_version()) > LooseVersion("1.6.0"): rank = len(u.ufl_shape) else: rank = u.rank() if rank == 0: self.f = Function(u.function_space()) elif rank >= 1: # Assume all subpaces are equal V = u.function_space().extract_sub_space([0]).collapse() mesh = V.mesh() el = V.ufl_element() self.f = Function(V) # Find out if we can operate directly on vectors, or if we have to use a projection # We can operate on vectors if all sub-dofmaps are ordered the same way # For simplicity, this is only tested for CG- or DG0-spaces # (this might always be true for these spaces, but better to be safe than sorry ) self.use_project = True if el.family() == "Lagrange" or (el.family() == "Discontinuous Lagrange" and el.degree() == 0): #dm = u.function_space().dofmap() dm0 = V.dofmap() self.use_project = False for i in xrange(u.function_space().num_sub_spaces()): Vi = u.function_space().extract_sub_space( [i]).collapse() dmi = Vi.dofmap() try: # For 1.6.0+ and newer diff = Vi.tabulate_dof_coordinates( ) - V.tabulate_dof_coordinates() except: # For 1.6.0 and older diff = dmi.tabulate_all_coordinates( mesh) - dm0.tabulate_all_coordinates(mesh) if len(diff) > 0: max_diff = max(abs(diff)) else: max_diff = 0.0 max_diff = MPI.max(mpi_comm_world(), max_diff) if max_diff > 1e-12: self.use_project = True break self.assigner = FunctionAssigner( [V] * u.function_space().num_sub_spaces(), u.function_space()) self.subfuncs = [ Function(V) for _ in range(u.function_space().num_sub_spaces()) ] # IF we have to use a projection, build projection matrix only once if self.use_project: self.v = TestFunction(V) M = assemble(inner(self.v, TrialFunction(V)) * dx) self.projection = KrylovSolver("cg", "default") self.projection.set_operator(M) elif isinstance(u, Iterable) and all( isinstance(_u, Number) for _u in u): pass elif isinstance(u, Number): pass else: # Don't know how to handle object cbc_warning( "Don't know how to calculate magnitude of object of type %s." % type(u)) def compute(self, get): u = get(self.valuename) if isinstance(u, Function): if not hasattr(self, "use_project"): self.before_first_compute(get) if LooseVersion(dolfin_version()) > LooseVersion("1.6.0"): rank = len(u.ufl_shape) else: rank = u.rank() if rank == 0: self.f.vector().zero() self.f.vector().axpy(1.0, u.vector()) self.f.vector().abs() return self.f elif rank >= 1: if self.use_project: b = assemble(sqrt(inner(u, u)) * self.v * dx(None)) self.projection.solve(self.f.vector(), b) else: self.assigner.assign(self.subfuncs, u) self.f.vector().zero() for i in xrange(u.function_space().num_sub_spaces()): vec = self.subfuncs[i].vector() vec.apply('') self.f.vector().axpy(1.0, vec * vec) try: sqrt_in_place(self.f.vector()) except: r = self.f.vector().local_range() self.f.vector()[r[0]:r[1]] = np.sqrt( self.f.vector()[r[0]:r[1]]) self.f.vector().apply('') return self.f elif isinstance(u, Iterable) and all( isinstance(_u, Number) for _u in u): return np.sqrt(sum(_u**2 for _u in u)) elif isinstance(u, Number): return abs(u) else: # Don't know how to handle object cbc_warning( "Don't know how to calculate magnitude of object of type %s. Returning object." % type(u)) return u
class AcousticWave(): """ Solution of forward and adjoint equations for acoustic inverse problem """ def __init__(self, functionspaces_V): """ Input: functionspaces_V = dict containing functionspaces for state/adj ('V') and med param ('Vl' for lambda and 'Vr' for rho) """ self.readV(functionspaces_V) self.verbose = False # print info self.lump = False # Lump the mass matrix self.lumpD = False # Lump the ABC matrix self.timestepper = None # 'backward', 'centered' self.exact = None # exact solution at final time self.u0init = None # provides u(t=t0) self.utinit = None # provides u_t(t=t0) self.u1init = None # provides u1 = u(t=t0+/-Dt) self.bc = None self.abc = False self.ftime = lambda t: 0.0 # ftime(tt) = src term @ t=tt (in np.array()) self.set_fwd() # default is forward problem def copy(self): """(hard) copy constructor""" newobj = self.__class__({'V':self.V, 'Vl':self.Vl, 'Vr':self.Vr}) newobj.lump = self.lump newobj.timestepper = self.timestepper newobj.exact = self.exact newobj.utinit = self.utinit newobj.u1init = self.u1init newobj.bc = self.bc if self.abc == True: newobj.set_abc(self.V.mesh(), self.class_bc_abc, self.lumpD) newobj.ftime = self.ftime newobj.update({'lambda':self.lam, 'rho':self.rho, \ 't0':self.t0, 'tf':self.tf, 'Dt':self.Dt, \ 'u0init':self.u0init, 'utinit':self.utinit, 'u1init':self.u1init}) return newobj def readV(self, functionspaces_V): # Solutions: self.V = functionspaces_V['V'] self.test = TestFunction(self.V) self.trial = TrialFunction(self.V) self.u0 = Function(self.V) # u(t-Dt) self.u1 = Function(self.V) # u(t) self.u2 = Function(self.V) # u(t+Dt) self.rhs = Function(self.V) self.sol = Function(self.V) # Parameters: self.Vl = functionspaces_V['Vl'] self.lam = Function(self.Vl) self.Vr = functionspaces_V['Vr'] self.rho = Function(self.Vr) if functionspaces_V.has_key('Vm'): self.Vm = functionspaces_V['Vm'] self.mu = Function(self.Vm) self.elastic = True assert(False) else: self.elastic = False self.weak_k = inner(self.lam*nabla_grad(self.trial), \ nabla_grad(self.test))*dx self.weak_m = inner(self.rho*self.trial,self.test)*dx def set_abc(self, mesh, class_bc_abc, lumpD=False): self.abc = True # False means zero-Neumann all-around if lumpD: self.lumpD = True abc_boundaryparts = FacetFunction("size_t", mesh) class_bc_abc.mark(abc_boundaryparts, 1) self.ds = Measure("ds")[abc_boundaryparts] self.weak_d = inner(sqrt(self.lam*self.rho)*self.trial, self.test)*self.ds(1) self.class_bc_abc = class_bc_abc # to make copies def set_fwd(self): self.fwdadj = 1.0 self.ftime = None def set_adj(self): self.fwdadj = -1.0 self.ftime = None def get_tt(self, nn): if self.fwdadj > 0.0: return self.times[nn] else: return self.times[-nn-1] def update(self, parameters_m): assert not self.timestepper == None, "You need to set a time stepping method" # Time options: if parameters_m.has_key('t0'): self.t0 = parameters_m['t0'] if parameters_m.has_key('tf'): self.tf = parameters_m['tf'] if parameters_m.has_key('Dt'): self.Dt = parameters_m['Dt'] if parameters_m.has_key('t0') or parameters_m.has_key('tf') or parameters_m.has_key('Dt'): self.Nt = int(round((self.tf-self.t0)/self.Dt)) self.Tf = self.t0 + self.Dt*self.Nt self.times = np.linspace(self.t0, self.Tf, self.Nt+1) assert isequal(self.times[1]-self.times[0], self.Dt, 1e-16), "Dt modified" self.Dt = self.times[1] - self.times[0] assert isequal(self.Tf, self.tf, 1e-2), "Final time differs by more than 1%" if not isequal(self.Tf, self.tf, 1e-12): print 'Final time modified from {} to {} ({}%)'.\ format(self.tf, self.Tf, abs(self.Tf-self.tf)/self.tf) # Initial conditions: if parameters_m.has_key('u0init'): self.u0init = parameters_m['u0init'] if parameters_m.has_key('utinit'): self.utinit = parameters_m['utinit'] if parameters_m.has_key('u1init'): self.u1init = parameters_m['u1init'] if parameters_m.has_key('um1init'): self.um1init = parameters_m['um1init'] # Medium parameters: setfct(self.lam, parameters_m['lambda']) if self.verbose: print 'lambda updated ' if self.elastic == True: setfct(self.mu, parameters_m['mu']) if self.verbose: print 'mu updated' if self.verbose: print 'assemble K', self.K = assemble(self.weak_k) if self.verbose: print ' -- K assembled' if parameters_m.has_key('rho'): setfct(self.rho, parameters_m['rho']) # Mass matrix: if self.verbose: print 'rho updated\nassemble M', Mfull = assemble(self.weak_m) if self.lump: self.solverM = LumpedMatrixSolverS(self.V) self.solverM.set_operator(Mfull, self.bc) self.M = self.solverM else: if mpisize == 1: self.solverM = LUSolver() self.solverM.parameters['reuse_factorization'] = True self.solverM.parameters['symmetric'] = True else: self.solverM = KrylovSolver('cg', 'amg') self.solverM.parameters['report'] = False self.M = Mfull if not self.bc == None: self.bc.apply(Mfull) self.solverM.set_operator(Mfull) if self.verbose: print ' -- M assembled' # Matrix D for abs BC if self.abc == True: if self.verbose: print 'assemble D', Mfull = assemble(self.weak_m) Dfull = assemble(self.weak_d) if self.lumpD: self.D = LumpedMatrixSolverS(self.V) self.D.set_operator(Dfull, None, False) if self.lump: self.solverMplD = LumpedMatrixSolverS(self.V) self.solverMplD.set_operators(Mfull, Dfull, .5*self.Dt, self.bc) self.MminD = LumpedMatrixSolverS(self.V) self.MminD.set_operators(Mfull, Dfull, -.5*self.Dt, self.bc) else: self.D = Dfull if self.verbose: print ' -- D assembled' else: self.D = 0.0 #@profile def solve(self): """ General solver method """ if self.timestepper == 'backward': def iterate(tt): self.iteration_backward(tt) elif self.timestepper == 'centered': def iterate(tt): self.iteration_centered(tt) else: print "Time stepper not implemented" sys.exit(1) if self.verbose: print 'Compute solution' solout = [] # Store computed solution # u0: tt = self.get_tt(0) if self.verbose: print 'Compute solution -- time {}'.format(tt) setfct(self.u0, self.u0init) solout.append([self.u0.vector().array(), tt]) # Compute u1: if not self.u1init == None: self.u1 = self.u1init else: assert(not self.utinit == None) setfct(self.rhs, self.ftime(tt)) self.rhs.vector().axpy(-self.fwdadj, self.D*self.utinit.vector()) self.rhs.vector().axpy(-1.0, self.K*self.u0.vector()) if not self.bc == None: self.bc.apply(self.rhs.vector()) self.solverM.solve(self.sol.vector(), self.rhs.vector()) setfct(self.u1, self.u0) self.u1.vector().axpy(self.fwdadj*self.Dt, self.utinit.vector()) self.u1.vector().axpy(0.5*self.Dt**2, self.sol.vector()) tt = self.get_tt(1) if self.verbose: print 'Compute solution -- time {}'.format(tt) solout.append([self.u1.vector().array(), tt]) # Iteration for nn in xrange(2, self.Nt+1): iterate(tt) # Advance to next time step setfct(self.u0, self.u1) setfct(self.u1, self.u2) tt = self.get_tt(nn) if self.verbose: print 'Compute solution -- time {}, rhs {}'.\ format(tt, np.max(np.abs(self.ftime(tt)))) solout.append([self.u1.vector().array(),tt]) if self.fwdadj > 0.0: assert isequal(tt, self.Tf, 1e-16), \ 'tt={}, Tf={}, reldiff={}'.format(tt, self.Tf, abs(tt-self.Tf)/self.Tf) else: assert isequal(tt, self.t0, 1e-16), \ 'tt={}, t0={}, reldiff={}'.format(tt, self.t0, abs(tt-self.t0)) return solout, self.computeerror() def iteration_centered(self, tt): setfct(self.rhs, (self.Dt**2)*self.ftime(tt)) self.rhs.vector().axpy(-1.0, self.MminD*self.u0.vector()) self.rhs.vector().axpy(2.0, self.M*self.u1.vector()) self.rhs.vector().axpy(-self.Dt**2, self.K*self.u1.vector()) if not self.bc == None: self.bc.apply(self.rhs.vector()) self.solverMplD.solve(self.u2.vector(), self.rhs.vector()) def iteration_backward(self, tt): setfct(self.rhs, self.Dt*self.ftime(tt)) self.rhs.vector().axpy(-self.Dt, self.K*self.u1.vector()) self.rhs.vector().axpy(-1.0, self.D*(self.u1.vector()-self.u0.vector())) if not self.bc == None: self.bc.apply(self.rhs.vector()) self.solverM.solve(self.sol.vector(), self.rhs.vector()) setfct(self.u2, 2.0*self.u1.vector()) self.u2.vector().axpy(-1.0, self.u0.vector()) self.u2.vector().axpy(self.Dt, self.sol.vector()) def computeerror(self): return self.computerelativeerror() def computerelativeerror(self): if not self.exact == None: #MM = assemble(inner(self.trial, self.test)*dx) MM = self.M norm_ex = np.sqrt(\ (MM*self.exact.vector()).inner(self.exact.vector())) diff = self.exact.vector() - self.u1.vector() if norm_ex > 1e-16: return np.sqrt((MM*diff).inner(diff))/norm_ex else: return np.sqrt((MM*diff).inner(diff)) else: return [] def computeabserror(self): if not self.exact == None: MM = self.M diff = self.exact.vector() - self.u1.vector() return np.sqrt((MM*diff).inner(diff)) else: return []