def isComplex(vals): """Check numpy or pg.Vector if have complex data type""" z = pg.CVector(2) if len(vals) > 0: if hasattr(vals, '__iter__'): if isinstance(vals[0], np.complex) or \ isinstance(vals[0], complex): return True return False
def test_CVectorToNumpy(self): """Implemented through hand_made_wrapper.py""" # check ob wirklich from array genommen wird! v = pg.CVector(10, 1.1 + 1j * 3) a = np.array(v) self.assertEqual(type(a), np.ndarray) self.assertEqual(a.dtype, np.complex) self.assertEqual(len(a), 10) self.assertEqual(a[0], 1.1 + 1j * 3)
def test_XVectorBasics(self): def testVector(v): for t in v: self.assertEqual(t, 1.0) self.assertEqual(sum(v), 5.0) self.assertFalse(pg.haveInfNaN(v)) v[1] = 0 self.assertEqual(v[1], 0) v[1] = 1 self.assertEqual(v[1], 1) #print(v/v) testVector(pg.RVector(5, 1.0)) testVector(pg.CVector(5, 1.0)) testVector(pg.BVector(5, True)) testVector(pg.IVector(5, 1))
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: 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) theta = kwargs.pop('theta', 1) if 'duB' not in kwargs: 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) return crankNicolson(times, theta, A, M, F, u0=u0, verbose=verbose) 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' in kwargs: if kwargs['progress']: print(("\t" + str(n) + "/" + str(len(times) - 1) + ": " + str(t_prep) + "/" + str(swatch.duration()))) if debug: print("Measure(" + str(len(times)) + "): ", measure, measure / len(times)) return U