예제 #1
0
    def __call__(self, pc):
        dm = pc.getDM()
        prefix = pc.getOptionsPrefix()
        sentinel = object()
        sweeps = PETSc.Options(prefix).getString(
            "pc_patch_construct_ps_sweeps", default=sentinel)
        if sweeps == sentinel:
            raise ValueError("Must set %spc_patch_construct_ps_sweeps" %
                             prefix)

        patches = []
        for sweep in sweeps.split(':'):
            axis = int(sweep[0])
            dir = {'+': +1, '-': -1}[sweep[1]]
            ndiv = int(sweep[2:])

            entities = self.sort_entities(dm, axis, dir, ndiv)
            for patch in entities:
                iset = PETSc.IS().createGeneral(patch, comm=PETSc.COMM_SELF)
                patches.append(iset)

        iterationSet = PETSc.IS().createStride(size=len(patches),
                                               first=0,
                                               step=1,
                                               comm=PETSc.COMM_SELF)
        return (patches, iterationSet)
예제 #2
0
    def __call__(self, pc):
        dm = pc.getDM()
        prefix = pc.getOptionsPrefix()
        opts = PETSc.Options(prefix)

        select = partial(self.select_entity, dm=dm, exclude="pyop2_ghost")
        (pStart, pEnd) = dm.getChart()
        owned_entities = set(filter(select, range(pStart, pEnd)))

        overlap = 1  # FIXME: get from options
        all_entities = set(owned_entities)
        for i in range(overlap):
            tmp_entities = set(all_entities)
            for entity in tmp_entities:
                for s in self.star(dm, entity):
                    for c in self.closure(dm, s):
                        all_entities.add(c)

        overlap_entities = list(all_entities.difference(owned_entities))

        iset = PETSc.IS().createGeneral(overlap_entities, comm=PETSc.COMM_SELF)
        patches = [iset]
        iterset = PETSc.IS().createStride(size=len(patches),
                                          first=0,
                                          step=1,
                                          comm=PETSc.COMM_SELF)

        #print("owned_entities: ", owned_entities)
        #print("all_entities: ", all_entities)
        #print("overlap_entities: ", overlap_entities)
        return (patches, iterset)
예제 #3
0
파일: relaxation.py 프로젝트: wence-/ssc
    def __call__(self, pc):
        dm = pc.getDM()
        prefix = pc.getOptionsPrefix()
        sentinel = object()
        opts = PETSc.Options(prefix)
        name = self.name
        assert self.name is not None

        self.set_options(dm, opts, name)

        select = partial(select_entity, dm=dm, exclude="pyop2_ghost")
        entities = list(filter(select, self.get_entities(opts, name, dm)))

        nclosure = opts.getInt("pc_patch_construction_%s_nclosures" % name, default=1)
        patches = []
        for entity in entities:
            subentities = self.callback(dm, entity, nclosure)
            iset = PETSc.IS().createGeneral(subentities, comm=PETSc.COMM_SELF)
            patches.append(iset)

        # Now make the iteration set.
        iterset = []

        sortorders = opts.getString("pc_patch_construction_%s_sort_order" % name, default=sentinel)
        if sortorders == sentinel:
            sortorders = "None"

        if sortorders == "None":
            piterset = PETSc.IS().createStride(size=len(patches), first=0, step=1, comm=PETSc.COMM_SELF)
            return (patches, piterset)

        coords = list(enumerate(self.coords(dm, p) for p in entities))
        for sortorder in sortorders.split("|"):
            sortdata = []
            for axis in sortorder.split(':'):
                ax = int(axis[0])
                if len(axis) > 1:
                    sgn = {'+': 1, '-': -1}[axis[1]]
                else:
                    sgn = 1
                sortdata.append((ax, sgn))

            def keyfunc(z):
                return tuple(sgn*z[1][ax] for (ax, sgn) in sortdata)

            iterset += [x[0] for x in sorted(coords, key=keyfunc)]

        piterset = PETSc.IS().createGeneral(iterset, comm=PETSc.COMM_SELF)
        return (patches, piterset)
