def _getNearestCellID(self, points): """ Test cases >>> from fipy import * >>> m0 = Grid2D(dx=(.1, 1., 10.), dy=(.1, 1., 10.)) >>> m1 = Grid2D(nx=2, ny=2, dx=5., dy=5.) >>> print m0._getNearestCellID(m1.getCellCenters().getGlobalValue()) [4 5 7 8] """ if self.globalNumberOfCells == 0: return numerix.arange(0) points = numerix.resize(points, (self.globalNumberOfCells, len(points), len(points[0]))).swapaxes(0,1) centers = self.getCellCenters().getGlobalValue()[...,numerix.newaxis] try: tmp = centers - points except TypeError: tmp = centers - PhysicalField(points) return numerix.argmin(numerix.dot(tmp, tmp, axis = 0), axis=0)
def _calcDistanceFunction(self, extensionVariable = None, narrowBandWidth = None, deleteIslands = False): if narrowBandWidth == None: narrowBandWidth = self.narrowBandWidth ## calculate interface values cellToCellIDs = self.mesh._getCellToCellIDs() if deleteIslands: adjVals = numerix.take(self.value, cellToCellIDs) adjInterfaceValues = MA.masked_array(adjVals, mask = (adjVals * self.value) > 0) masksum = numerix.sum(numerix.logical_not(MA.getmask(adjInterfaceValues)), 0) tmp = MA.logical_and(masksum == 4, self.value > 0) self.value = MA.where(tmp, -1, self.value) adjVals = numerix.take(self.value, cellToCellIDs) adjInterfaceValues = MA.masked_array(adjVals, mask = (adjVals * self.value) > 0) dAP = self.mesh._getCellToCellDistances() distances = abs(self.value * dAP / (self.value - adjInterfaceValues)) indices = MA.argsort(distances, 0) sign = (self.value > 0) * 2 - 1 s = distances[indices[0], numerix.arange(indices.shape[1])] if self.mesh.getDim() == 2: t = distances[indices[1], numerix.arange(indices.shape[1])] u = distances[indices[2], numerix.arange(indices.shape[1])] if indices.shape[1] > 0: ns = self.cellNormals[..., indices[0], numerix.arange(indices.shape[1])] nt = self.cellNormals[..., indices[1], numerix.arange(indices.shape[1])] else: ns = MA.zeros(self.cellNormals.shape[:-1] + (0,)) nt = MA.zeros(self.cellNormals.shape[:-1] + (0,)) signedDistance = MA.where(MA.getmask(s), self.value, MA.where(MA.getmask(t), sign * s, MA.where(abs(numerix.dot(ns,nt)) < 0.9, sign * s * t / MA.sqrt(s**2 + t**2), MA.where(MA.getmask(u), sign * s, sign * s * u / MA.sqrt(s**2 + u**2) ) ) ) ) else: signedDistance = MA.where(MA.getmask(s), self.value, sign * s) self.value = signedDistance ## calculate interface flag masksum = numerix.sum(numerix.logical_not(MA.getmask(distances)), 0) interfaceFlag = (masksum > 0).astype('l') ## spread the extensionVariable to the whole interface flag = True if extensionVariable is None: extensionVariable = numerix.zeros(self.mesh.getNumberOfCells(), 'd') flag = False ext = numerix.zeros(self.mesh.getNumberOfCells(), 'd') positiveInterfaceFlag = numerix.where(self.value > 0, interfaceFlag, 0) negativeInterfaceIDs = numerix.nonzero(numerix.where(self.value < 0, interfaceFlag, 0))[0] for id in negativeInterfaceIDs: tmp, extensionVariable[...,id] = self._calcTrialValue(id, positiveInterfaceFlag, extensionVariable) if flag: self.value = self.tmpValue.copy() ## evaluate the trialIDs adjInterfaceFlag = numerix.take(interfaceFlag, cellToCellIDs) hasAdjInterface = (numerix.sum(MA.filled(adjInterfaceFlag, 0), 0) > 0).astype('l') trialFlag = numerix.logical_and(numerix.logical_not(interfaceFlag), hasAdjInterface).astype('l') trialIDs = list(numerix.nonzero(trialFlag)[0]) evaluatedFlag = interfaceFlag for id in trialIDs: self.value[...,id], extensionVariable[id] = self._calcTrialValue(id, evaluatedFlag, extensionVariable) while len(trialIDs): id = trialIDs[numerix.argmin(abs(numerix.take(self.value, trialIDs)))] if abs(self.value[...,id]) > narrowBandWidth / 2: break trialIDs.remove(id) evaluatedFlag[...,id] = 1 for adjID in MA.filled(cellToCellIDs[...,id], -1): if adjID != -1: if not evaluatedFlag[...,adjID]: self.value[...,adjID], extensionVariable[...,adjID] = self._calcTrialValue(adjID, evaluatedFlag, extensionVariable) if adjID not in trialIDs: trialIDs.append(adjID) self.value = numerix.array(self.value)
def _getAddedMeshValues(self, other, resolution=1e-2): """Calculate the parameters to define a concatenation of `other` with `self` :Parameters: - `other`: The :class:`~fipy.meshes.numMesh.Mesh` to concatenate with `self` - `resolution`: How close vertices have to be (relative to the smallest cell-to-cell distance in either mesh) to be considered the same :Returns: A `dict` with 3 elements: the new mesh vertexCoords, faceVertexIDs, and cellFaceIDs. """ selfc = self._getConcatenableMesh() other = other._getConcatenableMesh() selfNumFaces = selfc.faceVertexIDs.shape[-1] selfNumVertices = selfc.vertexCoords.shape[-1] otherNumFaces = other.faceVertexIDs.shape[-1] otherNumVertices = other.vertexCoords.shape[-1] ## check dimensions if(selfc.vertexCoords.shape[0] != other.vertexCoords.shape[0]): raise MeshAdditionError, "Dimensions do not match" ## compute vertex correlates ## only try to match exterior (X) vertices self_Xvertices = numerix.unique(selfc._getFaceVertexIDs().filled()[..., selfc.getExteriorFaces().getValue()].flatten()) other_Xvertices = numerix.unique(other._getFaceVertexIDs().filled()[..., other.getExteriorFaces().getValue()].flatten()) self_XvertexCoords = selfc.vertexCoords[..., self_Xvertices] other_XvertexCoords = other.vertexCoords[..., other_Xvertices] # lifted from Mesh._getNearestCellID() other_vertexCoordMap = numerix.resize(other_XvertexCoords, (self_XvertexCoords.shape[-1], other_XvertexCoords.shape[0], other_XvertexCoords.shape[-1])).swapaxes(0,1) tmp = self_XvertexCoords[..., numerix.newaxis] - other_vertexCoordMap closest = numerix.argmin(numerix.dot(tmp, tmp), axis=0) # just because they're closest, doesn't mean they're close tmp = self_XvertexCoords[..., closest] - other_XvertexCoords distance = numerix.sqrtDot(tmp, tmp) # only want vertex pairs that are 100x closer than the smallest # cell-to-cell distance close = distance < resolution * min(selfc._getCellToCellDistances().min(), other._getCellToCellDistances().min()) vertexCorrelates = numerix.array((self_Xvertices[closest[close]], other_Xvertices[close])) # warn if meshes don't touch, but allow it if (selfc._getNumberOfVertices() > 0 and other._getNumberOfVertices() > 0 and vertexCorrelates.shape[-1] == 0): import warnings warnings.warn("Vertices are not aligned", UserWarning, stacklevel=4) ## compute face correlates # ensure that both sets of faceVertexIDs have the same maximum number of (masked) elements self_faceVertexIDs = selfc.faceVertexIDs other_faceVertexIDs = other.faceVertexIDs diff = self_faceVertexIDs.shape[0] - other_faceVertexIDs.shape[0] if diff > 0: other_faceVertexIDs = numerix.append(other_faceVertexIDs, -1 * numerix.ones((diff,) + other_faceVertexIDs.shape[1:]), axis=0) other_faceVertexIDs = MA.masked_values(other_faceVertexIDs, -1) elif diff < 0: self_faceVertexIDs = numerix.append(self_faceVertexIDs, -1 * numerix.ones((-diff,) + self_faceVertexIDs.shape[1:]), axis=0) self_faceVertexIDs = MA.masked_values(self_faceVertexIDs, -1) # want self's Faces for which all faceVertexIDs are in vertexCorrelates self_matchingFaces = numerix.in1d(self_faceVertexIDs, vertexCorrelates[0]).reshape(self_faceVertexIDs.shape).all(axis=0).nonzero()[0] # want other's Faces for which all faceVertexIDs are in vertexCorrelates other_matchingFaces = numerix.in1d(other_faceVertexIDs, vertexCorrelates[1]).reshape(other_faceVertexIDs.shape).all(axis=0).nonzero()[0] # map other's Vertex IDs to new Vertex IDs, # accounting for overlaps with self's Vertex IDs vertex_map = numerix.empty(otherNumVertices, dtype=int) verticesToAdd = numerix.delete(numerix.arange(otherNumVertices), vertexCorrelates[1]) vertex_map[verticesToAdd] = numerix.arange(otherNumVertices - len(vertexCorrelates[1])) + selfNumVertices vertex_map[vertexCorrelates[1]] = vertexCorrelates[0] # calculate hashes of faceVertexIDs for comparing Faces if self_matchingFaces.shape[-1] == 0: self_faceHash = numerix.empty(self_matchingFaces.shape[:-1] + (0,), dtype="str") else: # sort each of self's Face's vertexIDs for canonical comparison self_faceHash = numerix.sort(self_faceVertexIDs[..., self_matchingFaces], axis=0) # then hash the Faces for comparison (NumPy set operations are only for 1D arrays) self_faceHash = numerix.apply_along_axis(str, axis=0, arr=self_faceHash) face_sort = numerix.argsort(self_faceHash) self_faceHash = self_faceHash[face_sort] self_matchingFaces = self_matchingFaces[face_sort] if other_matchingFaces.shape[-1] == 0: other_faceHash = numerix.empty(other_matchingFaces.shape[:-1] + (0,), dtype="str") else: # convert each of other's Face's vertexIDs to new IDs other_faceHash = vertex_map[other_faceVertexIDs[..., other_matchingFaces]] # sort each of other's Face's vertexIDs for canonical comparison other_faceHash = numerix.sort(other_faceHash, axis=0) # then hash the Faces for comparison (NumPy set operations are only for 1D arrays) other_faceHash = numerix.apply_along_axis(str, axis=0, arr=other_faceHash) face_sort = numerix.argsort(other_faceHash) other_faceHash = other_faceHash[face_sort] other_matchingFaces = other_matchingFaces[face_sort] self_matchingFaces = self_matchingFaces[numerix.in1d(self_faceHash, other_faceHash)] other_matchingFaces = other_matchingFaces[numerix.in1d(other_faceHash, self_faceHash)] faceCorrelates = numerix.array((self_matchingFaces, other_matchingFaces)) # warn if meshes don't touch, but allow it if (selfc._getNumberOfFaces() > 0 and other._getNumberOfFaces() > 0 and faceCorrelates.shape[-1] == 0): import warnings warnings.warn("Faces are not aligned", UserWarning, stacklevel=4) # map other's Face IDs to new Face IDs, # accounting for overlaps with self's Face IDs face_map = numerix.empty(otherNumFaces, dtype=int) facesToAdd = numerix.delete(numerix.arange(otherNumFaces), faceCorrelates[1]) face_map[facesToAdd] = numerix.arange(otherNumFaces - len(faceCorrelates[1])) + selfNumFaces face_map[faceCorrelates[1]] = faceCorrelates[0] other_faceVertexIDs = vertex_map[other.faceVertexIDs[..., facesToAdd]] # ensure that both sets of cellFaceIDs have the same maximum number of (masked) elements self_cellFaceIDs = selfc.cellFaceIDs other_cellFaceIDs = face_map[other.cellFaceIDs] diff = self_cellFaceIDs.shape[0] - other_cellFaceIDs.shape[0] if diff > 0: other_cellFaceIDs = numerix.append(other_cellFaceIDs, -1 * numerix.ones((diff,) + other_cellFaceIDs.shape[1:]), axis=0) other_cellFaceIDs = MA.masked_values(other_cellFaceIDs, -1) elif diff < 0: self_cellFaceIDs = numerix.append(self_cellFaceIDs, -1 * numerix.ones((-diff,) + self_cellFaceIDs.shape[1:]), axis=0) self_cellFaceIDs = MA.masked_values(self_cellFaceIDs, -1) # concatenate everything and return return { 'vertexCoords': numerix.concatenate((selfc.vertexCoords, other.vertexCoords[..., verticesToAdd]), axis=1), 'faceVertexIDs': numerix.concatenate((self_faceVertexIDs, other_faceVertexIDs), axis=1), 'cellFaceIDs': MA.concatenate((self_cellFaceIDs, other_cellFaceIDs), axis=1) }