Ejemplo n.º 1
0
def find_corner(mesh):
    '''
    For 2D geometry
      find line (boundary between two bdr_attribute) and
      corner of lines
    '''
    use_parallel = hasattr(mesh, "GroupNVertices")

    if use_parallel:
        from mpi4py import MPI
        myid = MPI.COMM_WORLD.rank
        nprc = MPI.COMM_WORLD.size
        comm = MPI.COMM_WORLD

        from mfem.common.mpi_debug import nicePrint, niceCall
        from petram.helper.mpi_recipes import allgather, allgather_vector, gather_vector
        from petram.mesh.mesh_utils import distribute_shared_entity
        if not hasattr(mesh, "shared_info"):
            mesh.shared_info = distribute_shared_entity(mesh)
    else:
        myid = 0
        nprc = 1

    ndim = mesh.Dimension()
    sdim = mesh.SpaceDimension()
    ne = mesh.GetNEdges()
    assert ndim == 2, "find_edge_corner is for 3D mesh"

    get_edges = mesh.GetElementEdges
    get_attr = mesh.GetAttribute
    iattr = mesh.GetAttributeArray()  # min of this array is 1
    nattr = 0 if iattr.size == 0 else np.max(iattr)
    nb = mesh.GetNE()
    nbe = mesh.GetNBE()
    if use_parallel:
        nbe = sum(allgather(nbe))
    if nbe == 0:
        return {}, {}, {}

    if use_parallel:
        offset = np.hstack([0, np.cumsum(allgather(mesh.GetNEdges()))])
        offsetf = np.hstack([0, np.cumsum(allgather(mesh.GetNFaces()))])
        offsetv = np.hstack([0, np.cumsum(allgather(mesh.GetNV()))])
        myoffset = offset[myid]
        myoffsetf = offsetf[myid]
        myoffsetv = offsetv[myid]
        nattr = max(allgather(nattr))
        ne = sum(allgather(mesh.GetNEdges()))
    else:
        myoffset = np.array(0, dtype=int)
        myoffsetf = np.array(0, dtype=int)
        myoffsetv = np.array(0, dtype=int)

    if mesh.GetNBE() == 0:  # some parallel node may have zero boundary
        battrs = []
        iedges = np.array([], dtype=int)
    else:
        battrs = mesh.GetBdrAttributeArray()
        iedges = np.hstack([
            mesh.GetBdrElementEdgeIndex(ibdr) for ibdr in range(mesh.GetNBE())
        ]).astype(int, copy=False)

    line2edge = GlobalNamedList()
    line2edge.setlists(battrs, iedges)

    if use_parallel:
        ld, md = mesh.shared_info
        iedges = iedges + myoffset

    if use_parallel:
        for key2 in ld:
            if key2[0] == myid: continue
            iii = np.in1d(iedges, ld[key2][1], invert=True)
            if len(iii) == 0: continue
            iedges = iedges[iii]
            battrs = battrs[iii]

    line2realedge = GlobalNamedList()
    line2realedge.setlists(battrs, iedges)

    line2realvert = GlobalNamedList()
    for key in line2realedge:
        data = np.hstack([
            mesh.GetEdgeVertices(i - myoffset) + myoffsetv
            for i in line2realedge[key]
        ])
        if use_parallel:
            for key2 in ld:
                if key2[0] == myid: continue
                for lv, mv in zip(ld[key2][0], md[key2][0]):
                    iii = np.where(data == lv)[0]
                    data[iii] = mv
        line2realvert[key] = data

    line2realvert.sharekeys().gather(nprc, distribute=True)

    corners = GlobalNamedList()
    for key in line2realvert:
        seen = defaultdict(int)
        for iiv in line2realvert[key]:
            seen[iiv] += 1
        corners[key] = [kk for kk in seen if seen[kk] == 1]

    sorted_key = corners.sharekeys().globalkeys

    corners.allgather()
    u_own = np.unique(
        np.hstack([corners[key] for key in corners]).astype(int, copy=False))

    if use_parallel:
        idx = np.logical_and(u_own >= offsetv[myid], u_own < offsetv[myid + 1])
        u_own = u_own[idx]

    if len(u_own) > 0:
        vtx = np.hstack([mesh.GetVertexArray(i - myoffsetv) for i in u_own])
    else:
        vtx = np.atleast_1d([])
    if use_parallel:
        vtx = gather_vector(vtx)
        u_own = gather_vector(u_own)

    # sort vertex
    if myid == 0:
        vtx = vtx.reshape(-1, sdim)
        tmp = sorted([(k, tuple(x)) for k, x in enumerate(vtx)],
                     key=lambda x: x[1])
        if len(tmp) > 0:
            vtx = np.vstack([x[1] for x in tmp])
            u_own = np.hstack([[u_own[x[0]] for x in tmp]]).astype(int)
            ivert = np.arange(len(vtx), dtype=int) + 1
        else:
            u_own = np.atleast_1d([]).astype(int)
            ivert = np.atleast_1d([]).astype(int)

    if use_parallel:
        #if myid != 0:
        #    u_own = None; vtx = None
        u_own = comm.bcast(u_own)
        ivert = np.arange(len(u_own), dtype=int) + 1
        for key in ld:
            if key[0] == myid: continue
            for lv, mv in zip(ld[key][0], md[key][0]):
                iii = np.where(u_own == mv)[0]
                u_own[iii] = lv
        idx = np.logical_and(u_own >= offsetv[myid], u_own < offsetv[myid + 1])
        u_own = u_own[idx]
        vtx = comm.bcast(vtx)
        vtx = comm.bcast(vtx)[idx.flatten()]
        ivert = ivert[idx]

    vert2vert = {iv: iu - myoffsetv for iv, iu in zip(ivert, u_own)}

    # mapping line index to vertex index (not MFFEM vertex id)
    line2vert = {}
    #nicePrint(corners)
    corners.bcast(nprc, distributed=True)
    for j, key in enumerate(sorted_key):
        data = corners[key]
        if use_parallel:
            for key2 in ld:
                if key2[0] == myid: continue
                for lv, mv in zip(ld[key2][0], md[key2][0]):
                    iii = np.where(data == mv)[0]
                    data[iii] = lv
            idx = np.logical_and(data >= offsetv[myid],
                                 data < offsetv[myid + 1])
            data = data[idx]
        data = list(data - myoffsetv)
        line2vert[j + 1] = [k for k in vert2vert if vert2vert[k] in data]

    if debug:
        g = GlobalNamedList(line2vert)
        g.sharekeys()
        gg = g.gather(nprc, overwrite=False).unique()
        if myid == 0: print(gg)
        for i in range(nprc):
            if use_parallel: comm.barrier()
            if myid == i:
                for k in vert2vert:
                    print(myid, k, mesh.GetVertexArray(vert2vert[k]))
            if use_parallel: comm.barrier()
    return line2vert, line2edge, vert2vert