예제 #4
0
 def create_subdm(cls, dm, fields, *args, **kwargs):
     W = dm.getAttr('__fs__')()
     if len(fields) == 1:
         # Subspace is just a single FunctionSpace.
         field = fields[0]
         subdm = W[field]._dm
         iset = W._ises[field]
         return iset, subdm
     else:
         try:
             # Look up the subspace in the cache
             iset, subspace = W._subspaces[tuple(fields)]
             return iset, subspace._dm
         except KeyError:
             pass
         # Need to build an MFS for the subspace
         subspace = MixedFunctionSpace([W[f] for f in fields])
         # Index set mapping from W into subspace.
         iset = PETSc.IS().createGeneral(
             numpy.concatenate([W._ises[f].indices for f in fields]))
         # Keep hold of strong reference to created subspace (given we
         # only hold a weakref in the shell DM), and so we can
         # reuse it later.
         W._subspaces[tuple(fields)] = iset, subspace
         return iset, subspace._dm
예제 #5
0
    def get_patches(self, V):
        mesh = V._mesh
        assert mesh.cell_set._extruded
        dm = mesh._topology_dm
        section = V.dm.getDefaultSection()
        # Obtain the codimensions to loop over from options, if present
        codim_list = PETSc.Options().getString(self.prefix + "codims", "0, 1")
        codim_list = [int(ii) for ii in codim_list.split(",")]

        # Build index sets for the patches
        ises = []
        for codim in codim_list:
            for p in range(*dm.getHeightStratum(codim)):
                # Only want to build patches over owned faces
                if dm.getLabelValue("pyop2_ghost", p) != -1:
                    continue
                dof = section.getDof(p)
                if dof <= 0:
                    continue
                off = section.getOffset(p)
                indices = numpy.arange(off * V.value_size,
                                       V.value_size * (off + dof),
                                       dtype='int32')
                iset = PETSc.IS().createGeneral(indices, comm=COMM_SELF)
                ises.append(iset)

        return ises
예제 #6
0
def create_subdm(dm, fields, *args, **kwargs):
    """Callback to create a sub-DM describing the specified fields.

    :arg DM: The DM.
    :arg fields: The fields in the new sub-DM.
    """
    W = get_function_space(dm)
    ctx = get_appctx(dm)
    coarsen = get_ctx_coarsener(dm)
    if len(fields) == 1:
        # Subspace is just a single FunctionSpace.
        idx, = fields
        subdm = W[idx].dm
        iset = W._ises[idx]
        if ctx is not None:
            ctx, = ctx.split([(idx, )])
            push_appctx(subdm, ctx)
            push_ctx_coarsener(subdm, coarsen)
        return iset, subdm
    else:
        # Need to build an MFS for the subspace
        subspace = firedrake.MixedFunctionSpace([W[f] for f in fields])
        # Index set mapping from W into subspace.
        iset = PETSc.IS().createGeneral(numpy.concatenate(
            [W._ises[f].indices for f in fields]),
                                        comm=W.comm)
        if ctx is not None:
            ctx, = ctx.split([fields])
            push_appctx(subspace.dm, ctx)
            push_ctx_coarsener(subspace.dm, coarsen)
        return iset, subspace.dm
