Пример #1
0
def slice_indices(A, slice_ind):
    """
    Function for slicing sparse matrix along rows or columns.
    If A is a csc_matrix A will be sliced along columns, while if A is a
    csr_matrix A will be sliced along the rows.

    Parameters
    ----------
    A (scipy.sparse.csc/csr_matrix): A sparse matrix.
    slice_ind (np.array): Array containing indices to be sliced

    Returns
    -------
    indices (np.array): If A is csc_matrix:
                            The nonzero row indices or columns slice_ind
                        If A is csr_matrix:
                            The nonzero columns indices or rows slice_ind
    Examples
    --------
    A = sps.csc_matrix(np.eye(10))
    rows = slice_indices(A, np.array([0,2,3]))
    """
    assert A.getformat() == "csc" or A.getformat() == "csr"
    if isinstance(slice_ind, int):
        indices = A.indices[slice(A.indptr[int(slice_ind)],
                                  A.indptr[int(slice_ind + 1)])]
    elif slice_ind.size == 1:
        indices = A.indices[slice(A.indptr[int(slice_ind)],
                                  A.indptr[int(slice_ind + 1)])]
    else:
        indices = A.indices[mcolon(A.indptr[slice_ind],
                                   A.indptr[slice_ind + 1])]
    return indices
Пример #2
0
def duplicate_faces(gh, face_cells):
    """
    Duplicate all faces that are connected to a lower-dim cell

    Parameters
    ----------
    gh         - Higher-dim grid
    face_cells - A list of connection matrices. Each matrix gives
                 the mapping from the cells of a lower-dim grid
                 to the faces of the higher diim grid.
    """
    # We find the indices of the higher-dim faces to be duplicated.
    # Each of these faces are duplicated, and the duplication is
    # attached to the same nodes. We do not attach the faces to
    # any cells as this connection will have to be undone later
    # anyway.
    frac_id = face_cells.nonzero()[1]
    frac_id = np.unique(frac_id)
    rem = gh.has_face_tag(FaceTag.BOUNDARY)[frac_id]
    gh.add_face_tag(frac_id[rem], FaceTag.FRACTURE)
    gh.remove_face_tag(frac_id, FaceTag.TIP)

    frac_id = frac_id[~rem]
    if frac_id.size == 0:
        return frac_id

    node_start = gh.face_nodes.indptr[frac_id]
    node_end = gh.face_nodes.indptr[frac_id + 1]
    nodes = gh.face_nodes.indices[mcolon(node_start, node_end)]
    added_node_pos = np.cumsum(node_end - node_start) + \
        gh.face_nodes.indptr[-1]
    assert (added_node_pos.size == frac_id.size)
    assert (added_node_pos[-1] - gh.face_nodes.indptr[-1] == nodes.size)
    gh.face_nodes.indices = np.hstack((gh.face_nodes.indices, nodes))
    gh.face_nodes.indptr = np.hstack((gh.face_nodes.indptr, added_node_pos))
    gh.face_nodes.data = np.hstack(
        (gh.face_nodes.data, np.ones(nodes.size, dtype=bool)))
    gh.face_nodes._shape = (gh.num_nodes,
                            gh.face_nodes.shape[1] + frac_id.size)
    assert (gh.face_nodes.indices.size == gh.face_nodes.indptr[-1])

    node_start = gh.face_nodes.indptr[frac_id]
    node_end = gh.face_nodes.indptr[frac_id + 1]

    #frac_nodes = gh.face_nodes[:, frac_id]

    #gh.face_nodes = sps.hstack((gh.face_nodes, frac_nodes))
    # We also copy the attributes of the original faces.
    gh.num_faces += frac_id.size
    gh.face_normals = np.hstack((gh.face_normals, gh.face_normals[:, frac_id]))
    gh.face_areas = np.append(gh.face_areas, gh.face_areas[frac_id])
    gh.face_centers = np.hstack((gh.face_centers, gh.face_centers[:, frac_id]))

    # Not sure if this still does the correct thing. Might have to
    # send in a logical array instead of frac_id.
    gh.add_face_tag(frac_id, FaceTag.FRACTURE | FaceTag.BOUNDARY)
    gh.remove_face_tag(frac_id, FaceTag.TIP)
    gh.face_tags = np.append(gh.face_tags, gh.face_tags[frac_id])

    return frac_id
Пример #3
0
def zero_out_sparse_rows(A, rows, diag=None):
    """
    zeros out given rows from sparse csr matrix. Optionally also set values on
    the diagonal.

    Parameters:
        A: Sparse matrix
        rows (np.ndarray of int): Indices of rows to be eliminated.
        diag (np.ndarray, double, optional): Values to be set to the diagonal
            on the eliminated rows.

    """

    # If matrix is not csr, it will be converted to csr, then the rows will be
    # zeroed, and the matrix converted back.
    flag = False
    if not A.getformat() == 'csr':
        mat_format = A.getformat()
        A = A.tocsr()
        flag = True

    ip = A.indptr
    row_indices = mcolon.mcolon(ip[rows], ip[rows + 1])
    A.data[row_indices] = 0
    if diag is not None:
        # now we set the diagonal
        diag_vals = np.zeros(A.shape[1])
        diag_vals[rows] = diag
        A += sps.dia_matrix((diag_vals, 0), shape=A.shape)

    if flag:
        # Convert matrix back
        A = A.astype(mat_format)

    return A