Ejemplo n.º 2
0
def find_edge_corner(mesh):
    '''
    For 3D geometry
      find line (boundary between two bdr_attribute) and
      corner of lines
    '''
    use_parallel = hasattr(mesh, "GroupNVertices")

    if use_parallel:
        from mpi4py import MPI
        myid = MPI.COMM_WORLD.rank
        nprc = MPI.COMM_WORLD.size
        comm = MPI.COMM_WORLD

        from mfem.common.mpi_debug import nicePrint, niceCall
        from petram.helper.mpi_recipes import allgather, allgather_vector, gather_vector
        from petram.mesh.mesh_utils import distribute_shared_entity
        if not hasattr(mesh, "shared_info"):
            mesh.shared_info = distribute_shared_entity(mesh)
    else:
        myid = 0
        nprc = 1
    ndim = mesh.Dimension()
    sdim = mesh.SpaceDimension()
    ne = mesh.GetNEdges()
    assert ndim == 3, "find_edge_corner is for 3D mesh"

    # 3D mesh
    get_edges = mesh.GetBdrElementEdges
    get_attr = mesh.GetBdrAttribute
    iattr = mesh.GetBdrAttributeArray()  # min of this array is 1
    nattr = 0 if iattr.size == 0 else np.max(iattr)
    nb = mesh.GetNBE()
    if mesh.GetNBE() == 0 and nprc == 1:
        return {}, {}, {}, {}

    if use_parallel:
        offset = np.hstack([0, np.cumsum(allgather(mesh.GetNEdges()))])
        offsetf = np.hstack([0, np.cumsum(allgather(mesh.GetNFaces()))])
        offsetv = np.hstack([0, np.cumsum(allgather(mesh.GetNV()))])
        myoffset = offset[myid]
        myoffsetf = offsetf[myid]
        myoffsetv = offsetv[myid]
        nattr = max(allgather(nattr))
        ne = sum(allgather(mesh.GetNEdges()))
    else:
        myoffset = np.array(0, dtype=int)
        myoffsetf = np.array(0, dtype=int)
        myoffsetv = np.array(0, dtype=int)

    edges = defaultdict(list)
    iedges = np.arange(nb, dtype=int)

    if use_parallel:
        # eliminate slave faces from consideration
        iface = np.array([mesh.GetBdrElementEdgeIndex(i) for i in iedges],
                         dtype=int) + myoffsetf
        mask = np.array([True] * len(iface), dtype=bool)
        ld, md = mesh.shared_info
        for key in ld.keys():
            mid, g_in_master = key
            if mid == myid: continue
            iii = np.in1d(iedges, ld[key][2], invert=True)
            mask = np.logical_and(mask, iii)
        iedges = iedges[mask]

    # nicePrint(len(iedges)) np 1,2,4 gives 900... ok

    for i in iedges:
        ie, io = get_edges(i)
        ie += myoffset
        iattr = get_attr(i)
        edges[iattr].extend(list(ie))

    if use_parallel:
        # collect edges using master edge number
        # and gather it to a node.
        edgesc = {}
        ld, md = mesh.shared_info
        for j in range(1, nattr + 1):
            if j in edges:
                data = np.array(edges[j], dtype=int)
                for key in ld.keys():
                    mid, g_in_master = key
                    if mid == myid: continue
                    for le, me in zip(ld[key][1], md[key][1]):
                        iii = np.where(data == le)[0]
                        data[iii] = me
            else:
                data = np.atleast_1d([]).astype(int)
            data = gather_vector(data, root=j % nprc)
            if data is not None: edgesc[j] = data
        edges = edgesc

    # for each iattr real edge appears only once
    for key in edges.keys():
        seen = defaultdict(int)
        for x in edges[key]:
            seen[x] += 1
        edges[key] = [k for k in seen if seen[k] == 1]

    #nicePrint('Num edges',
    nedge = sum([len(edges[k]) for k in edges])
    if nedge != 0:
        N = np.hstack(
            [np.zeros(len(edges[k]), dtype=int) + k - 1 for k in edges.keys()])
        M = np.hstack([np.array(edges[k]) for k in edges.keys()])
    else:
        N = np.atleast_1d([]).astype(int)
        M = np.atleast_1d([]).astype(int)
    M = M.astype(int, copy=False)
    N = N.astype(int, copy=False)

    if use_parallel:
        # send attribute to owner of edges
        for j in range(nprc):
            idx = np.logical_and(M >= offset[j], M < offset[j + 1])
            Mpart = M[idx]
            Npart = N[idx]
            Mpart = gather_vector(Mpart, root=j)
            Npart = gather_vector(Npart, root=j)
            if j == myid: M2, N2 = Mpart, Npart
        M, N = M2, N2

    #nicePrint('unique edge', len(np.unique(M)))
    #nicePrint('N', len(N))
    data = M * 0 + 1
    table1 = coo_matrix((data, (M, N)), shape=(ne, nattr), dtype=int)

    csr = table1.tocsr()

    #embeded surface only touches to one iattr
    idx = np.where(np.diff(csr.indptr) >= 1)[0]
    csr = csr[idx, :]

    # this is true bdr edges.
    bb_edges = defaultdict(list)
    indptr = csr.indptr
    indices = csr.indices

    for i in range(csr.shape[0]):
        idxs = tuple(sorted(indices[indptr[i]:indptr[i + 1]] + 1))
        bb_edges[idxs].append(idx[i])
    bb_edges.default_factory = None

    # sort keys (= attribute set)
    keys = list(bb_edges)
    if use_parallel:
        keys = comm.gather(keys)
        if myid == 0: keys = sum(keys, [])
    sorted_key = None
    if myid == 0:
        sorted_key = list(set(keys))
        sorted_key.sort(key=lambda x: (len(x), x))

    if use_parallel: sorted_key = comm.bcast(sorted_key, root=0)
    bb_edgess = OrderedDict()
    for k in sorted_key:
        if k in bb_edges:
            bb_edgess[k] = bb_edges[k]
        else:
            bb_edgess[k] = [
            ]  # in parallel, put empty so that key order is kept
    bb_edges = bb_edgess
    '''
    res = []
    for key in sorted_key:
        tmp = allgather(len(bb_edges[key]))
        if myid == 0:
            res.append((key, sum(tmp)))
    if myid == 0: print res
    '''
    # at this point each node has its own edges populated in bb_edges (no shadow)
    ivert = {}

    for k in sorted_key:
        if len(bb_edges[k]) > 0:
            ivert[k] = np.hstack([
                mesh.GetEdgeVertices(i - myoffset) + myoffsetv
                for i in np.unique(bb_edges[k])
            ]).astype(int)
        else:
            ivert[k] = np.atleast_1d([]).astype(int)

    if use_parallel:
        # convert shadow vertex to real
        for k in sorted_key:
            data = ivert[k]
            for key in ld:
                if key[0] == myid: continue
                for le, me in zip(ld[key][0], md[key][0]):
                    iii = np.where(data == le)[0]
                    data[iii] = me
            ivert[k] = data
        ivertc = {}
        for j, k in enumerate(sorted_key):
            data = gather_vector(ivert[k], root=j % nprc)
            if data is not None:
                ivertc[k] = data
        ivert = ivertc

    corners = {}
    for key in ivert:
        seen = defaultdict(int)
        for iiv in ivert[key]:
            seen[iiv] += 1
        corners[key] = [kk for kk in seen if seen[kk] == 1]

    if len(corners) == 0:
        u = np.atleast_1d([]).astype(int)
    else:
        u = np.unique(np.hstack([corners[key]
                                 for key in corners])).astype(int, copy=False)

    # collect vertex on each node and gather to node 0
    u_own = u
    if use_parallel:
        u = np.unique(allgather_vector(u))
        u_own = u.copy()
        for key in ld:
            if key[0] == myid: continue
            for lv, mv in zip(ld[key][0], md[key][0]):
                iii = np.where(u == mv)[0]
                u[iii] = lv
        idx = np.logical_and(u >= offsetv[myid], u < offsetv[myid + 1])
        u = u[idx]  # u include shared vertex
        idx = np.logical_and(u_own >= offsetv[myid], u_own < offsetv[myid + 1])
        u_own = u_own[idx]  # u_own is only owned vertex

    #nicePrint('u_own',mesh.GetNV(),",",  u_own)
    if len(u_own) > 0:
        vtx = np.vstack([mesh.GetVertexArray(i - myoffsetv) for i in u_own])
    else:
        vtx = np.atleast_1d([]).reshape(-1, sdim)
    if use_parallel:
        u_own = gather_vector(u_own)
        vtx = gather_vector(vtx.flatten())

    # sort vertex
    if myid == 0:
        vtx = vtx.reshape(-1, sdim)
        #print('vtx shape', vtx.shape)
        tmp = sorted([(k, tuple(x)) for k, x in enumerate(vtx)],
                     key=lambda x: x[1])
        if len(tmp) > 0:
            vtx = np.vstack([x[1] for x in tmp])
            u_own = np.hstack([[u_own[x[0]] for x in tmp]]).astype(int)
            ivert = np.arange(len(vtx), dtype=int) + 1
        else:
            vtx = np.atleast_1d([]).astype(float)
            u_own = np.atleast_1d([]).astype(int)
            u_own = np.atleast_1d([]).astype(int)
    if use_parallel:
        #if myid != 0:
        #    u_own = None; vtx = None
        u_own = comm.bcast(u_own)
        ivert = np.arange(len(u_own), dtype=int) + 1
        for key in ld:
            if key[0] == myid: continue
            for lv, mv in zip(ld[key][0], md[key][0]):
                iii = np.where(u_own == mv)[0]
                u_own[iii] = lv
        idx = np.logical_and(u_own >= offsetv[myid], u_own < offsetv[myid + 1])
        u_own = u_own[idx]
        ivert = ivert[idx]
        #vtx = comm.bcast(vtx)
        #vtx = comm.bcast(vtx)[idx.flatten()]

    vert2vert = {iv: iu - myoffsetv for iv, iu in zip(ivert, u_own)}
    #nicePrint('vert2vert', vert2vert)

    # mapping line index to vertex index (not MFFEM vertex id)
    line2vert = {}
    #nicePrint(corners)
    for j, key in enumerate(sorted_key):
        data = corners[key] if key in corners else None
        if use_parallel:
            data = comm.bcast(data, root=j % nprc)
            data = np.array(data, dtype=int)
            for key2 in ld:
                if key2[0] == myid: continue
                for lv, mv in zip(ld[key2][0], md[key2][0]):
                    iii = np.where(data == mv)[0]
                    data[iii] = lv
            idx = np.logical_and(data >= offsetv[myid],
                                 data < offsetv[myid + 1])
            data = data[idx]
        else:
            data = np.array(data, dtype=int)
        data = list(data - myoffsetv)
        line2vert[j + 1] = [k for k in vert2vert if vert2vert[k] in data]

    # finish-up edge data
    if use_parallel:
        # distribute edges, convert (add) from master to local
        # number
        for attr_set in bb_edges:
            data = sum(allgather(bb_edges[attr_set]), [])
            data = np.array(data, dtype=int)
            for key in ld:
                if key[0] == myid: continue
                for le, me in zip(ld[key][1], md[key][1]):
                    iii = np.where(data == me)[0]
                    data[iii] = le

            idx = np.logical_and(data >= offset[myid], data < offset[myid + 1])
            data = data[idx]
            bb_edges[attr_set] = list(data - myoffset)

        attrs = list(edges)
        attrsa = np.unique(sum(allgather(attrs), []))

        for a in attrsa:
            if a in attrs:
                data = np.array(edges[a], dtype=int)
            else:
                data = np.atleast_1d([]).astype(int)
            data = allgather_vector(data)

            for key in ld:
                if key[0] == myid: continue
                for le, me in zip(ld[key][1], md[key][1]):
                    iii = np.where(data == me)[0]
                    data[iii] = le
            idx = np.logical_and(data >= offset[myid], data < offset[myid + 1])
            data = data[idx]
            edges[a] = list(data - myoffset)

    line2edge = {}
    for k, attr_set in enumerate(sorted_key):
        if attr_set in bb_edges:
            line2edge[k + 1] = bb_edges[attr_set]
        else:
            line2edge[k + 1] = []
    '''
    # debug find true (non-shadow) edges
    line2edge_true = {}
    for k, attr_set in enumerate(sorted_key):
        if attr_set in bb_edges:
            data = np.array(bb_edges[attr_set], dtype=int)
            for key in ld:
                if key[0] == myid: continue                                
                iii = np.in1d(data+myoffset, ld[key][1], invert = True)
                data = data[iii]
            line2edge_true[k+1] = data
        else:
            line2edge_true[k+1] = []
    nicePrint([sum(allgather(len(line2edge_true[key]))) for key in line2edge])
    '''
    surf2line = {k + 1: [] for k in range(nattr)}
    for k, attr_set in enumerate(sorted_key):
        for a in attr_set:
            surf2line[a].append(k + 1)

    if debug:
        g = GlobalNamedList(line2vert)
        g.sharekeys()
        gg = g.gather(nprc, overwrite=False).unique()
        if myid == 0: print("debug (gathered line2vert)", gg)

    return surf2line, line2vert, line2edge, vert2vert