예제 #7
0
    def get_patches(self, V):
        mesh = V._mesh
        mesh_dm = mesh._topology_dm

        # Obtain the topological entities to use to construct the stars
        depth = PETSc.Options().getInt(self.prefix + "construct_dim",
                                       default=-1)
        height = PETSc.Options().getInt(self.prefix + "construct_codim",
                                        default=-1)
        if (depth == -1 and height == -1) or (depth != -1 and height != -1):
            raise ValueError(
                f"Must set exactly one of {self.prefix}construct_dim or {self.prefix}construct_codim"
            )

        # Accessing .indices causes the allocation of a global array,
        # so we need to cache these for efficiency
        V_local_ises_indices = []
        for (i, W) in enumerate(V):
            V_local_ises_indices.append(V.dof_dset.local_ises[i].indices)

        # Build index sets for the patches
        ises = []
        if depth != -1:
            (start, end) = mesh_dm.getDepthStratum(depth)
        else:
            (start, end) = mesh_dm.getHeightStratum(height)

        for seed in range(start, end):
            # Only build patches over owned DoFs
            if mesh_dm.getLabelValue("pyop2_ghost", seed) != -1:
                continue

            # Create point list from mesh DM
            star, _ = mesh_dm.getTransitiveClosure(seed, useCone=False)
            pt_array = set()
            for pt in star.tolist():
                closure, _ = mesh_dm.getTransitiveClosure(seed, useCone=True)
                pt_array.update(closure.tolist())

            # Get DoF indices for patch
            indices = []
            for (i, W) in enumerate(V):
                section = W.dm.getDefaultSection()
                for p in pt_array:
                    dof = section.getDof(p)
                    if dof <= 0:
                        continue
                    off = section.getOffset(p)
                    # Local indices within W
                    W_indices = numpy.arange(off * W.value_size,
                                             W.value_size * (off + dof),
                                             dtype='int32')
                    indices.extend(V_local_ises_indices[i][W_indices])
            iset = PETSc.IS().createGeneral(indices, comm=COMM_SELF)
            ises.append(iset)

        return ises
예제 #8
0
    def __init__(self, Q, free_bids=["on_boundary"]):
        (V, I_interp) = Q.get_space_for_inner()

        self.free_bids = free_bids

        u = fd.TrialFunction(V)
        v = fd.TestFunction(V)

        n = fd.FacetNormal(V.mesh())

        def surf_grad(u):
            return fd.sym(fd.grad(u) - fd.outer(fd.grad(u) * n, n))

        a = (fd.inner(surf_grad(u), surf_grad(v)) + fd.inner(u, v)) * fd.ds
        # petsc doesn't like matrices with zero rows
        a += 1e-10 * fd.inner(u, v) * fd.dx
        A = fd.assemble(a, mat_type="aij")
        A = A.petscmat
        tdim = V.mesh().topological_dimension()

        lsize = fd.Function(V).vector().local_size()

        def get_nodes_bc(bc):
            nodes = bc.nodes
            return nodes[nodes < lsize]

        def get_nodes_bid(bid):
            bc = fd.DirichletBC(V, fd.Constant(tdim * (0, )), bid)
            return get_nodes_bc(bc)

        free_nodes = np.concatenate(
            [get_nodes_bid(bid) for bid in self.free_bids])
        free_dofs = np.concatenate(
            [tdim * free_nodes + i for i in range(tdim)])
        free_dofs = np.unique(np.sort(free_dofs))
        self.free_is = PETSc.IS().createGeneral(free_dofs)
        lgr, lgc = A.getLGMap()
        self.global_free_is_row = lgr.applyIS(self.free_is)
        self.global_free_is_col = lgc.applyIS(self.free_is)
        A = A.createSubMatrix(self.global_free_is_row, self.global_free_is_col)
        # A.view()
        A.assemble()
        self.A = A
        Aksp = PETSc.KSP().create()
        Aksp.setOperators(self.A)
        Aksp.setOptionsPrefix("A_")
        opts = PETSc.Options()
        opts["A_ksp_type"] = "cg"
        opts["A_pc_type"] = "hypre"
        opts["A_ksp_atol"] = 1e-10
        opts["A_ksp_rtol"] = 1e-10
        Aksp.setUp()
        Aksp.setFromOptions()
        self.Aksp = Aksp