Пример #4
0
def grid_sequence_fixed_lines(
    basedim, num_levels, grid_type, pert=0, ref_rate=2, domain=None, subdom_func=None
):
    dim = basedim.shape[0]

    if domain is None:
        domain = np.ones(dim)

    for iter1 in range(num_levels):
        nx = basedim * ref_rate ** iter1
        g = make_grid(grid_type, nx, domain, dim)

        if pert > 0:
            g.compute_geometry()
            old_nodes = g.nodes.copy()
            dx = np.max(domain / nx)
            g = perturb(g, pert, dx)
            if subdom_func is not None:
                # Characteristic function for all cell centers
                xc = g.cell_centers
                chi = subdom_func(xc[0], xc[1])
                #
                chi_face = np.abs(g.cell_faces * chi)
                bnd_face = np.argwhere(chi_face > 0).squeeze(1)
                node_ptr = g.face_nodes.indptr
                node_ind = g.face_nodes.indices
                bnd_nodes = node_ind[
                    mcolon.mcolon(node_ptr[bnd_face], node_ptr[bnd_face + 1])
                ]
                g.nodes[:, bnd_nodes] = old_nodes[:, bnd_nodes]

        g.compute_geometry()
        yield g
Пример #5
0
    def test_mcolon_last_equal(self):
        # Motivated by Github issue #11
        indPtr = np.array([1, 5, 5, 6])
        select = np.array([0, 1])

        c = mcolon.mcolon(indPtr[select], indPtr[select + 1])
        c_known = np.array([1, 2, 3, 4])
        assert np.allclose(c, c_known)
Пример #6
0
    def test_mcolon_middle_equal(self):
        # Motivated by Github issue #11
        indPtr = np.array([1, 5, 5, 6])
        select = np.array([0, 1, 2])

        c = mcolon.mcolon(indPtr[select], indPtr[select + 1])
        c_known = np.array([1, 2, 3, 4, 5])
        self.assertTrue(np.allclose(c, c_known))
Пример #7
0
def _tag_faces(grids, check_highest_dim=True):
    """
    Tag faces of grids. Three different tags are given to different types of
    faces:
        NONE: None of the below (i.e. an internal face)
        DOMAIN_BOUNDARY: All faces that lie on the domain boundary
            (i.e. should be given a boundary condition).
        FRACTURE: All faces that are split (i.e. has a connection to a
            lower dim grid).
        TIP: A boundary face that is not on the domain boundary, nor
            coupled to a lower domentional domain.

    Parameters:
        grids (list): List of grids to be tagged. Sorted per dimension.
        check_highest_dim (boolean, default=True): If true, we require there is
            a single mesh in the highest dimension. The test is useful, but
            should be waived for dfn meshes.

    """

    # Assume only one grid of highest dimension
    if check_highest_dim:
        assert len(grids[0]) == 1, 'Must be exactly'\
            '1 grid of dim: ' + str(len(grids))

    for g_h in np.atleast_1d(grids[0]):
        bnd_faces = g_h.get_all_boundary_faces()
        domain_boundary_tags = np.zeros(
            g_h.num_faces, dtype=bool)
        domain_boundary_tags[bnd_faces] = True
        g_h.tags['domain_boundary_faces'] = domain_boundary_tags
        bnd_nodes, _, _ = sps.find(g_h.face_nodes[:, bnd_faces])
        bnd_nodes = np.unique(bnd_nodes)
        for g_dim in grids[1:-1]:
            for g in g_dim:
                # We find the global nodes of all boundary faces
                bnd_faces_l = g.get_all_boundary_faces()
                indptr = g.face_nodes.indptr
                fn_loc = mcolon.mcolon(
                    indptr[bnd_faces_l], indptr[bnd_faces_l + 1])
                nodes_loc = g.face_nodes.indices[fn_loc]
                # Convert to global numbering
                nodes_glb = g.global_point_ind[nodes_loc]
                # We then tag each node as a tip node if it is not a global
                # boundary node
                is_tip = np.in1d(nodes_glb, bnd_nodes, invert=True)
                # We reshape the nodes such that each column equals the nodes of
                # one face. If a face only contains global boundary nodes, the
                # local face is also a boundary face. Otherwise, we add a TIP tag.
                n_per_face = _nodes_per_face(g)
                is_tip = np.any(is_tip.reshape(
                    (n_per_face, bnd_faces_l.size), order='F'), axis=0)

                g.tags['tip_faces'][bnd_faces_l[is_tip]] = True
                domain_boundary_tags = np.zeros(g.num_faces, dtype=bool)
                domain_boundary_tags[bnd_faces_l[is_tip == False]] = True
                g.tags['domain_boundary_faces'] = domain_boundary_tags
Пример #8
0
    def get_boundary_nodes(self):
        """
        Get nodes on the boundary

        Returns:
            np.ndarray (1d), index of nodes on the boundary

        """
        b_faces = self.get_boundary_faces()
        first = self.face_nodes.indptr[b_faces]
        second = self.face_nodes.indptr[b_faces + 1]
        return np.unique(self.face_nodes.indices[mcolon.mcolon(first, second)])