Ejemplo n.º 3
0
def volume(mesh, in_attr, filename='', precision=8):
    '''
    make a new mesh which contains only spedified attributes.

    note: 
       1) boundary elements are also copied and bdr_attributes
          are maintained
       2) in parallel, new mesh must be geometrically continuous.
          this routine does not check it
         
    mesh must have sdim == 3:
    in_attr : domain attribute
    filename : an option to save the file 

    return new volume mesh
    '''
    in_attr = np.atleast_1d(in_attr)
    sdim = mesh.SpaceDimension()
    dim = mesh.Dimension()
    Nodal = mesh.GetNodalFESpace()
    hasNodal = (Nodal is not None)

    if sdim != 3: assert False, "sdim must be three for volume mesh"
    if dim != 3: assert False, "sdim must be three for volume mesh"

    idx, attrs, ivert, nverts, base = _collect_data(in_attr, mesh, 'dom')

    v2s = mesh.extended_connectivity['vol2surf']
    in_battr = np.unique(np.hstack([v2s[k]
                                    for k in in_attr])).astype(int, copy=False)
    if isParMesh(mesh):
        in_battr = np.unique(allgather_vector(in_battr))

    bidx, battrs, bivert, nbverts, bbase = _collect_data(in_battr, mesh, 'bdr')
    iface = np.array([mesh.GetBdrElementEdgeIndex(i) for i in bidx], dtype=int)

    # note u is sorted unique
    u, indices = np.unique(np.hstack((ivert, bivert)), return_inverse=True)

    kbelem = np.array([True] * len(bidx), dtype=bool)
    u_own = u

    if isParMesh(mesh):
        shared_info = distribute_shared_entity(mesh)
        u_own, ivert, bivert = _gather_shared_vertex(mesh, u, shared_info,
                                                     ivert, bivert)

    if len(u_own) > 0:
        vtx = np.vstack([mesh.GetVertexArray(i) for i in u_own])
    else:
        vtx = np.array([]).reshape((-1, sdim))

    if isParMesh(mesh):
        #
        # distribute vertex/element data
        #
        base = allgather_vector(base)
        nverts = allgather_vector(nverts)
        attrs = allgather_vector(attrs)

        ivert = allgather_vector(ivert)
        bivert = allgather_vector(bivert)

        vtx = allgather_vector(vtx.flatten()).reshape(-1, sdim)

        u, indices = np.unique(np.hstack([ivert, bivert]), return_inverse=True)

        #
        # take care of shared boundary (face)
        #
        #  2018.11.28
        #   skip_adding is on. This basically skip shared_element
        #   processing. Check em3d_TEwg7 if you need to remov this.
        #
        kbelem, battrs, nbverts, bbase, bivert = (_gather_shared_element(
            mesh,
            'face',
            shared_info,
            iface,
            kbelem,
            battrs,
            nbverts,
            bbase,
            bivert,
            skip_adding=True))

    #indices0  = np.array([np.where(u == biv)[0][0] for biv in ivert])
    #bindices0 = np.array([np.where(u == biv)[0][0] for biv in bivert])

    iv, ivi = np.unique(ivert, return_inverse=True)
    tmp = np.where(np.in1d(u, ivert, assume_unique=True))[0]
    indices = tmp[ivi]
    iv, ivi = np.unique(bivert, return_inverse=True)
    tmp = np.where(np.in1d(u, bivert, assume_unique=True))[0]
    bindices = tmp[ivi]

    #print('check', np.sum(np.abs(indices - indices0)))

    Nvert = len(vtx)
    Nelem = len(attrs)
    Nbelem = np.sum(kbelem)  #len(battrs)

    if myid == 0:
        print("NV, NBE, NE: " +
              ",".join([str(x) for x in (Nvert, Nbelem, Nelem)]))

    omesh = mfem.Mesh(3, Nvert, Nelem, Nbelem, sdim)
    #omesh = mfem.Mesh(3, Nvert, Nelem, 0, sdim)

    _fill_mesh_elements(omesh, vtx, indices, nverts, attrs, base)
    _fill_mesh_bdr_elements(omesh, vtx, bindices, nbverts, battrs, bbase,
                            kbelem)

    omesh.FinalizeTopology()
    omesh.Finalize(refine=True, fix_orientation=True)

    if hasNodal:
        odim = omesh.Dimension()
        fec = Nodal.FEColl()
        dNodal = mfem.FiniteElementSpace(omesh, fec, sdim)
        omesh.SetNodalFESpace(dNodal)
        omesh._nodal = dNodal

        GetXDofs = Nodal.GetElementDofs
        GetNX = Nodal.GetNE
        dGetXDofs = dNodal.GetElementDofs
        dGetNX = dNodal.GetNE

        DofToVDof = Nodal.DofToVDof
        dDofToVDof = dNodal.DofToVDof

        #nicePrint(dGetNX(),',', GetNX())
        nodes = mesh.GetNodes()
        node_ptx1 = nodes.GetDataArray()

        onodes = omesh.GetNodes()
        node_ptx2 = onodes.GetDataArray()
        #nicePrint(len(idx), idx)

        if len(idx) > 0:
            dof1_idx = np.hstack([[DofToVDof(i, d) for d in range(sdim)]
                                  for j in idx for i in GetXDofs(j)])
            data = node_ptx1[dof1_idx]
        else:
            dof1_idx = np.array([])
            data = np.array([])
        if isParMesh(mesh): data = allgather_vector(data)
        if isParMesh(mesh): idx = allgather_vector(idx)
        #nicePrint(len(data), ',', len(idx))

        dof2_idx = np.hstack([[dDofToVDof(i, d) for d in range(sdim)]
                              for j in range(len(idx)) for i in dGetXDofs(j)])
        node_ptx2[dof2_idx] = data
        #nicePrint(len(dof2_idx))

    if isParMesh(mesh):
        omesh = mfem.ParMesh(comm, omesh)

    if filename != '':
        if isParMesh(mesh):
            smyid = '{:0>6d}'.format(myid)
            filename = filename + '.' + smyid
        omesh.PrintToFile(filename, precision)

    return omesh