예제 #9
0
    def copy_massmatrix_into_dat(self):

        # Copy the velocity mass matrix into a Dat
        vmat = self.imass_velocity.handle
        dofs_per_entity = self.U.fiat_element.entity_dofs()
        dofs_per_entity = sum(
            self.mesh.make_dofs_per_plex_entity(dofs_per_entity))
        arity = dofs_per_entity * self.U.topological.dim
        self.velocity_mass_asdat = Dat(DataSet(self.mesh.cell_set,
                                               arity * arity),
                                       dtype='double')
        istart, iend = vmat.getOwnershipRange()
        idxs = [
            PETSc.IS().createGeneral(np.arange(i, i + arity, dtype=np.int32),
                                     comm=PETSc.COMM_SELF)
            for i in range(istart, iend, arity)
        ]
        submats = vmat.getSubMatrices(idxs, idxs)
        for i, m in enumerate(submats):
            self.velocity_mass_asdat.data[i] = m[:, :].flatten()
        info("Computed velocity mass matrix")

        # Copy the stress mass matrix into a Dat
        smat = self.imass_stress.handle
        dofs_per_entity = self.S.fiat_element.entity_dofs()
        dofs_per_entity = sum(
            self.mesh.make_dofs_per_plex_entity(dofs_per_entity))
        arity = dofs_per_entity * self.S.topological.dim
        self.stress_mass_asdat = Dat(DataSet(self.mesh.cell_set,
                                             arity * arity),
                                     dtype='double')
        istart, iend = smat.getOwnershipRange()
        idxs = [
            PETSc.IS().createGeneral(np.arange(i, i + arity, dtype=np.int32),
                                     comm=PETSc.COMM_SELF)
            for i in range(istart, iend, arity)
        ]
        submats = smat.getSubMatrices(idxs, idxs)
        for i, m in enumerate(submats):
            self.stress_mass_asdat.data[i] = m[:, :].flatten()
        info("Computed stress mass matrix")
예제 #10
0
    def gather(self, global_indices=None):
        """Gather a :class:`Vector` to all processes

        :arg global_indices: the globally numbered indices to gather
                            (should be the same on all processes).  If
                            `None`, gather the entire :class:`Vector`."""
        if global_indices is None:
            N = self.size()
            v = PETSc.Vec().createSeq(N, comm=PETSc.COMM_SELF)
            is_ = PETSc.IS().createStride(N, 0, 1, comm=PETSc.COMM_SELF)
        else:
            global_indices = np.asarray(global_indices, dtype=np.int32)
            N = len(global_indices)
            v = PETSc.Vec().createSeq(N, comm=PETSc.COMM_SELF)
            is_ = PETSc.IS().createGeneral(global_indices, comm=PETSc.COMM_SELF)

        with self.dat.vec_ro as vec:
            vscat = PETSc.Scatter().create(vec, is_, v, None)
            vscat.scatterBegin(vec, v, addv=PETSc.InsertMode.INSERT_VALUES)
            vscat.scatterEnd(vec, v, addv=PETSc.InsertMode.INSERT_VALUES)
        return v.array
예제 #11
0
파일: elastic.py 프로젝트: shineusn/seigen
 def matrix_to_dat(self, massmatrix, functionspace):
     r""" Copy the clock diagonal entries of the velocity mass
     matrix into a pyop2.Dat"""
     arity = sum(functionspace.topological.dofs_per_entity)*functionspace.topological.dim
     dat = Dat(DataSet(self.mesh.cell_set, arity*arity), dtype='double')
     istart, iend = massmatrix.handle.getOwnershipRange()
     idxs = [PETSc.IS().createGeneral(np.arange(i, i+arity, dtype=np.int32),
                                      comm=PETSc.COMM_SELF)
             for i in range(istart, iend, arity)]
     submats = massmatrix.handle.getSubMatrices(idxs, idxs)
     for i, m in enumerate(submats):
         dat.data[i] = m[:, :].flatten()
     return dat
예제 #12
0
def matern(V, variance=1, smoothness=1, correlation_length=1, rg=None):
    d = V.mesh().topological_dimension()
    nu = smoothness
    sigma = np.sqrt(variance)
    lambd = correlation_length

    k = (nu + d / 2) / 2
    assert k == 1
    kappa = np.sqrt(8) / lambd
    sigma_hat = gamma(nu) * nu**(d / 2)
    sigma_hat /= gamma(nu + d / 2)
    sigma_hat *= (2 / np.pi)**(d / 2)
    sigma_hat *= lambd**(-d)
    eta = sigma / sigma_hat

    if rg is None:
        pcg = PCG64(seed=100)
        rg = RandomGenerator(pcg)

    Wnoise = rg.normal(V, 0.0, 1.0)
    HWnoise = Function(V)

    u = TrialFunction(V)
    v = TestFunction(V)

    #bcs = DirichletBC(V, 0, (1,2,3,4))

    mass = inner(u, v) * dx

    M = assemble(mass).M
    print(M.handle.type)
    iset = PETSc.IS().createGeneral(range(M.nrows), comm=COMM_SELF)
    M.handle.factorCholesky(iset)

    M = get_np_matrix(mass)
    H = cholesky(M)
    HWnoise.dat.data[:] = H @ Wnoise.dat.data

    a = (inner(u, v) + Constant(1 / (kappa**2)) * inner(grad(u), grad(v))) * dx
    l = Constant(eta) * inner(HWnoise, v) * dx

    u_h = Function(V)
    solve(a == l,
          u_h,
          bcs=bcs,
          solver_parameters={
              'ksp_type': 'cg',
              'pc_type': 'gamg'
          })

    return u_h