Пример #9
0
def slice_mat(A, ind):
    """
    Function for slicing sparse matrix along rows or columns.
    If A is a csc_matrix A will be sliced along columns, while if A is a
    csr_matrix A will be sliced along the rows.

    Parameters
    ----------
    A (scipy.sparse.csc/csr_matrix): A sparse matrix.
    ind (np.array): Array containing indices to be sliced.

    Returns
    -------
    A_sliced (scipy.sparse.csc/csr_matrix): The sliced matrix
        if A is a csc_matrix A_sliced = A[:, ind]
        if A is a csr_matrix A_slice = A[ind, :]

    Examples
    --------
    A = sps.csc_matrix(np.eye(10))
    rows = slice_mat(A, np.array([0,2,3]))
    """
    assert A.getformat() == "csc" or A.getformat() == "csr"

    if np.asarray(ind).dtype == "bool":
        # convert to indices.
        # First check for dimension
        if ind.size != A.indptr.size - 1:
            raise IndexError("boolean index did not match indexed array")
        ind = np.where(ind)[0]

    if isinstance(ind, int):
        N = 1
        indptr = np.zeros(2)
        ind_slice = slice(A.indptr[int(ind)], A.indptr[int(ind + 1)])
    elif ind.size == 1:
        N = 1
        indptr = np.zeros(2)
        ind_slice = slice(A.indptr[int(ind)], A.indptr[int(ind + 1)])
    else:
        N = ind.size
        indptr = np.zeros(ind.size + 1)
        ind_slice = mcolon(A.indptr[ind], A.indptr[ind + 1])

    indices = A.indices[ind_slice]
    indptr[1:] = np.cumsum(A.indptr[ind + 1] - A.indptr[ind])
    data = A.data[ind_slice]

    if A.getformat() == "csc":
        return sps.csc_matrix((data, indices, indptr), shape=(A.shape[0], N))
    elif A.getformat() == "csr":
        return sps.csr_matrix((data, indices, indptr), shape=(N, A.shape[1]))
Пример #10
0
def slice_mat(A, ind):
    """
    Function for slicing sparse matrix along rows or columns.
    If A is a csc_matrix A will be sliced along columns, while if A is a
    csr_matrix A will be sliced along the rows.

    Parameters
    ----------
    A (scipy.sparse.csc/csr_matrix): A sparse matrix.
    ind (np.array): Array containing indices to be sliced.

    Returns
    -------
    A_sliced (scipy.sparse.csc/csr_matrix): The sliced matrix
        if A is a csc_matrix A_sliced = A[:, ind]
        if A is a csr_matrix A_slice = A[ind, :]

    Examples
    --------
    A = sps.csc_matrix(np.eye(10))
    rows = slice_mat(A, np.array([0,2,3]))
    """
    assert A.getformat() == 'csc' or A.getformat() == 'csr'

    if isinstance(ind, int):
        N = 1
        indptr = np.zeros(2)
        ind_slice = slice(A.indptr[int(ind)], A.indptr[int(ind + 1)])
    elif ind.size == 1:
        N = 1
        indptr = np.zeros(2)
        ind_slice = slice(A.indptr[int(ind)], A.indptr[int(ind + 1)])
    else:
        N = ind.size
        indptr = np.zeros(ind.size + 1)
        ind_slice = mcolon(A.indptr[ind], A.indptr[ind + 1])

    indices = A.indices[ind_slice]
    indptr[1:] = np.cumsum(A.indptr[ind + 1] - A.indptr[ind])
    data = A.data[ind_slice]

    if A.getformat() == 'csc':
        return sps.csc_matrix((data, indices, indptr), shape=(A.shape[0], N))
    elif A.getformat() == 'csr':
        return sps.csr_matrix((data, indices, indptr), shape=(N, A.shape[1]))
Пример #11
0
    def update_boundary_node_tag(self) -> None:
        """Tag nodes on the boundary of the grid with boundary tag."""

        mask = {
            "domain_boundary_faces": "domain_boundary_nodes",
            "fracture_faces": "fracture_nodes",
            "tip_faces": "tip_nodes",
        }
        zeros = np.zeros(self.num_nodes, dtype=bool)

        for face_tag, node_tag in mask.items():
            self.tags[node_tag] = zeros.copy()
            faces = np.where(self.tags[face_tag])[0]
            if faces.size > 0:
                first = self.face_nodes.indptr[faces]
                second = self.face_nodes.indptr[faces + 1]
                nodes = self.face_nodes.indices[mcolon.mcolon(first, second)]
                self.tags[node_tag][nodes] = True
Пример #12
0
def slice_indices(A, slice_ind, return_array_ind=False):
    """
    Function for slicing sparse matrix along rows or columns.
    If A is a csc_matrix A will be sliced along columns, while if A is a
    csr_matrix A will be sliced along the rows.

    Parameters
    ----------
    A (scipy.sparse.csc/csr_matrix): A sparse matrix.
    slice_ind (np.array): Array containing indices to be sliced

    Returns
    -------
    indices (np.array): If A is csc_matrix:
                            The nonzero row indices or columns slice_ind
                        If A is csr_matrix:
                            The nonzero columns indices or rows slice_ind
    Examples
    --------
    A = sps.csc_matrix(np.eye(10))
    rows = slice_indices(A, np.array([0,2,3]))
    """
    assert A.getformat() == "csc" or A.getformat() == "csr"
    if np.asarray(slice_ind).dtype == "bool":
        # convert to indices.
        # First check for dimension
        if slice_ind.size != A.indptr.size - 1:
            raise IndexError("boolean index did not match indexed array")
        slice_ind = np.where(slice_ind)[0]

    if isinstance(slice_ind, int):
        array_ind = slice(A.indptr[int(slice_ind)], A.indptr[int(slice_ind + 1)])
        indices = A.indices[array_ind]
    elif slice_ind.size == 1:
        array_ind = slice(A.indptr[int(slice_ind)], A.indptr[int(slice_ind + 1)])
        indices = A.indices[array_ind]
    else:
        array_ind = mcolon(A.indptr[slice_ind], A.indptr[slice_ind + 1])
        indices = A.indices[array_ind]
    if return_array_ind:
        return indices, array_ind
    else:
        return indices