Ejemplo n.º 4
0
def surface(mesh, in_attr, filename='', precision=8):
    '''
    mesh must be 
    if sdim == 3:
       a domain of   2D mesh
       a boundary of 3D mesh
    if sdim == 2:
       a domain  in 2D mesh

    in_attr : eihter
    filename : an option to save the file 
    return new surface mesh

    '''
    sdim = mesh.SpaceDimension()
    dim = mesh.Dimension()
    Nodal = mesh.GetNodalFESpace()
    hasNodal = (Nodal is not None)

    if sdim == 3 and dim == 3:
        mode = 'bdr', 'edge'
    elif sdim == 3 and dim == 2:
        mode = 'dom', 'bdr'
    elif sdim == 2 and dim == 2:
        mode = 'dom', 'bdr'
    else:
        assert False, "unsupported mdoe"

    idx, attrs, ivert, nverts, base = _collect_data(in_attr, mesh, mode[0])

    s2l = mesh.extended_connectivity['surf2line']
    in_eattr = np.unique(np.hstack([s2l[k]
                                    for k in in_attr])).astype(int, copy=False)
    if isParMesh(mesh):
        in_eattr = np.unique(allgather_vector(in_eattr))

    eidx, eattrs, eivert, neverts, ebase = _collect_data(
        in_eattr, mesh, mode[1])

    u, indices = np.unique(np.hstack((ivert, eivert)), return_inverse=True)
    keelem = np.array([True] * len(eidx), dtype=bool)
    u_own = u

    if isParMesh(mesh):
        shared_info = distribute_shared_entity(mesh)
        u_own, ivert, eivert = _gather_shared_vertex(mesh, u, shared_info,
                                                     ivert, eivert)
    Nvert = len(u)
    if len(u_own) > 0:
        vtx = np.vstack([mesh.GetVertexArray(i) for i in u_own])
    else:
        vtx = np.array([]).reshape((-1, sdim))

    if isParMesh(mesh):
        #
        # distribute vertex/element data
        #
        base = allgather_vector(base)
        nverts = allgather_vector(nverts)
        attrs = allgather_vector(attrs)

        ivert = allgather_vector(ivert)
        eivert = allgather_vector(eivert)

        vtx = allgather_vector(vtx.flatten()).reshape(-1, sdim)

        u, indices = np.unique(np.hstack([ivert, eivert]), return_inverse=True)

        #
        # take care of shared boundary (edge)
        #
        keelem, eattrs, neverts, ebase, eivert = (_gather_shared_element(
            mesh,
            'edge',
            shared_info,
            eidx,
            keelem,
            eattrs,
            neverts,
            ebase,
            eivert,
            skip_adding=True))

    #indices  = np.array([np.where(u == biv)[0][0] for biv in ivert])
    #eindices = np.array([np.where(u == biv)[0][0] for biv in eivert])

    iv, ivi = np.unique(ivert, return_inverse=True)
    tmp = np.where(np.in1d(u, ivert, assume_unique=True))[0]
    indices = tmp[ivi]
    iv, ivi = np.unique(eivert, return_inverse=True)
    tmp = np.where(np.in1d(u, eivert, assume_unique=True))[0]
    eindices = tmp[ivi]

    Nvert = len(vtx)
    Nelem = len(attrs)
    Nbelem = len(eattrs)

    if myid == 0:
        print("NV, NBE, NE: " +
              ",".join([str(x) for x in (Nvert, Nbelem, Nelem)]))

    omesh = mfem.Mesh(2, Nvert, Nelem, Nbelem, sdim)

    _fill_mesh_elements(omesh, vtx, indices, nverts, attrs, base)
    _fill_mesh_bdr_elements(omesh, vtx, eindices, neverts, eattrs, ebase,
                            keelem)

    omesh.FinalizeTopology()
    omesh.Finalize(refine=True, fix_orientation=True)

    if hasNodal:
        odim = omesh.Dimension()
        print("odim, dim, sdim", odim, " ", dim, " ", sdim)
        fec = Nodal.FEColl()
        dNodal = mfem.FiniteElementSpace(omesh, fec, sdim)
        omesh.SetNodalFESpace(dNodal)
        omesh._nodal = dNodal

        if sdim == 3:
            if dim == 3:
                GetXDofs = Nodal.GetBdrElementDofs
                GetNX = Nodal.GetNBE
            elif dim == 2:
                GetXDofs = Nodal.GetElementDofs
                GetNX = Nodal.GetNE
            else:
                assert False, "not supported ndim 1"
            if odim == 3:
                dGetXDofs = dNodal.GetBdrElementDofs
                dGetNX = dNodal.GetNBE
            elif odim == 2:
                dGetXDofs = dNodal.GetElementDofs
                dGetNX = dNodal.GetNE
            else:
                assert False, "not supported ndim (3->1)"
        elif sdim == 2:
            GetNX = Nodal.GetNE
            dGetNX = dNodal.GetNE
            GetXDofs = Nodal.GetElementDofs
            dGetXDofs = dNodal.GetElementDofs

        DofToVDof = Nodal.DofToVDof
        dDofToVDof = dNodal.DofToVDof

        #nicePrint(dGetNX(),',', GetNX())
        nodes = mesh.GetNodes()
        node_ptx1 = nodes.GetDataArray()

        onodes = omesh.GetNodes()
        node_ptx2 = onodes.GetDataArray()
        #nicePrint(len(idx), idx)

        if len(idx) > 0:
            dof1_idx = np.hstack([[DofToVDof(i, d) for d in range(sdim)]
                                  for j in idx for i in GetXDofs(j)])
            data = node_ptx1[dof1_idx]
        else:
            dof1_idx = np.array([])
            data = np.array([])
        if isParMesh(mesh): data = allgather_vector(data)
        if isParMesh(mesh): idx = allgather_vector(idx)
        #nicePrint(len(data), ',', len(idx))

        dof2_idx = np.hstack([[dDofToVDof(i, d) for d in range(sdim)]
                              for j in range(len(idx)) for i in dGetXDofs(j)])
        node_ptx2[dof2_idx] = data
        #nicePrint(len(dof2_idx))

    if isParMesh(mesh):
        omesh = mfem.ParMesh(comm, omesh)

    if filename != '':
        if isParMesh(mesh):
            smyid = '{:0>6d}'.format(myid)
            filename = filename + '.' + smyid
        omesh.PrintToFile(filename, precision)

    return omesh