예제 #13
0
def find_sub_block(iset, ises):
    """Determine if iset comes from a concatenation of some subset of
    ises.

    :arg iset: a PETSc IS to find in ``ises``.
    :arg ises: An iterable of PETSc ISes.

    :returns: The indices into ``ises`` that when concatenated
        together produces ``iset``.

    :raises LookupError: if ``iset`` could not be found in
        ``ises``.
    """
    found = []
    sfound = set()
    comm = iset.comm
    while True:
        match = False
        for i, iset_ in enumerate(ises):
            if i in sfound:
                continue
            lsize = iset_.getLocalSize()
            if lsize > iset.getLocalSize():
                continue
            indices = iset.indices
            tmp = PETSc.IS().createGeneral(indices[:lsize], comm=comm)
            if tmp.equal(iset_):
                found.append(i)
                sfound.add(i)
                iset = PETSc.IS().createGeneral(indices[lsize:], comm=comm)
                match = True
                continue
        if not match:
            break
    if iset.getSize() > 0:
        raise LookupError("Unable to find %s in %s" % (iset, ises))
    return found
예제 #14
0
 def create_subdm(cls, dm, fields, *args, **kwargs):
     W = dm.getAttr('__fs__')()
     if len(fields) == 1:
         # Subspace is just a single FunctionSpace.
         subspace = W[fields[0]]
     else:
         # Need to build an MFS for the subspace
         subspace = MixedFunctionSpace([W[f] for f in fields])
     # Sub-DM is just the DM belonging to the subspace.
     subdm = subspace._dm
     # Keep hold of strong reference, to created subspace (given we
     # only hold a weakref in the shell DM)
     W._subspaces.append(subspace)
     # Index set mapping from W into subspace.
     iset = PETSc.IS().createGeneral(
         np.concatenate([W._ises[f].indices for f in fields]))
     return iset, subdm
예제 #15
0
    def get_patches(self, V):
        mesh = V._mesh
        mesh_dm = mesh.topology_dm
        if mesh.layers:
            warning("applying ASMStarPC on an extruded mesh")

        # Obtain the topological entities to use to construct the stars
        depth = PETSc.Options().getInt(self.prefix+"construct_dim", default=0)

        # Accessing .indices causes the allocation of a global array,
        # so we need to cache these for efficiency
        V_local_ises_indices = []
        for (i, W) in enumerate(V):
            V_local_ises_indices.append(V.dof_dset.local_ises[i].indices)

        # Build index sets for the patches
        ises = []
        (start, end) = mesh_dm.getDepthStratum(depth)
        for seed in range(start, end):
            # Only build patches over owned DoFs
            if mesh_dm.getLabelValue("pyop2_ghost", seed) != -1:
                continue

            # Create point list from mesh DM
            pt_array, _ = mesh_dm.getTransitiveClosure(seed, useCone=False)

            # Get DoF indices for patch
            indices = []
            for (i, W) in enumerate(V):
                section = W.dm.getDefaultSection()
                for p in pt_array.tolist():
                    dof = section.getDof(p)
                    if dof <= 0:
                        continue
                    off = section.getOffset(p)
                    # Local indices within W
                    W_indices = numpy.arange(off*W.value_size, W.value_size * (off + dof), dtype=IntType)
                    indices.extend(V_local_ises_indices[i][W_indices])
            iset = PETSc.IS().createGeneral(indices, comm=COMM_SELF)
            ises.append(iset)

        return ises