Пример #13
0
def zero_columns(A, cols):
    """
    Function to zero out columns in matrix A. Note that this function does not
    change the sparcity structure of the matrix, it only changes the column
    values to 0

    Parameter
    ---------
    A (scipy.sparse.spmatrix): A sparce matrix
    cols (ndarray): A numpy array of columns that should be zeroed
    Return
    ------
    None


    """

    if A.getformat() != "csc":
        raise ValueError("Need a csc matrix")

    indptr = A.indptr
    col_indptr = mcolon(indptr[cols], indptr[cols + 1])
    A.data[col_indptr] = 0
Пример #14
0
def block_diag_index(m, n=None):
    """
    Get row and column indices for block diagonal matrix

    This is intended as the equivalent of the corresponding method in MRST.

    Examples:
    >>> m = np.array([2, 3])
    >>> n = np.array([1, 2])
    >>> i, j = block_diag_index(m, n)
    >>> i, j
    (array([0, 1, 2, 3, 4, 2, 3, 4]), array([0, 0, 1, 1, 1, 2, 2, 2]))
    >>> a = np.array([1, 3])
    >>> i, j = block_diag_index(a)
    >>> i, j
    (array([0, 1, 2, 3, 1, 2, 3, 1, 2, 3]), array([0, 1, 1, 1, 2, 2, 2, 3, 3, 3]))

    Parameters:
        m - ndarray, dimension 1
        n - ndarray, dimension 1, defaults to m

    """
    if n is None:
        n = m

    start = np.hstack((np.zeros(1, dtype='int'), m))
    pos = np.cumsum(start)
    p1 = pos[0:-1]
    p2 = pos[1:] - 1
    p1_full = matrix_compression.rldecode(p1, n)
    p2_full = matrix_compression.rldecode(p2, n)

    i = mcolon.mcolon(p1_full, p2_full + 1)
    sumn = np.arange(np.sum(n))
    m_n_full = matrix_compression.rldecode(m, n)
    j = matrix_compression.rldecode(sumn, m_n_full)
    return i, j
