def assembleNeumannBC(S, boundaryPairs, time=0.0, userData=None): r"""Apply Neumann condition to the system matrix S. Apply Neumann condition to the system matrix S. .. math:: \frac{\partial u(\arr{r}, t)}{\partial\textbf{n}} = \textbf{n}\grad u(\arr{r}, t) = g \quad\text{with}\quad\arr{r} \quad\text{on}\quad \partial\Omega Parameters ---------- S : :gimliapi:`GIMLI::RSparseMatrix` System matrix of the system equation. boundaryPair : list() List of pairs [ :gimliapi:`GIMLI::Boundary`, g ]. The value g will assigned to the nodes of the boundaries. Later assignment overwrites prior. :math:`g` need to be a scalar value (float or int) or a value generator function that will be executed at runtime. See :py:mod:`pygimli.solver.solver.parseArgToBoundaries` See tutorial section for an example, e.g., Modelling with Boundary Conditions time : float Will be forwarded to value generator. userData : class Will be forwarded to value generator. """ Se = pg.ElementMatrix() if not hasattr(boundaryPairs, '__getitem__'): raise BaseException("Boundary pairs need to be a list of " "[boundary, value]") for pair in boundaryPairs: boundary = pair[0] val = pair[1] g = generateBoundaryValue(boundary, val, time, userData) if g is not 0.0 and g is not None: Se.u2(boundary) Se *= g S += Se
def createStiffnessMatrix(mesh, a=None): """Create the Stiffness matrix. Calculates the scaled stiffness matrix for the given mesh scaled with the per cell values a. ..math:: ... Parameters ---------- mesh : :gimliapi:`GIMLI::Mesh` Arbitrary mesh to calculate the stiffness for. Type of base and shape functions depends on the cell types. a : array, either complex or real, callable Per cell values., e.g., physical parameter. If None given default is 1. Returns ------- A : :gimliapi:`GIMLI::RSparseMatrix` Stiffness matrix """ if a is None: a = pg.RVector(mesh.cellCount(), 1.0) A = None if isinstance(a[0], float) or isinstance(a[0], np.float64): A = pg.RSparseMatrix() A.fillStiffnessMatrix(mesh, a) return A else: A = pg.CSparseMatrix() # create matrix structure regarding the mesh A.buildSparsityPattern(mesh) # define a local element matrix A_l = pg.ElementMatrix() for c in mesh.cells(): A_l.ux2uy2uz2(c) A.add(A_l, scale=a[c.id()]) # if c.id() == 0: # print(c.id(), A_l) return A
def createJacobian(self, model): """Create Jacobian matrix.""" if self.subPotentials is None: self.response(model) J = self.jacobian() J.resize(self.data.size(), self.regionManager().parameterCount()) cells = self.mesh().findCellByMarker(0, -1) Si = pg.ElementMatrix() St = pg.ElementMatrix() u = self.subPotentials pg.tic() if self.verbose(): print("Calculate sensitivity matrix for model: ", min(model), max(model)) Jt = pg.RMatrix(self.data.size(), self.regionManager().parameterCount()) for kIdx, w in enumerate(self.w): k = self.k[kIdx] w = self.w[kIdx] Jt *= 0. A = pg.ElementMatrixMap() for i, c in enumerate(cells): modelIdx = c.marker() # 2.5D Si.u2(c) Si *= k * k Si += St.ux2uy2uz2(c) # 3D # Si.ux2uy2uz2(c); w = w* 2 A.add(modelIdx, Si) for dataIdx in range(self.data.size()): a = int(self.data('a')[dataIdx]) b = int(self.data('b')[dataIdx]) m = int(self.data('m')[dataIdx]) n = int(self.data('n')[dataIdx]) Jt[dataIdx] = A.mult(u[kIdx][a] - u[kIdx][b], u[kIdx][m] - u[kIdx][n]) J += w * Jt m2 = model * model k = self.data('k') for i in range(J.rows()): J[i] /= (m2 / k[i]) if self.verbose(): sumsens = np.zeros(J.rows()) for i in range(J.rows()): sumsens[i] = pg.sum(J[i]) print("sens sum: median = ", pg.median(sumsens), " min = ", pg.min(sumsens), " max = ", pg.max(sumsens))
def assembleForceVector(mesh, f, userData=None): """Create right hand side vector based on the given mesh and force values. Create right hand side vector based on the given mesh and force values. Parameters ---------- f: float, array, callable(cell, [userData]), [...] Force Values float -> ones(mesh.nodeCount()) * vals, for each node [0 .. mesh.nodeCount()] for each cell [0 .. mesh.cellCount()] """ if isinstance(f, list) or hasattr(f, 'ndim'): if isinstance(f, list): rhs = np.zeros((len(f), mesh.nodeCount())) for i, fi in enumerate(f): userData['i'] = i rhs[i] = assembleForceVector(mesh, fi, userData) return rhs elif f.ndim == 2: # assume rhs [n, nNodes] array is already a valid return f rhs = pg.RVector(mesh.nodeCount(), 0) if hasattr(f, '__call__') and not isinstance(f, pg.RVector): for c in mesh.cells(): if userData is not None: f(c, rhs, userData) else: f(c, rhs) else: fArray = None if hasattr(f, '__len__'): if len(f) == mesh.cellCount() or len(f) == mesh.nodeCount(): fArray = f if fArray is None: fArray = parseArgToArray(f, mesh.cellCount(), mesh, userData) if len(fArray) == mesh.cellCount(): b_l = pg.ElementMatrix() for c in mesh.cells(): if fArray[c.id()] != 0.0: b_l.u(c) rhs.add(b_l, fArray[c.id()]) # print("test reference solution:") # rhsRef = pg.RVector(mesh.nodeCount(), 0) # for c in mesh.cells(): # b_l.u(c) # for i, idx in enumerate(b_l.idx()): # rhsRef[idx] += b_l.row(0)[i] * fArray[c.id()] # np.testing.assert_allclose(rhs, rhsRef) # print("Remove revtest in assembleForceVector after check") elif len(fArray) == mesh.nodeCount(): b_l = pg.ElementMatrix() for c in mesh.cells(): b_l.u(c) # rhs.addVal(b_l.row(0) * fArray[b_l.idx()], b_l.idx()) rhs.add(b_l, fArray) print("test reference solution:") rhsRef = pg.RVector(mesh.nodeCount(), 0) for c in mesh.cells(): b_l.u(c) for i, idx in enumerate(b_l.idx()): rhsRef[idx] += b_l.row(0)[i] * fArray[idx] np.testing.assert_allclose(rhs, rhsRef) print("Remove revtest in assembleForceVector after check") # rhs = pg.RVector(fArray) else: raise Exception("Forcevector have the wrong size: " + str(len(fArray))) return rhs