예제 #16
0
def create_subdm(dm, fields, *args, **kwargs):
    """Callback to create a sub-DM describing the specified fields.

    :arg DM: The DM.
    :arg fields: The fields in the new sub-DM.

    .. note::

       This should, but currently does not, transfer appropriately
       split application contexts onto the sub-DMs.
    """
    W = get_function_space(dm)
    # TODO: Correct splitting of SNESContext for len(fields) > 1 case
    if len(fields) == 1:
        # Subspace is just a single FunctionSpace.
        idx = fields[0]
        subdm = W[idx].dm
        iset = W._ises[idx]
        return iset, subdm
    else:
        try:
            # Look up the subspace in the cache
            iset, subspace = W._subspaces[tuple(fields)]
            return iset, subspace.dm
        except KeyError:
            pass
        # Need to build an MFS for the subspace
        subspace = firedrake.MixedFunctionSpace([W[f] for f in fields])
        # Index set mapping from W into subspace.
        iset = PETSc.IS().createGeneral(numpy.concatenate([W._ises[f].indices
                                                           for f in fields]),
                                        comm=W.comm)
        # Keep hold of strong reference to created subspace (given we
        # only hold a weakref in the shell DM), and so we can
        # reuse it later.
        W._subspaces[tuple(fields)] = iset, subspace
        return iset, subspace.dm
예제 #17
0
    def initialize(self, pc):
        _, P = pc.getOperators()
        assert P.type == "python"
        context = P.getPythonContext()
        (self.J, self.bcs) = (context.a, context.row_bcs)

        test, trial = self.J.arguments()
        if test.function_space() != trial.function_space():
            raise NotImplementedError("test and trial spaces must be the same")

        Pk = test.function_space()
        element = Pk.ufl_element()
        shape = element.value_shape()
        mesh = Pk.ufl_domain()
        if len(shape) == 0:
            P1 = firedrake.FunctionSpace(mesh, "CG", 1)
        elif len(shape) == 2:
            P1 = firedrake.VectorFunctionSpace(mesh, "CG", 1, dim=shape[0])
        else:
            P1 = firedrake.TensorFunctionSpace(mesh,
                                               "CG",
                                               1,
                                               shape=shape,
                                               symmetry=element.symmetry())

        # TODO: A smarter low-order operator would also interpolate
        # any coefficients to the coarse space.
        mapper = ArgumentReplacer({
            test: firedrake.TestFunction(P1),
            trial: firedrake.TrialFunction(P1)
        })
        self.lo_J = map_integrands.map_integrand_dags(mapper, self.J)

        lo_bcs = []
        for bc in self.bcs:
            # Don't actually need the value, since it's only used for
            # killing parts of the restriction matrix.
            lo_bcs.append(
                firedrake.DirichletBC(P1,
                                      firedrake.zero(P1.shape),
                                      bc.sub_domain,
                                      method=bc.method))

        self.lo_bcs = tuple(lo_bcs)

        mat_type = PETSc.Options().getString(
            pc.getOptionsPrefix() + "lo_mat_type",
            firedrake.parameters["default_matrix_type"])
        self.lo_op = firedrake.assemble(self.lo_J,
                                        bcs=self.lo_bcs,
                                        mat_type=mat_type)
        self.lo_op.force_evaluation()
        A, P = pc.getOperators()
        nearnullsp = P.getNearNullSpace()
        if nearnullsp.handle != 0:
            # Actually have a near nullspace
            tmp = firedrake.Function(Pk)
            low = firedrake.Function(P1)
            vecs = []
            for vec in nearnullsp.getVecs():
                with tmp.dat.vec as v:
                    vec.copy(v)
                low.interpolate(tmp)
                with low.dat.vec_ro as v:
                    vecs.append(v.copy())
            nullsp = PETSc.NullSpace().create(vectors=vecs, comm=pc.comm)
            self.lo_op.petscmat.setNearNullSpace(nullsp)
        lo = PETSc.PC().create(comm=pc.comm)
        lo.incrementTabLevel(1, parent=pc)
        lo.setOperators(self.lo_op.petscmat, self.lo_op.petscmat)
        lo.setOptionsPrefix(pc.getOptionsPrefix() + "lo_")
        lo.setFromOptions()
        self.lo = lo
        self.restriction = restriction_matrix(Pk, P1, self.bcs, self.lo_bcs)

        self.work = self.lo_op.petscmat.createVecs()
        if len(self.bcs) > 0:
            bc_nodes = numpy.unique(
                numpy.concatenate([bc.nodes for bc in self.bcs]))
            bc_nodes = bc_nodes[bc_nodes < Pk.dof_dset.size]
            bc_iset = PETSc.IS().createBlock(numpy.prod(shape),
                                             bc_nodes,
                                             comm=PETSc.COMM_SELF)
            self.bc_indices = bc_iset.getIndices()
            bc_iset.destroy()
        else:
            self.bc_indices = numpy.empty(0, dtype=numpy.int32)
