M = solver.createMassMatrix(grid, np.ones(grid.cellCount())) ut = pg.RVector(dof, 0.0) rhs = pg.RVector(dof, 0.0) b = pg.RVector(dof, 0.0) theta = 0 boundUdir = solver.parseArgToBoundaries(dirichletBC, grid) for n in range(1, len(times)): b = (M - A * dt) * u[n - 1] + rhs * dt S = M solver.assembleDirichletBC(S, boundUdir, rhs=b) solve = pg.LinSolver(S) solve.solve(b, ut) u[n, :] = ut plt.plot(times, u[:, probeID], label='Explicit Euler') theta = 1 for n in range(1, len(times)): b = (M + A * (dt*(theta - 1.0))) * u[n-1] + \ rhs * (dt*(1.0 - theta)) + \ rhs * dt * theta b = M * u[n - 1] + rhs * dt
def solveFiniteVolume(mesh, a=1.0, b=0.0, f=0.0, fn=0.0, vel=None, u0=0.0, times=None, uL=None, relax=1.0, ws=None, scheme='CDS', **kwargs): """Calculate for u. NOTE works only for steady boundary conditions!!! !!Refactor with solver class and Runga-Kutte solver!! Parameters ---------- mesh : :gimliapi:`GIMLI::Mesh` Mesh represents spatial discretization of the calculation domain a : value | array | callable(cell, userData) cell values b : value | array | callable(cell, userData) TODO What is b f : iterable(cell) TODO What is f fn : iterable(cell) TODO What is fn vel : ndarray (N,dim) | RMatrix(N,dim) velocity field [[v_j,]_i,] with i=[1..3] for the mesh dimension and j = [0 .. N-1] with N either Amount of Cells, Nodes or Boundaries. Velocity per boundary is preferred. u0 : value | array | callable(cell, userData) Starting field ws : Workspace This can be an empty class that will used as an Workspace to store and cache data. If ws is given: The system matrix is taken from ws or calculated once and stored in ws for further usage. The system matrix is cached in this Workspace as ws.S The LinearSolver with the factorized matrix is cached in this Workspace as ws.solver The rhs vector is only stored in this Workspace as ws.rhs scheme : str [CDS] Finite volume scheme: :py:mod:`pygimli.solver.diffusionConvectionKernel` **kwargs: * uB : Dirichlet boundary conditions * duB : Neumann boundary conditions Returns ------- u : ndarray(nTimes, nCells) solution field for all time steps """ verbose = kwargs.pop('verbose', False) # The Workspace is to hold temporary data or preserve matrix rebuild # swatch = pg.Stopwatch(True) sparse = True workspace = pg.solver.WorkSpace() if ws: workspace = ws a = pg.solver.parseArgToArray(a, [mesh.cellCount(), mesh.boundaryCount()]) f = pg.solver.parseArgToArray(f, mesh.cellCount()) fn = pg.solver.parseArgToArray(fn, mesh.cellCount()) boundsDirichlet = None boundsNeumann = None # BEGIN check for Courant-Friedrichs-Lewy if vel is not None: if isinstance(vel, float): print("Warning! .. velocity is float and no vector field") if len(vel) != mesh.cellCount() and \ len(vel) != mesh.nodeCount() and \ len(vel) != mesh.boundaryCount(): print("mesh:", mesh) print("vel:", vel.shape) raise BaseException("Velocity field has wrong dimension.") if len(vel) is not mesh.nodeCount(): vel = pg.meshtools.nodeDataToBoundaryData(mesh, vel) vmax = 0 if mesh.dimension() == 3: vmax = np.max(np.sqrt(vel[:, 0]**2 + vel[:, 1]**2 + vel[:, 2]**2)) else: vmax = np.max(np.sqrt(vel[:, 0]**2 + vel[:, 1]**2)) dt = times[1] - times[0] dx = min(mesh.boundarySizes()) c = vmax * dt / dx if c > 1: print( "Courant-Friedrichs-Lewy Number:", c, "but sould be lower 1 to ensure movement inside a cell " "per timestep. (" "vmax =", vmax, "dt =", dt, "dx =", dx, "dt <", dx / vmax, " | N > ", int((times[-1] - times[0]) / (dx / vmax)) + 1, ")") # END check for Courant-Friedrichs-Lewy if not hasattr(workspace, 'S'): if 'uB' in kwargs: boundsDirichlet = pg.solver.parseArgToBoundaries( kwargs['uB'], mesh) if 'duB' in kwargs: boundsNeumann = pg.solver.parseArgToBoundaries(kwargs['duB'], mesh) workspace.S, workspace.rhsBCScales = diffusionConvectionKernel( mesh=mesh, a=a, b=b, uB=boundsDirichlet, duB=boundsNeumann, # u0=u0, fn=fn, vel=vel, scheme=scheme, sparse=sparse, userData=kwargs.pop('userData', None)) dof = len(workspace.rhsBCScales) workspace.ap = np.zeros(dof) # for nonlinears if uL is not None: for i in range(dof): val = 0.0 if sparse: val = workspace.S.getVal(i, i) / relax workspace.S.setVal(i, i, val) else: val = workspace.S[i, i] / relax workspace.S[i, i] = val workspace.ap[i] = val # print('FVM kernel 2:', swatch.duration(True)) # endif: not hasattr(workspace, 'S'): workspace.rhs = np.zeros(len(workspace.rhsBCScales)) workspace.rhs[0:mesh.cellCount()] = f # * mesh.cellSizes() workspace.rhs += workspace.rhsBCScales # for nonlinear: relax progress with scaled last result if uL is not None: workspace.rhs += (1. - relax) * workspace.ap * uL # print('FVM: Prep:', swatch.duration(True)) if not hasattr(times, '__len__'): if sparse and not hasattr(workspace, 'solver'): Sm = pg.RSparseMatrix(workspace.S) # hold Sm until we have reference counting, # loosing Sm here will kill LinSolver later workspace.Sm = Sm workspace.solver = pg.LinSolver(Sm, verbose=verbose) u = None if sparse: u = workspace.solver.solve(workspace.rhs) else: u = np.linalg.solve(workspace.S, workspace.rhs) # print('FVM solve:', swatch.duration(True)) return u[0:mesh.cellCount():1] else: theta = kwargs.pop('theta', 0.5 + 1e-6) if sparse: I = pg.solver.identity(len(workspace.rhs)) else: I = np.diag(np.ones(len(workspace.rhs))) if verbose: print("Solve timesteps with Crank-Nicolson.") return pg.solver.crankNicolson(times, theta, workspace.S, I, f=workspace.rhs, u0=pg.solver.parseArgToArray( u0, mesh.cellCount(), mesh), verbose=verbose)
def solveFiniteElements(mesh, a=1.0, b=0.0, f=0.0, times=None, userData=None, verbose=False, stats=None, **kwargs): r"""Solve partial differential equation with Finite Elements. This function is a syntactic sugar proxy for using the Finite Element functionality of the library core to solve elliptic and parabolic partial differential of the following type: .. math:: \frac{\partial u}{\partial t} & = \nabla\cdot(a \nabla u) + b u + f(\mathbf{r},t) \\ u(\mathbf{r}, t) & = u_B \quad\mathbf{r}\in\Gamma_{\text{Dirichlet}}\\ \frac{\partial u(\mathbf{r}, t)}{\partial \mathbf{n}} & = u_{\partial \text{B}} \quad\mathbf{r}\in\Gamma_{\text{Neumann}}\\ u(\mathbf{r}, t=0) & = u_0 \quad\text{with} \quad\mathbf{r}\in\Omega\quad\text{for}\quad t\neq 0 The Domain :math:`\Omega` and the Boundary :math:`\Gamma` are defined through the given mesh with appropriate boundary marker. The solution :math:`u(\mathbf{r}, t)` is given for each node in mesh. TODO: * unsteady ub and dub Parameters ---------- mesh : :gimliapi:`GIMLI::Mesh` Mesh represents spatial discretization of the calculation domain a : value | array | callable(cell, userData) Cell values b : value | array | callable(cell, userData) Cell values u0 : value | array | callable(pos, userData) Node values uB : value | array | callable(pos, userData) Dirichlet values for u at the boundary uN : list([node, value]) Dirichlet values for u at given nodes duB : value | array | callable(pos, userData) Neumann values for du/dn at the boundary f : value | array(cells) | array(nodes) | callable(args, kwargs) force values times : array [None] Solve as time dependent problem for the given times. theta : float [1] - :math:`theta = 0` means explicit Euler, maybe stable for :math:`\Delta t \quad\text{near}\quad h` - :math:`theta = 0.5`, Crank-Nicolsen, maybe instable - :math:`theta = 1`, implicit Euler If unsure choose :math:`\theta = 0.5 + \epsilon`, which is probably stable. progress : bool Give some calculation progress. ret : Workspace for results so no new memory will be allocated. Returns ------- u : array Returns the solution u either 1,n array for stationary problems or for m,n array for m time steps Examples -------- >>> import pygimli as pg >>> from pygimli.meshtools import polytools as plc >>> from pygimli.mplviewer import drawField, drawMesh >>> import matplotlib.pyplot as plt >>> world = plc.createWorld(start=[-10, 0], end=[10, -10], ... marker=1, worldMarker=False) >>> c1 = plc.createCircle(pos=[0.0, -5.0], radius=3.0, area=.1, marker=2) >>> mesh = pg.meshtools.createMesh([world, c1], quality=34.3) >>> u = pg.solver.solveFiniteElements(mesh, a=[[1, 100], [2, 1]], ... uB=[[4, 1.0], [2, 0.0]]) >>> fig, ax = plt.subplots() >>> pc = drawField(ax, mesh, u) >>> drawMesh(ax, mesh) >>> plt.show() See Also -------- other solver TODO """ if 'uDirichlet' in kwargs or 'uBoundary' in kwargs: raise BaseException("use uB instead") if 'uBoundary' in kwargs: raise BaseException("use duB instead") debug = kwargs.pop('debug', False) mesh.createNeighbourInfos() if verbose: print("Mesh: ", str(mesh)) dof = mesh.nodeCount() swatch = pg.Stopwatch(True) swatch2 = pg.Stopwatch(True) # check for material parameter a = parseArgToArray(a, ndof=mesh.cellCount(), mesh=mesh, userData=userData) b = parseArgToArray(b, ndof=mesh.cellCount(), mesh=mesh, userData=userData) if debug: print("2: ", swatch2.duration(True)) # assemble the stiffness matrix A = createStiffnessMatrix(mesh, a) if debug: print("3: ", swatch2.duration(True)) M = createMassMatrix(mesh, b) if debug: print("4: ", swatch2.duration(True)) S = A + M if debug: print("5: ", swatch2.duration(True)) if times is None: rhs = assembleForceVector(mesh, f, userData=userData) if debug: print("6a: ", swatch2.duration(True)) if 'duB' in kwargs: assembleNeumannBC(S, parseArgToBoundaries(kwargs['duB'], mesh), time=0.0, userData=userData) if debug: print("6b: ", swatch2.duration(True)) if 'uB' in kwargs: assembleDirichletBC(S, parseArgToBoundaries(kwargs['uB'], mesh), rhs, time=0.0, userData=userData) if 'uN' in kwargs: assembleDirichletBC(S, [], nodePairs=kwargs['uN'], rhs=rhs, time=0.0, userData=userData) if debug: print("6c: ", swatch2.duration(True)) # create result array u = kwargs.pop('ret', None) singleForce = True if hasattr(rhs, 'ndim'): if rhs.ndim == 2: singleForce = False if u is None: u = np.zeros(rhs.shape) else: if isinstance(a[0], complex): if u is None: u = pg.CVector(rhs.size(), 0.0) rhs = pg.toComplex(rhs) else: if u is None: u = pg.RVector(rhs.size(), 0.0) assembleTime = swatch.duration(True) if stats: stats.assembleTime = assembleTime if verbose: print(("Asssemblation time: ", assembleTime)) # showSparseMatrix(S) solver = pg.LinSolver(False) solver.setMatrix(S, 0) if singleForce: u = solver.solve(rhs) else: for i, r in enumerate(rhs): solver.solve(r, u[i]) solverTime = swatch.duration(True) if verbose: if stats: stats.solverTime = solverTime print(("Solving time: ", solverTime)) return u else: # times given pg.solver.checkCFL(times, mesh, max(a)) if debug: print("start TL", swatch.duration()) M = createMassMatrix(mesh) F = assembleForceVector(mesh, f) if 'u0' in kwargs: u0 = parseArgToArray(kwargs['u0'], dof, mesh, userData) progress = None if 'progress' in kwargs: from pygimli.utils import ProgressBar progress = ProgressBar(its=len(times), width=40, sign='+') theta = kwargs.pop('theta', 1.0) dynamic = kwargs.pop('dynamic', False) if not dynamic: A = createStiffnessMatrix(mesh, a) if 'uB' in kwargs: assembleDirichletBC(A, parseArgToBoundaries(kwargs['uB'], mesh), rhs=F) if 'uN' in kwargs: assembleDirichletBC(A, [], nodePairs=kwargs['uN'], rhs=F) if 'duB' in kwargs: assembleNeumannBC(A, parseArgToBoundaries(kwargs['duB'], mesh), time=None, userData=userData) return crankNicolson(times, theta, A, M, F, u0=u0, progress=progress) rhs = np.zeros((len(times), dof)) # rhs kann zeitabhängig sein ..wird hier nicht berücksichtigt rhs[:] = F # this is slow: optimize if debug: print("rhs", swatch.duration()) U = np.zeros((len(times), dof)) U[0, :] = u0 # init state u = pg.RVector(dof, 0.0) if debug: print("u0", swatch.duration()) measure = 0. for n in range(1, len(times)): swatch.reset() dt = times[n] - times[n - 1] # previous timestep # print "i: ", i, dt, U[i - 1] if 'duB' in kwargs: # aufschreiben und checken ob neumann auf A oder auf S mit # skaliertem val*dt angewendet wird A = createStiffnessMatrix(mesh, a) assembleNeumannBC(A, parseArgToBoundaries(kwargs['duB'], mesh), time=times[n], userData=userData) swatch.reset() # (A + a*B)u is fastest, # followed by A*u + (B*u)*a and finally A*u + a*B*u and b = (M + (dt * (theta - 1.)) * A) * U[n - 1] + \ dt * ((1.0 - theta) * rhs[n - 1] + theta * rhs[n]) # print ('a',swatch.duration(True)) # b = M * U[n - 1] - (A * U[n - 1]) * (dt*(1.0 - theta)) + \ # dt * ((1.0 - theta) * rhs[n - 1] + theta * rhs[n]) # print ('b',swatch.duration(True)) # b = M * U[n - 1] - (dt*(1.0 - theta)) * A * U[n - 1] + \ # dt * ((1.0 - theta) * rhs[n - 1] + theta * rhs[n]) # print ('c',swatch.duration(True)) measure += swatch.duration() S = M + A * dt * theta if 'uB' in kwargs: assembleDirichletBC(S, parseArgToBoundaries(kwargs['uB'], mesh), rhs=b, time=times[n], userData=userData) if 'uN' in kwargs: assembleDirichletBC(S, [], nodePairs=kwargs['uN'], rhs=b, time=times[n], userData=userData) # u = S/b t_prep = swatch.duration(True) solver = pg.LinSolver(S, verbose) solver.solve(b, u) if 'plotTimeStep' in kwargs: kwargs['plotTimeStep'](u, times[n]) U[n, :] = np.asarray(u) if progress: progress.update(n, ' t_prep: ' + str(round(t_prep, 5)) + 's' + \ ' t_step: ' + str(round(swatch.duration(),5)) + 's') if debug: print("Measure(" + str(len(times)) + "): ", measure, measure / len(times)) return U
def crankNicolson(times, theta, S, I, f, u0=None, progress=None, debug=None): """ S = constant over time f = constant over time """ if len(times) < 2: raise BaseException("We need at least 2 times for Crank " "Nicolsen time discretization." + str(len(times))) sw = pg.Stopwatch(True) if u0 is None: u0 = np.zeros(len(f)) u = np.zeros((len(times), len(f))) u[0, :] = u0 dt = times[1] - times[0] rhs = np.zeros((len(times), len(f))) rhs[:] = f A = I + S * dt * theta solver = pg.LinSolver(A, verbose=False) timeAssemble = [] timeSolve = [] # print('0', min(u[0]), max(u[0]), min(f), max(f)) timeMeasure = False if progress: timeMeasure = True for n in range(1, len(times)): if timeMeasure: pg.tic() # pg.tic() #bRef = (I + (dt * (theta - 1.)) * S) * u[n - 1] + \ #dt * ((1.0 - theta) * rhs[n - 1] + theta * rhs[n]) # pg.toc() # # pg.tic() #b = I * u[n - 1] + ((dt * (theta - 1.)) * S) * u[n - 1] + \ #dt * ((1.0 - theta) * rhs[n - 1] + theta * rhs[n]) # pg.toc() # # pg.tic() b = I * u[n - 1] + S.mult(dt * (theta - 1.) * u[n - 1]) + \ dt * ((1.0 - theta) * rhs[n - 1] + theta * rhs[n]) # pg.toc() # # print(np.linalg.norm(b-b1)) #np.testing.assert_allclose(bRef, b) if timeMeasure: timeAssemble.append(pg.dur()) if timeMeasure: pg.tic() u[n, :] = solver.solve(b) if timeMeasure: timeSolve.append(pg.dur()) # A = (I + dt * theta * S) # u[n, : ] = linsolve(A, b) if progress: progress.update(n, ' t_prep: ' + str(round(timeAssemble[-1], 5)) + 's' + \ ' t_step: ' + str(round(timeSolve[-1], 5)) + 's') #if verbose and (n % verbose == 0): ## print(min(u[n]), max(u[n])) #print("timesteps:", n, "/", len(times), #'runtime:', sw.duration(), "s", #'assemble:', np.mean(timeAssemble), #'solve:', np.mean(timeSolve)) return u
def solvePressureWave(mesh, velocities, times, sourcePos, uSource, verbose): r""" Solve pressure wave equation. Solve pressure wave for a given source function .. math:: \frac{\partial^2 u}{\partial t^2} & = \diverg(a\grad u) + f\\ finalize equation Parameters ---------- mesh : :gimliapi:`GIMLI::Mesh` Mesh to solve on velocities : array velocities for each cell of the mesh time : array Time base definition sourcePos : RVector3 Source position uSource : array u(t, sourcePos) source movement of length(times) Usually a Ricker wavelet of the desired seismic signal frequency. Returns ------- u : RMatrix Return Examples -------- See TODO write example """ A = pg.RSparseMatrix() M = pg.RSparseMatrix() # F = pg.RVector(mesh.nodeCount(), 0.0) rhs = pg.RVector(mesh.nodeCount(), 0.0) u = pg.RMatrix(len(times), mesh.nodeCount()) v = pg.RMatrix(len(times), mesh.nodeCount()) sourceID = mesh.findNearestNode(sourcePos) if len(uSource) != len(times): raise Exception("length of uSource does not fit length of times: " + str(uSource) + " != " + len(times)) A.fillStiffnessMatrix(mesh, velocities * velocities) M.fillMassMatrix(mesh) # M.fillMassMatrix(mesh, velocities) FV = 0 if FV: A, rhs = pygimli.solver.diffusionConvectionKernel(mesh, velocities * velocities, sparse=1) M = pygimli.solver.identity(len(rhs)) u = pg.RMatrix(len(times), mesh.cellCount()) v = pg.RMatrix(len(times), mesh.cellCount()) sourceID = mesh.findCell(sourcePos).id() dt = times[1] - times[0] theta = 0.51 S1 = M + dt * dt * theta * theta * A S2 = M solver1 = pg.LinSolver(S1, verbose=False) solver2 = pg.LinSolver(S2, verbose=False) swatch = pg.Stopwatch(True) # ut = pg.RVector(mesh.nodeCount(), .0) # vt = pg.RVector(mesh.nodeCount(), .0) timeIter1 = np.zeros(len(times)) timeIter2 = np.zeros(len(times)) timeIter3 = np.zeros(len(times)) timeIter4 = np.zeros(len(times)) for n in range(1, len(times)): u[n - 1, sourceID] = uSource[n - 1] # solve for u tic = time.time() # + * dt*dt * F rhs = dt * M * v[n - 1] + (M - dt * dt * theta * (1. - theta) * A) * u[n - 1] timeIter1[n - 1] = time.time() - tic tic = time.time() u[n] = solver1.solve(rhs) timeIter2[n - 1] = time.time() - tic # solve for v tic = time.time() rhs = M * v[n - 1] - dt * \ ((1. - theta) * A * u[n - 1] + theta * A * u[n]) # + dt * F timeIter3[n - 1] = time.time() - tic tic = time.time() v[n] = solver2.solve(rhs) timeIter4[n - 1] = time.time() - tic # same as above # rhs = M * v[n-1] - dt * A * u[n-1] + dt * F # v[n] = solver1.solve(rhs) t1 = swatch.duration(True) if verbose and (n % verbose == 0): print( str(n) + "/" + str(len(times)), times[n], dt, t1, min(u[n]), max(u[n])) # plt.figure() # plt.plot(timeIter1, label='Ass:1') # plt.plot(timeIter2, label='Sol:1') # plt.plot(timeIter3, label='Ass:2') # plt.plot(timeIter4, label='Sol:2') # plt.legend() # plt.figure() return u