Ejemplo n.º 5
0
def edge(mesh, in_attr, filename='', precision=8):
    '''
    make a new mesh which contains only spedified edges.

    in_attr : eihter
    filename : an option to save the file 
    return new surface mesh
    '''
    sdim = mesh.SpaceDimension()
    dim = mesh.Dimension()
    Nodal = mesh.GetNodalFESpace()
    hasNodal = (Nodal is not None)

    if sdim == 3 and dim == 3:
        mode = 'edge', 'vertex'
    elif sdim == 3 and dim == 2:
        mode = 'bdr', 'vertex'
    elif sdim == 2 and dim == 2:
        mode = 'bdr', 'vertex'
    elif sdim == 2 and dim == 1:
        mode = 'dom', 'vertex'
    else:
        assert False, "unsupported mdoe"

    idx, attrs, ivert, nverts, base = _collect_data(in_attr, mesh, mode[0])

    l2v = mesh.extended_connectivity['line2vert']
    in_eattr = np.unique(np.hstack([l2v[k]
                                    for k in in_attr])).astype(int, copy=False)
    if isParMesh(mesh):
        in_eattr = np.unique(allgather_vector(in_eattr))
    eidx, eattrs, eivert, neverts, ebase = _collect_data(
        in_eattr, mesh, mode[1])

    u, indices = np.unique(np.hstack((ivert, eivert)), return_inverse=True)
    keelem = np.array([True] * len(eidx), dtype=bool)
    u_own = u

    if isParMesh(mesh):
        shared_info = distribute_shared_entity(mesh)
        u_own, ivert, eivert = _gather_shared_vertex(mesh, u, shared_info,
                                                     ivert, eivert)
    Nvert = len(u)
    if len(u_own) > 0:
        vtx = np.vstack([mesh.GetVertexArray(i) for i in u_own])
    else:
        vtx = np.array([]).reshape((-1, sdim))

    if isParMesh(mesh):
        #
        # distribute vertex/element data
        #
        base = allgather_vector(base)
        nverts = allgather_vector(nverts)
        attrs = allgather_vector(attrs)

        ivert = allgather_vector(ivert)
        eivert = allgather_vector(eivert)

        vtx = allgather_vector(vtx.flatten()).reshape(-1, sdim)

        u, indices = np.unique(np.hstack([ivert, eivert]), return_inverse=True)

        #
        # take care of shared boundary (edge)
        #
        keelem, eattrs, neverts, ebase, eivert = (_gather_shared_element(
            mesh, 'vertex', shared_info, eidx, keelem, eattrs, neverts, ebase,
            eivert))

    indices = np.array([np.where(u == biv)[0][0] for biv in ivert])
    eindices = np.array([np.where(u == biv)[0][0] for biv in eivert])

    Nvert = len(vtx)
    Nelem = len(attrs)
    Nbelem = len(eattrs)

    dprint1("NV, NBE, NE: " +
            ",".join([str(x) for x in (Nvert, Nbelem, Nelem)]))

    omesh = mfem.Mesh(1, Nvert, Nelem, Nbelem, sdim)

    _fill_mesh_elements(omesh, vtx, indices, nverts, attrs, base)
    _fill_mesh_bdr_elements(omesh, vtx, eindices, neverts, eattrs, ebase,
                            keelem)

    omesh.FinalizeTopology()

    if hasNodal:
        odim = omesh.Dimension()

        dprint1("odim, dim, sdim", odim, " ", dim, " ", sdim)
        fec = Nodal.FEColl()
        dNodal = mfem.FiniteElementSpace(omesh, fec, sdim)
        omesh.SetNodalFESpace(dNodal)
        omesh._nodal = dNodal

        GetXDofs = Nodal.GetElementDofs
        if dim == 3:
            GetXDofs = Nodal.GetEdgeDofs
        elif dim == 2:
            GetXDofs = Nodal.GetBdrElementDofs
        elif dim == 1:
            GetXDofs = Nodal.GetElementDofs

        dGetXDofs = dNodal.GetElementDofs

        DofToVDof = Nodal.DofToVDof
        dDofToVDof = dNodal.DofToVDof

        #nicePrint(dGetNX(),',', GetNX())
        nodes = mesh.GetNodes()
        node_ptx1 = nodes.GetDataArray()

        onodes = omesh.GetNodes()
        node_ptx2 = onodes.GetDataArray()
        #nicePrint(len(idx), idx)

        if len(idx) > 0:
            dof1_idx = np.hstack([[DofToVDof(i, d) for d in range(sdim)]
                                  for j in idx for i in GetXDofs(j)])
            data = node_ptx1[dof1_idx]
        else:
            dof1_idx = np.array([])
            data = np.array([])
        if isParMesh(mesh): data = allgather_vector(data)
        if isParMesh(mesh): idx = allgather_vector(idx)
        #nicePrint(len(data), ',', len(idx))

        dof2_idx = np.hstack([[dDofToVDof(i, d) for d in range(sdim)]
                              for j in range(len(idx)) for i in dGetXDofs(j)])
        node_ptx2[dof2_idx] = data
        #nicePrint(len(dof2_idx))

    # this should be after setting HO nodals...
    omesh.Finalize(refine=True, fix_orientation=True)

    if isParMesh(mesh):
        if omesh.GetNE() < nprc * 3:
            parts = omesh.GeneratePartitioning(1, 1)
        else:
            parts = None
        omesh = mfem.ParMesh(comm, omesh, parts)

    if filename != '':
        if isParMesh(mesh):
            smyid = '{:0>6d}'.format(myid)
            filename = filename + '.' + smyid
        omesh.PrintToFile(filename, precision)

    return omesh