예제 #18
0
def create_subdm(dm, fields, *args, **kwargs):
    """Callback to create a sub-DM describing the specified fields.

    :arg DM: The DM.
    :arg fields: The fields in the new sub-DM.
    """
    W = get_function_space(dm)
    ctx = get_appctx(dm)
    coarsen = get_ctx_coarsener(dm)
    parent = get_parent(dm)
    if len(fields) == 1:
        # Subspace is just a single FunctionSpace.
        idx, = fields
        subdm = W[idx].dm
        iset = W._ises[idx]
        add_hook(parent,
                 setup=partial(push_parent, subdm, parent),
                 teardown=partial(pop_parent, subdm, parent),
                 call_setup=True)

        if ctx is not None:
            ctx, = ctx.split([(idx, )])
            add_hook(parent,
                     setup=partial(push_appctx, subdm, ctx),
                     teardown=partial(pop_appctx, subdm, ctx),
                     call_setup=True)
            add_hook(parent,
                     setup=partial(push_ctx_coarsener, subdm, coarsen),
                     teardown=partial(pop_ctx_coarsener, subdm, coarsen),
                     call_setup=True)
        return iset, subdm
    else:
        # Need to build an MFS for the subspace
        subspace = firedrake.MixedFunctionSpace([W[f] for f in fields])

        # Pass any transfer operators over
        transfer = get_transfer_operators(dm)

        add_hook(parent,
                 setup=partial(push_transfer_operators, subspace.dm, transfer),
                 teardown=partial(pop_transfer_operators, subspace.dm,
                                  transfer),
                 call_setup=True)
        add_hook(parent,
                 setup=partial(push_parent, subspace.dm, parent),
                 teardown=partial(pop_parent, subspace.dm, parent),
                 call_setup=True)
        # Index set mapping from W into subspace.
        iset = PETSc.IS().createGeneral(numpy.concatenate(
            [W._ises[f].indices for f in fields]),
                                        comm=W.comm)
        if ctx is not None:
            ctx, = ctx.split([fields])
            add_hook(parent,
                     setup=partial(push_appctx, subspace.dm, ctx),
                     teardown=partial(pop_appctx, subspace.dm, ctx),
                     call_setup=True)
            add_hook(parent,
                     setup=partial(push_ctx_coarsener, subspace.dm, coarsen),
                     teardown=partial(pop_ctx_coarsener, subspace.dm, coarsen),
                     call_setup=True)
        return iset, subspace.dm