Пример #15
0
def _tag_faces(grids, check_highest_dim=True):
    """
        Tag faces of grids. Three different tags are given to different types of
        faces:
            NONE: None of the below (i.e. an internal face)
            DOMAIN_BOUNDARY: All faces that lie on the domain boundary
                (i.e. should be given a boundary condition).
            FRACTURE: All faces that are split (i.e. has a connection to a
                lower dim grid).
            TIP: A boundary face that is not on the domain boundary, nor
                coupled to a lower domentional domain.

        Parameters:
            grids (list): List of grids to be tagged. Sorted per dimension.
            check_highest_dim (boolean, default=True): If true, we require there is
                a single mesh in the highest dimension. The test is useful, but
                should be waived for dfn meshes.
            tag_tip_node (bool, default=True): If True, nodes in the highest-dimensional grid
    are tagged
    """

    # Assume only one grid of highest dimension
    if check_highest_dim:
        assert len(grids[0]) == 1, "Must be exactly" "1 grid of dim: " + str(len(grids))

    for g_h in np.atleast_1d(grids[0]):
        bnd_faces = g_h.get_all_boundary_faces()
        domain_boundary_tags = np.zeros(g_h.num_faces, dtype=bool)
        domain_boundary_tags[bnd_faces] = True
        g_h.tags["domain_boundary_faces"] = domain_boundary_tags
        bnd_nodes, _, _ = sps.find(g_h.face_nodes[:, bnd_faces])

        # Boundary nodes of g_h in terms of global indices
        bnd_nodes_glb = g_h.global_point_ind[np.unique(bnd_nodes)]

        # Find the nodes in g_h that are on the tip of a fracture (not counting
        # fracture endings on the global boundary). Do this by identifying tip nodes
        # among the grids of dimension dim_max - 1. Exclude tips that also occur on
        # other fractures (this will correspond to T or L-intersection, which look
        # like tips from individual fractures).

        # IMPLEMENTATION NOTE: To account for cases where not all global nodes are
        # present in g_h (e.g. for DFN-type grids), and avoid issues relating to nodes
        # in lower-dimensional grids that are not present in g_h, we store node numbers
        # instead of using booleans arrays on the nodes in g_h.

        # Keep track of nodes in g_h that correspond to tip nodes of a fracture.
        global_node_as_fracture_tip = np.array([], dtype=int)
        # Also count the number of occurences of nodes on fractures
        num_occ_nodes = np.array([], dtype=int)

        for g_dim in grids[1:-1]:
            for g in g_dim:
                # We find the global nodes of all boundary faces
                bnd_faces_l = g.get_all_boundary_faces()
                indptr = g.face_nodes.indptr
                fn_loc = mcolon.mcolon(indptr[bnd_faces_l], indptr[bnd_faces_l + 1])
                nodes_loc = g.face_nodes.indices[fn_loc]
                # Convert to global numbering
                nodes_glb = g.global_point_ind[nodes_loc]
                # We then tag each node as a tip node if it is not a global
                # boundary node
                node_on_tip_not_global_bnd = np.in1d(
                    nodes_glb, bnd_nodes_glb, invert=True
                )

                # We reshape the nodes such that each column equals the nodes of
                # one face. If a face only contains global boundary nodes, the
                # local face is also a boundary face. Otherwise, we add a TIP tag.
                # Note that we only consider boundary faces, hence any is okay.
                n_per_face = _nodes_per_face(g)
                is_tip_face = np.any(
                    node_on_tip_not_global_bnd.reshape(
                        (n_per_face, bnd_faces_l.size), order="F"
                    ),
                    axis=0,
                )

                # Tag faces on tips and boundaries
                g.tags["tip_faces"][bnd_faces_l[is_tip_face]] = True
                domain_boundary_tags = np.zeros(g.num_faces, dtype=bool)
                domain_boundary_tags[bnd_faces_l[np.logical_not(is_tip_face)]] = True
                g.tags["domain_boundary_faces"] = domain_boundary_tags

                # Also tag the nodes of the lower-dimensional grid that are on a tip
                is_tip_node = np.zeros(g.num_nodes, dtype=bool)
                is_tip_node[nodes_loc[node_on_tip_not_global_bnd]] = True
                g.tags["tip_nodes"] = is_tip_node

                if g.dim == g_h.dim - 1:
                    # For co-dimension 1, we also register those nodes in the host grid which
                    # are correspond to the tip of a fracture. We use a slightly wider
                    # definition of a fracture tip in this context: Nodes that are on the
                    # domain boundary, but also part of a tip face (on the fracture) which
                    # extends into the domain are also considered to be tip nodes. Filtering
                    # away these will be simple, using the domain_boundary_nodes tag, if
                    # necessary.
                    nodes_on_fracture_tip = np.unique(
                        nodes_glb.reshape((n_per_face, bnd_faces_l.size), order="F")[
                            :, is_tip_face
                        ]
                    )

                    global_node_as_fracture_tip = np.hstack(
                        (global_node_as_fracture_tip, nodes_on_fracture_tip)
                    )
                    # Count all global nodes used in this
                    num_occ_nodes = np.hstack((num_occ_nodes, g.global_point_ind))

        # The tip nodes should both be on the tip of a fracture, and not be present
        # on other fractures.
        may_be_tip = np.where(np.bincount(global_node_as_fracture_tip) == 1)[0]
        occurs_once = np.where(np.bincount(num_occ_nodes) == 1)[0]
        true_tip = np.intersect1d(may_be_tip, occurs_once)

        # Take the intersection between tip nodes and the nodes in this fracture.
        _, local_true_tip = pp.utils.setmembership.ismember_rows(
            true_tip, g_h.global_point_ind
        )
        tip_tag = np.zeros(g_h.num_nodes, dtype=bool)
        tip_tag[local_true_tip] = True
        # Tag nodes that are on the tip of a fracture, and not involved in other fractures
        g_h.tags["node_is_fracture_tip"] = tip_tag

        on_any_tip = np.where(np.bincount(global_node_as_fracture_tip) > 0)[0]
        _, local_any_tip = pp.utils.setmembership.ismember_rows(
            on_any_tip, g_h.global_point_ind
        )
        tip_of_a_fracture = np.zeros_like(tip_tag)
        tip_of_a_fracture[local_any_tip] = True
        # Tag nodes that are on the tip of a fracture, independent of whether it is
        g_h.tags["node_is_tip_of_some_fracture"] = tip_of_a_fracture
Пример #16
0
 def test_type_conversion(self):
     a = np.array([1, 3], dtype=np.int8)
     b = 1 + np.array([1, 3], dtype=np.int16)
     c = mcolon.mcolon(a, b)
     assert c.dtype == np.int64
Пример #17
0
 def test_mcolon_simple(self):
     a = np.array([1, 2])
     b = np.array([3, 4])
     c = mcolon.mcolon(a, b)
     assert np.all((c - np.array([1, 2, 2, 3])) == 0)
