def compute_dist_face_cell( g, subcell_topology, eta, eta_at_bnd=False, return_paired=True ): """ Compute vectors from cell centers continuity points on each sub-face. The location of the continuity point is given by x_cp = (1-eta) * x_facecenter + eta * x_vertex On the boundary, eta is set to zero, thus the continuity point is at the face center Parameters ---------- g: Grid subcell_topology: Of class subcell topology in this module eta: [0,1), eta = 0 gives cont. pt. at face midpoint, eta = 1 means at the vertex Returns ------- sps.csr() matrix representation of vectors. Size g.nf x (g.nc * g.nd) """ _, blocksz = matrix_compression.rlencode( np.vstack((subcell_topology.cno, subcell_topology.nno)) ) dims = g.dim _, cols = np.meshgrid(subcell_topology.subhfno, np.arange(dims)) cols += matrix_compression.rldecode(np.cumsum(blocksz) - blocksz[0], blocksz) eta_vec = eta * np.ones(subcell_topology.fno.size) # Set eta values to zero at the boundary if not eta_at_bnd: bnd = np.in1d(subcell_topology.fno, g.get_all_boundary_faces()) eta_vec[bnd] = 0 cp = g.face_centers[:, subcell_topology.fno] + eta_vec * ( g.nodes[:, subcell_topology.nno] - g.face_centers[:, subcell_topology.fno] ) dist = cp - g.cell_centers[:, subcell_topology.cno] ind_ptr = np.hstack((np.arange(0, cols.size, dims), cols.size)) mat = sps.csr_matrix((dist.ravel("F"), cols.ravel("F"), ind_ptr)) if return_paired: return subcell_topology.pair_over_subfaces(mat) else: return mat
def _tensor_vector_prod(g, k, subcell_topology, apertures=None): """ Compute product of normal vectors and tensors on a sub-cell level. This is essentially defining Darcy's law for each sub-face in terms of sub-cell gradients. Thus, we also implicitly define the global ordering of sub-cell gradient variables (via the interpretation of the columns in nk). NOTE: In the local numbering below, in particular in the variables i and j, it is tacitly assumed that g.dim == g.nodes.shape[0] == g.face_normals.shape[0] etc. See implementation note in main method. Parameters: g (core.grids.grid): Discretization grid k (core.constit.second_order_tensor): The permeability tensor subcell_topology (fvutils.SubcellTopology): Wrapper class containing subcell numbering. Returns: nk: sub-face wise product of normal vector and permeability tensor. cell_node_blocks pairings of node and cell indices, which together define a sub-cell. sub_cell_ind: index of all subcells """ # Stack cell and nodes, and remove duplicate rows. Since subcell_mapping # defines cno and nno (and others) working cell-wise, this will # correspond to a unique rows (Matlab-style) from what I understand. # This also means that the pairs in cell_node_blocks uniquely defines # subcells, and can be used to index gradients etc. cell_node_blocks, blocksz = matrix_compression.rlencode( np.vstack((subcell_topology.cno, subcell_topology.nno))) nd = g.dim # Duplicates in [cno, nno] corresponds to different faces meeting at the # same node. There should be exactly nd of these. This test will fail # for pyramids in 3D assert np.all(blocksz == nd) # Define row and column indices to be used for normal_vectors * perm. # Rows are based on sub-face numbers. # Columns have nd elements for each sub-cell (to store a gradient) and # is adjusted according to block sizes _, j = np.meshgrid(subcell_topology.subhfno, np.arange(nd)) sum_blocksz = np.cumsum(blocksz) j += matrix_compression.rldecode(sum_blocksz - blocksz[0], blocksz) # Distribute faces equally on the sub-faces num_nodes = np.diff(g.face_nodes.indptr) normals = g.face_normals[:, subcell_topology.fno] / num_nodes[ subcell_topology.fno] if apertures is not None: normals = normals * apertures[subcell_topology.cno] # Represent normals and permeability on matrix form ind_ptr = np.hstack((np.arange(0, j.size, nd), j.size)) normals_mat = sps.csr_matrix((normals.ravel('F'), j.ravel('F'), ind_ptr)) k_mat = sps.csr_matrix( (k.perm[::, ::, cell_node_blocks[0]].ravel('F'), j.ravel('F'), ind_ptr)) nk = normals_mat * k_mat # Unique sub-cell indexes are pulled from column indices, we only need # every nd column (since nd faces of the cell meet at each vertex) sub_cell_ind = j[::, 0::nd] return nk, cell_node_blocks, sub_cell_ind