예제 #19
0
    def get_patches(self, V):
        mesh = V.mesh()
        nlayers = mesh.layers
        mesh_dm = mesh.topology_dm

        # Obtain the topological entities to use to construct the stars
        depth = PETSc.Options().getInt(self.prefix+"construct_dim",
                                       default=0)

        # Accessing .indices causes the allocation of a global array,
        # so we need to cache these for efficiency
        V_ises = tuple(iset.indices for iset in V.dof_dset.local_ises)
        basemeshoff = []
        basemeshdof = []
        basemeshlayeroffsets = []
        for (i, W) in enumerate(V):
            boff, bdof, blayer_offsets = get_basemesh_nodes(W)
            basemeshoff.append(boff)
            basemeshdof.append(bdof)
            basemeshlayeroffsets.append(blayer_offsets)

        # Build index sets for the patches
        ises = []
        start, end = mesh_dm.getDepthStratum(depth)
        pstart, _ = mesh_dm.getChart()
        for seed in range(start, end):
            # Only build patches over owned DoFs
            if mesh_dm.getLabelValue("pyop2_ghost", seed) != -1:
                continue

            # Create point list from mesh DM
            points, _ = mesh_dm.getTransitiveClosure(seed, useCone=False)
            points -= pstart  # offset by chart start
            for k in range(nlayers):
                indices = []
                # Get DoF indices for patch
                for i, W in enumerate(V):
                    iset = V_ises[i]
                    for p in points:
                        # How to walk up one layer
                        blayer_offset = basemeshlayeroffsets[i][p]
                        if blayer_offset <= 0:
                            # In this case we don't have any dofs on
                            # this entity.
                            continue
                        # Offset in the global array for the bottom of
                        # the column
                        off = basemeshoff[i][p]
                        # Number of dofs in the interior of the
                        # vertical interval cell on top of this base
                        # entity
                        dof = basemeshdof[i][p]
                        # Hard-code taking the star
                        if k == 0:
                            begin = off
                            end = off + blayer_offset
                        elif k < nlayers - 1:
                            begin = off + (k-1) * blayer_offset + dof
                            end = off + (k+1) * blayer_offset
                        else:  # k == nlayers - 1:
                            begin = off + (k-1) * blayer_offset + dof
                            end = off + k * blayer_offset + dof
                        zlice = slice(W.value_size * begin,
                                      W.value_size * end)
                        indices.extend(iset[zlice])
                iset = PETSc.IS().createGeneral(indices, comm=COMM_SELF)
                ises.append(iset)
        return ises
예제 #20
0
    def construct_1d_interpolation_matrices(self, V):
        """
        Create a list of sparse matrices (one per geometric dimension).

        Each matrix has size (M, n[dim]), where M is the dimension of the
        self.V_r.sub(0), and n[dim] is the dimension of the univariate
        spline space associated to the dimth-geometric coordinate.
        The ith column of such a matrix is computed by evaluating the ith
        univariate B-spline on the dimth-geometric coordinate of the dofs
        of self.V_r(0)
        """
        interp_1d = []

        # this code is correct but can be made more beautiful
        # by replacing x_fct with self.id
        x_fct = fd.SpatialCoordinate(self.mesh_r)  # used for x_int
        # compute self.M, x_int will be overwritten below
        x_int = fd.interpolate(x_fct[0], V.sub(0))
        self.M = x_int.vector().size()

        comm = self.comm

        u, v = fd.TrialFunction(V.sub(0)), fd.TestFunction(V.sub(0))
        mass_temp = fd.assemble(u * v * fd.dx)
        self.lg_map_fe = mass_temp.petscmat.getLGMap()[0]

        for dim in range(self.dim):

            order = self.orders[dim]
            knots = self.knots[dim]
            n = self.n[dim]

            # owned part of global problem
            local_n = n // comm.size + int(comm.rank < (n % comm.size))
            I = PETSc.Mat().create(comm=self.comm)
            I.setType(PETSc.Mat.Type.AIJ)
            lsize = x_int.vector().local_size()
            gsize = x_int.vector().size()
            I.setSizes(((lsize, gsize), (local_n, n)))

            I.setUp()
            x_int = fd.interpolate(x_fct[dim], V.sub(0))
            x = x_int.vector().get_local()
            for idx in range(n):
                coeffs = np.zeros(knots.shape, dtype=float)

                # impose boundary regularity
                coeffs[idx + self.boundary_regularities[dim]] = 1
                degree = order - 1  # splev uses degree, not order
                tck = (knots, coeffs, degree)

                values = splev(x, tck, der=0, ext=1)
                rows = np.where(values != 0)[0].astype(np.int32)
                values = values[rows]
                rows_is = PETSc.IS().createGeneral(rows)
                global_rows_is = self.lg_map_fe.applyIS(rows_is)
                rows = global_rows_is.array
                I.setValues(rows, [idx], values)

            I.assemble()  # lazy strategy for kron
            interp_1d.append(I)

        # from IPython import embed; embed()
        return interp_1d