Ejemplo n.º 6
0
def find_loop_par(pmesh, *face):
    '''
    find_loop_ser(mesh, 1)          # loop around boundary index = 1
    find_loop_ser(mesh, [1, 2, 3])    # loop around boundary made by union of index = 1,2,3
    '''

    import mfem.par as mfem
    from mpi4py import MPI

    myid = MPI.COMM_WORLD.rank
    nprc = MPI.COMM_WORLD.size
    comm = MPI.COMM_WORLD

    if nprc == 1:
        return find_loop_ser(pmesh, face)

    from mfem.common.mpi_debug import nicePrint
    from petram.helper.mpi_recipes import allgather, allgather_vector, gather_vector

    battrs = pmesh.GetBdrAttributeArray()

    face = face[0]
    faces = np.atleast_1d(face)
    bidx = np.where(np.in1d(battrs, faces))[0]

    offset_e = np.hstack([0, np.cumsum(allgather(pmesh.GetNEdges()))])

    edges = [pmesh.GetBdrElementEdges(i) for i in bidx]
    iedges = np.array(sum([e1[0]
                           for e1 in edges], []), dtype=int) + offset_e[myid]
    dirs = np.array(sum([e1[1] for e1 in edges], []), dtype=int)

    from petram.mesh.mesh_utils import distribute_shared_entity

    shared_info = distribute_shared_entity(pmesh)

    keys = shared_info[0].keys()
    local_edges = np.hstack([shared_info[0][key][1] for key in keys])
    global_edges = np.hstack([shared_info[1][key][1] for key in keys])

    own_edges = []
    for k, ie in enumerate(iedges):
        iii = np.where(local_edges == ie)[0]
        if len(iii) != 0:
            if ie == global_edges[iii[0]]:
                own_edges.append(ie)
            iedges[k] = global_edges[iii[0]]
        else:
            own_edges.append(ie)

    #nicePrint("iedges", iedges)
    iedges_all = allgather_vector(iedges)
    dirs = allgather_vector(dirs)

    seens = defaultdict(int)
    seendirs = defaultdict(int)
    for k, ie in enumerate(iedges_all):
        seens[ie] = seens[ie] + 1
        seendirs[ie] = dirs[k]
    seens.default_factory = None

    idx = []
    signs = []
    for k in seens.keys():
        if seens[k] == 1:
            idx.append(k)
            signs.append(seendirs[k])
    iedges_g = np.array(idx)
    # here idx is global numbering

    #nicePrint("global_index", idx, signs)
    #nicePrint("own_edges", own_edges)

    iedges_l = []
    signs_l = []
    for k, ie in enumerate(iedges_g):
        iii = np.where(own_edges == ie)[0]
        if len(iii) != 0:
            iedges_l.append(ie)
            signs_l.append(signs[k])
    iedges_l = np.array(iedges_l) - offset_e[myid]
    signs_l = np.array(signs_l)

    #nicePrint("local_index", iedges_l, signs_l)

    return iedges_l, signs_l