Пример #18
0
def merge_matrices(A, B, lines):
    """
    Replace rows/coloms of matrix A with rows/cols of matrix B.
    If A and B are csc matrices this function is equivalent with
    A[:, lines] = B
    If A and B are csr matrices this funciton is equivalent iwth
    A[lines, :] = B

    Parameter
    ---------
    A (scipy.sparse.spmatrix): A sparce matrix
    B (scipy.sparse.spmatrix): A sparce matrix
    lines (ndarray): Lines of A to be replaced by B.

    Return
    ------
    None


    """
    if A.getformat() != "csc" and A.getformat() != "csr":
        raise ValueError("Need a csc or csr matrix")
    elif A.getformat() != B.getformat():
        raise ValueError("A and B must be of same matrix type")
    if A.getformat() == "csc":
        if A.shape[0] != B.shape[0]:
            raise ValueError("A.shape[0] must equal B.shape[0]")
    if A.getformat() == "csr":
        if A.shape[1] != B.shape[1]:
            raise ValueError("A.shape[0] must equal B.shape[0]")

    if B.getformat() == "csc":
        if lines.size != B.shape[1]:
            raise ValueError("B.shape[1] must equal size of lines")
    if B.getformat() == "csr":
        if lines.size != B.shape[0]:
            raise ValueError("B.shape[0] must equal size of lines")

    if np.unique(lines).shape != lines.shape:
        raise ValueError("Can only merge unique lines")

    indptr = A.indptr
    indices = A.indices
    data = A.data

    ind_ix = mcolon(indptr[lines], indptr[lines + 1])

    # First we remove the old data
    num_rem = np.zeros(indptr.size, dtype=np.int32)
    num_rem[lines + 1] = indptr[lines + 1] - indptr[lines]
    num_rem = np.cumsum(num_rem, dtype=num_rem.dtype)

    indptr = indptr - num_rem

    keep = np.ones(A.data.size, dtype=np.bool)
    keep[ind_ix] = False
    indices = indices[keep]
    data = data[keep]

    # Then we add the new
    b_indptr = B.indptr
    b_indices = B.indices
    b_data = B.data

    num_added = np.zeros(indptr.size, dtype=np.int32)
    num_added[lines + 1] = b_indptr[1:] - b_indptr[:-1]
    num_added = np.cumsum(num_added, dtype=num_added.dtype)

    rep = np.diff(b_indptr)
    indPos = np.repeat(indptr[lines], rep)

    A.indices = np.insert(indices, indPos, b_indices)
    A.data = np.insert(data, indPos, b_data)
    A.indptr = indptr + num_added
Пример #19
0
 def test_mcolon_zero_output(self):
     a = np.array([1, 2])
     b = np.array([1, 2])
     c = mcolon.mcolon(a, b)
     assert c.size == 0
Пример #20
0
 def test_mcolon_one_missing(self):
     a = np.array([1, 2])
     b = np.array([3, 1])
     c = mcolon.mcolon(a, b)
     assert np.all((c - np.array([1, 2])) == 0)
Пример #21
0
def create_partition(A, seeds=None, **kwargs):
    """
    Create the partition based on an input matrix using the algebraic multigrid
    method coarse/fine-splittings based on direct couplings. The standard values
    for cdepth and epsilon are taken from the following reference.

    For more information see: U. Trottenberg, C. W. Oosterlee, and A. Schuller.
    Multigrid. Academic press, 2000.

    Parameters
    ----------
    A: sparse matrix used for the agglomeration
    cdepth: the greather is the more intense the aggregation will be, e.g. less
        cells if it is used combined with generate_coarse_grid
    epsilon: weight for the off-diagonal entries to define the "strong
        negatively cupling"
    seeds: (optional) to define a-priori coarse cells

    Returns
    -------
    out: agglomeration indices

    How to use
    ----------
    part = create_partition(tpfa_matrix(g))
    g = generate_coarse_grid(g, part)

    """

    cdepth = int(kwargs.get("cdepth", 2))
    epsilon = kwargs.get("epsilon", 0.25)

    if A.size == 0:
        return np.zeros(1)
    Nc = A.shape[0]

    # For each node, which other nodes are strongly connected to it
    ST = sps.lil_matrix((Nc, Nc), dtype=np.bool)

    # In the first instance, all cells are strongly connected to each other
    At = A.T

    for i in np.arange(Nc):
        loc = slice(At.indptr[i], At.indptr[i + 1])
        ci, vals = At.indices[loc], At.data[loc]
        neg = vals < 0.
        nvals = vals[neg]
        nci = ci[neg]
        minId = np.argmin(nvals)
        ind = -nvals >= epsilon * np.abs(nvals[minId])
        ST[nci[ind], i] = True

    # Temporary field, will store connections of depth 1
    for _ in np.arange(2, cdepth + 1):
        STold = ST.copy()
        for j in np.arange(Nc):
            rowj = np.array(STold.rows[j])
            if rowj.size == 0:
                continue
            row = np.hstack([STold.rows[r] for r in rowj])
            ST[j, np.concatenate((rowj, row))] = True

    del STold

    ST.setdiag(False)
    lmbda = np.array([len(s) for s in ST.rows])

    # Define coarse nodes
    candidate = np.ones(Nc, dtype=np.bool)
    is_fine = np.zeros(Nc, dtype=np.bool)
    is_coarse = np.zeros(Nc, dtype=np.bool)

    # cells that are not important for any other cells are on the fine scale.
    for row_id, row in enumerate(ST.rows):
        if not row:
            is_fine[row_id] = True
            candidate[row_id] = False

    ST = ST.tocsr()
    it = 0
    while np.any(candidate):
        i = np.argmax(lmbda)
        is_coarse[i] = True
        j = ST.indices[ST.indptr[i]:ST.indptr[i + 1]]
        jf = j[candidate[j]]
        is_fine[jf] = True
        candidate[np.r_[i, jf]] = False
        loop = ST.indices[mcolon.mcolon(ST.indptr[jf], ST.indptr[jf + 1])]
        for row in np.unique(loop):
            s = ST.indices[ST.indptr[row]:ST.indptr[row + 1]]
            lmbda[row] = s[candidate[s]].size + 2 * s[is_fine[s]].size
        lmbda[np.logical_not(candidate)] = -1
        it = it + 1

        # Something went wrong during aggregation
        assert it <= Nc

    del lmbda, ST

    if seeds is not None:
        is_coarse[seeds] = True
        is_fine[seeds] = False

    # If two neighbors are coarse, eliminate one of them without touching the
    # seeds
    c2c = np.abs(A) > 0
    c2c_rows, _, _ = sps.find(c2c)

    pairs = np.empty((0, 2), dtype=np.int)
    for idx, it in enumerate(np.where(is_coarse)[0]):
        loc = slice(c2c.indptr[it], c2c.indptr[it + 1])
        ind = np.setdiff1d(c2c_rows[loc], it)
        cind = ind[is_coarse[ind]]
        new_pair = np.stack((np.repeat(it, cind.size), cind), axis=-1)
        pairs = np.append(pairs, new_pair, axis=0)

    # Remove one of the neighbors cells
    if pairs.size:
        pairs = setmembership.unique_rows(np.sort(pairs, axis=1))[0]
        for ij in pairs:
            A_val = np.array(A[ij, ij]).ravel()
            ids = ij[np.argsort(A_val)]
            ids = np.setdiff1d(ids, seeds, assume_unique=True)
            if ids.size:
                is_coarse[ids[0]] = False
                is_fine[ids[0]] = True

    coarse = np.where(is_coarse)[0]

    # Primal grid
    NC = coarse.size
    primal = sps.lil_matrix((NC, Nc), dtype=np.bool)
    primal[np.arange(NC), coarse[np.arange(NC)]] = True

    connection = sps.lil_matrix((Nc, Nc), dtype=np.double)
    for it in np.arange(Nc):
        n = np.setdiff1d(c2c_rows[c2c.indptr[it]:c2c.indptr[it + 1]], it)
        loc = slice(A.indptr[it], A.indptr[it + 1])
        A_idx, A_row = A.indices[loc], A.data[loc]
        mask = A_idx != it
        connection[it, n] = np.abs(A_row[mask] / A_row[np.logical_not(mask)])

    connection = connection.tocsr()

    candidates_rep = np.ediff1d(connection.indptr)
    candidates_idx = np.repeat(is_coarse, candidates_rep)
    candidates = np.stack(
        (
            connection.indices[candidates_idx],
            np.repeat(np.arange(NC), candidates_rep[is_coarse]),
        ),
        axis=-1,
    )

    connection_idx = mcolon.mcolon(connection.indptr[coarse],
                                   connection.indptr[coarse + 1])
    vals = sps.csr_matrix(
        accumarray.accum(candidates,
                         connection.data[connection_idx],
                         size=[Nc, NC]))
    del candidates_rep, candidates_idx, connection_idx

    it = NC
    not_found = np.logical_not(is_coarse)
    # Process the strongest connection globally
    while np.any(not_found):

        np.argmax(vals.data)
        vals.argmax(axis=0)
        mcind = np.atleast_1d(np.squeeze(np.asarray(vals.argmax(axis=0))))
        mcval = -np.inf * np.ones(mcind.size)
        for c, r in enumerate(mcind):
            loc = slice(vals.indptr[r], vals.indptr[r + 1])
            vals_idx, vals_data = vals.indices[loc], vals.data[loc]
            mask = vals_idx == c
            if vals_idx.size == 0 or not np.any(mask):
                continue
            mcval[c] = vals_data[mask]

        mi = np.argmax(mcval)
        nadd = mcind[mi]

        primal[mi, nadd] = True
        it = it + 1
        if it > Nc + 5:
            break

        not_found[nadd] = False
        vals.data[vals.indptr[nadd]:vals.indptr[nadd + 1]] = 0

        loc = slice(connection.indptr[nadd], connection.indptr[nadd + 1])
        nc = connection.indices[loc]
        af = not_found[nc]
        nc = nc[af]
        nv = mcval[mi] * connection[nadd, :]
        nv = nv.data[af]
        if len(nc) > 0:
            vals += sps.csr_matrix((nv, (nc, np.repeat(mi, len(nc)))),
                                   shape=(Nc, NC))

    coarse, fine = primal.tocsr().nonzero()
    return coarse[np.argsort(fine)]
Пример #22
0
def _duplicate_specific_faces(gh: pp.Grid, frac_id: np.ndarray) -> np.ndarray:
    """
    Duplicate faces of gh specified by frac_id.
    """

    # Find which of the faces to split are tagged with a standard face tag,
    # that is, as fracture, tip or domain_boundary
    rem = tags.all_face_tags(gh.tags)[frac_id]

    # Set the faces to be split to fracture faces
    # Q: Why only if the face already had a tag (e.g., why [rem])?
    # Possible answer: We wil not split them (see redefinition of frac_id below),
    # but want them to be tagged as fracture_faces
    gh.tags["fracture_faces"][frac_id[rem]] = True
    # Faces to be split should not be tip
    gh.tags["tip_faces"][frac_id] = False

    # Only consider previously untagged faces for splitting
    frac_id = frac_id[~rem]
    if frac_id.size == 0:
        return frac_id

    # Expand the face-node relation to include duplicated nodes
    # Do this by directly manipulating the CSC-format of the matrix
    # Nodes of the target faces
    node_start = gh.face_nodes.indptr[frac_id]
    node_end = gh.face_nodes.indptr[frac_id + 1]
    nodes = gh.face_nodes.indices[mcolon(node_start, node_end)]

    # Start point for the new columns. They will be appended to the matrix, thus
    # the offset of the previous size of gh.face_nodes
    added_node_pos = np.cumsum(node_end -
                               node_start) + gh.face_nodes.indptr[-1]
    # Sanity checks
    assert added_node_pos.size == frac_id.size
    assert added_node_pos[-1] - gh.face_nodes.indptr[-1] == nodes.size
    # Expand row-data by adding node indices
    gh.face_nodes.indices = np.hstack((gh.face_nodes.indices, nodes))
    # Expand column pointers
    gh.face_nodes.indptr = np.hstack((gh.face_nodes.indptr, added_node_pos))
    # Expand data array
    gh.face_nodes.data = np.hstack(
        (gh.face_nodes.data, np.ones(nodes.size, dtype=bool)))
    # Update matrix shape
    gh.face_nodes._shape = (gh.num_nodes,
                            gh.face_nodes.shape[1] + frac_id.size)
    assert gh.face_nodes.indices.size == gh.face_nodes.indptr[-1]

    # We also copy the attributes of the original faces.
    gh.num_faces += frac_id.size
    gh.face_normals = np.hstack((gh.face_normals, gh.face_normals[:, frac_id]))
    gh.face_areas = np.append(gh.face_areas, gh.face_areas[frac_id])
    gh.face_centers = np.hstack((gh.face_centers, gh.face_centers[:, frac_id]))

    # Not sure if this still does the correct thing. Might have to
    # send in a logical array instead of frac_id.
    gh.tags["fracture_faces"][frac_id] = True
    gh.tags["tip_faces"][frac_id] = False
    update_fields = gh.tags.keys()
    update_values: List[List[np.ndarray]] = [[]] * len(update_fields)
    for i, key in enumerate(update_fields):
        # faces related tags are doubled and the value is inherit from the original
        if key.endswith("_faces"):
            update_values[i] = gh.tags[key][frac_id]
    tags.append_tags(gh.tags, update_fields, update_values)

    return frac_id
Пример #23
0
def duplicate_faces(gh, face_cells):
    """
    Duplicate all faces that are connected to a lower-dim cell

    Parameters
    ----------
    gh         - Higher-dim grid
    face_cells - A list of connection matrices. Each matrix gives
                 the mapping from the cells of a lower-dim grid
                 to the faces of the higher diim grid.
    """
    # We find the indices of the higher-dim faces to be duplicated.
    # Each of these faces are duplicated, and the duplication is
    # attached to the same nodes. We do not attach the faces to
    # any cells as this connection will have to be undone later
    # anyway.
    frac_id = face_cells.nonzero()[1]
    frac_id = np.unique(frac_id)
    rem = tags.all_face_tags(gh.tags)[frac_id]
    gh.tags["fracture_faces"][frac_id[rem]] = True
    gh.tags["tip_faces"][frac_id] = False

    frac_id = frac_id[~rem]
    if frac_id.size == 0:
        return frac_id

    node_start = gh.face_nodes.indptr[frac_id]
    node_end = gh.face_nodes.indptr[frac_id + 1]
    nodes = gh.face_nodes.indices[mcolon(node_start, node_end)]
    added_node_pos = np.cumsum(node_end -
                               node_start) + gh.face_nodes.indptr[-1]
    assert added_node_pos.size == frac_id.size
    assert added_node_pos[-1] - gh.face_nodes.indptr[-1] == nodes.size
    gh.face_nodes.indices = np.hstack((gh.face_nodes.indices, nodes))
    gh.face_nodes.indptr = np.hstack((gh.face_nodes.indptr, added_node_pos))
    gh.face_nodes.data = np.hstack(
        (gh.face_nodes.data, np.ones(nodes.size, dtype=bool)))
    gh.face_nodes._shape = (gh.num_nodes,
                            gh.face_nodes.shape[1] + frac_id.size)
    assert gh.face_nodes.indices.size == gh.face_nodes.indptr[-1]

    node_start = gh.face_nodes.indptr[frac_id]
    node_end = gh.face_nodes.indptr[frac_id + 1]

    # frac_nodes = gh.face_nodes[:, frac_id]

    # gh.face_nodes = sps.hstack((gh.face_nodes, frac_nodes))
    # We also copy the attributes of the original faces.
    gh.num_faces += frac_id.size
    gh.face_normals = np.hstack((gh.face_normals, gh.face_normals[:, frac_id]))
    gh.face_areas = np.append(gh.face_areas, gh.face_areas[frac_id])
    gh.face_centers = np.hstack((gh.face_centers, gh.face_centers[:, frac_id]))

    # Not sure if this still does the correct thing. Might have to
    # send in a logical array instead of frac_id.
    gh.tags["fracture_faces"][frac_id] = True
    gh.tags["tip_faces"][frac_id] = False
    update_fields = gh.tags.keys()
    update_values = [[]] * len(update_fields)
    for i, key in enumerate(update_fields):
        # faces related tags are doubled and the value is inherit from the original
        if key.endswith("_faces"):
            update_values[i] = gh.tags[key][frac_id]
    tags.append_tags(gh.tags, update_fields, update_values)

    return frac_id