def create_fs_dm(V, problem=None):
    comm = V.mesh().mpi_comm()

    # this way the DM knows the function space it comes from
    dm = PETSc.DMShell().create(comm=comm)
    dm.setAttr('__fs__', V)
    dm.setAttr('__problem__', problem)

    # this gives the DM a template to create vectors inside snes
    dm.setGlobalVector(as_backend_type(Function(V).vector()).vec())

    if backend.__version__ > "2017.1.0":
        # this tells the DM how to interpolate from mesh to mesh
        # it depends on DOLFIN > 2017.1.0
        dm.setCreateInterpolation(create_interpolation)

    # if we have a mixed function space, then we need to tell PETSc
    # how to divvy up the different parts of the function space.
    # This is not needed for non-mixed elements.
    try:
        ufl_el = V.ufl_element()
        if isinstance(ufl_el, (MixedElement, VectorElement)):
            dm.setCreateSubDM(create_subdm)
            dm.setCreateFieldDecomposition(create_field_decomp)
    except AttributeError:
        # version of petsc4py is too old
        pass

    return dm
示例#2
0
    def setup_fields(self):

        if self.log_level >= 1:
            timer = df.Timer("pFibs: Setup fields")

        ## Default if empty ##
        if not self.block_field:
            for i in range(self.V.num_sub_spaces()):
                self.block_field.update({i: [i, i, {}]})
            self.num_fields = len(self.block_field)
        if self.num_fields == 0:
            self.num_fields += 1
        ## Create PetscSection ##
        self.section = PETSc.Section().create()
        self.section.setNumFields(self.num_fields)
        self.section.setChart(0, len(self.V.dofmap().dofs()))

        if self.log_level >= 2:
            timer_iterBlockFields = df.Timer(
                "pFibs: Setup fields - Iterate through block fields")

        ## Iterate through all the block fields ##
        for key in self.block_field:
            self.section.setFieldName(self.block_field[key][0], str(key))

            ## Extract dofs ##
            (dofs, ndof) = self.extract_dofs(key)

            ## Record dof count for each field ##
            self.field_size.update({self.block_field[key][0]: ndof})

            if self.log_level >= 3:
                timer_assignDof = df.Timer(
                    "pFibs: Setup fields - Iterate through block fields - assign dof"
                )

            assign_dof(self.section, dofs, self.goffset,
                       self.block_field[key][0])

            if self.log_level >= 3:
                timer_assignDof.stop()

        if self.log_level >= 2:
            timer_iterBlockFields.stop()

        ## Create DM and assign PetscSection ##
        self.section.setUp()
        self.dm = PETSc.DMShell().create()
        self.dm.setDefaultSection(self.section)
        self.dm.setUp()

        ## Prevent any further modification to block_field ##
        self.finalize_field = True

        if self.log_level >= 1:
            timer.stop()
示例#3
0
 def testCoarsenRefine(self):
     cdm = PETSc.DMShell().create(comm=self.COMM)
     def coarsen(dm, comm):
         return cdm
     def refine(dm, comm):
         return self.dm
     cdm.setRefine(refine)
     self.dm.setCoarsen(coarsen)
     coarsened = self.dm.coarsen()
     self.assertEqual(coarsened, cdm)
     refined = coarsened.refine()
     self.assertEqual(refined, self.dm)
示例#4
0
    def _setup_snes(self, pressure_mesh, residual_size):
        # Creates the Jacobian matrix structure.
        j_structure = _calculate_jacobian_mask(pressure_mesh.nx, pressure_mesh.ny, 3)
        logger.info("Jacobian NNZ=%s" % (j_structure.nnz,))
        csr = (j_structure.indptr, j_structure.indices, j_structure.data)
        self._petsc_jacobian = PETSc.Mat().createAIJWithArrays(j_structure.shape, csr)
        self._petsc_jacobian.assemble(assembly=self._petsc_jacobian.AssemblyType.FINAL_ASSEMBLY)

        self._comm = PETSc.COMM_WORLD
        self._dm = PETSc.DMShell().create(comm=self._comm)
        self._dm.setMatrix(self._petsc_jacobian)

        # residual vector
        self._r = PETSc.Vec().createSeq(residual_size)

        self._snes = PETSc.SNES().create(comm=self._comm)
        self._snes.setFunction(self.residual_function_for_petsc, self._r)
        self._snes.setDM(self._dm)
        self._snes.setConvergenceHistory()
        self._snes.setFromOptions()
        self._snes.setUseFD(True)
        self._snes.setMonitor(self._solver_monitor)
        self._snes.ksp.setMonitor(self._linear_solver_monitor)
        self._snes.setTolerances(rtol=1e-4, atol=1e-4, stol=1e-4, max_it=50)
示例#5
0
 def dm(self):
     dm = PETSc.DMShell().create(comm=self.comm)
     dm.setGlobalVector(self.layout_vec)
     return dm
示例#6
0
 def setUp(self):
     self.dm = PETSc.DMShell().create(comm=self.COMM)
示例#7
0
    for i in range(N):
        f[i] = (x[i] ** 2 - 9.0).item()
    f[N // 2] = x[N // 2]
    f.assemble()

COMM = PETSc.COMM_WORLD
N = 4
J = PETSc.Mat().createAIJ(N, comm=COMM)
J.setPreallocationNNZ(N)
for i in range(N):
    for j in range(N):
        J.setValue(i, j, 0.0)
J.setUp()
J.assemble()

dm = PETSc.DMShell().create(comm=COMM)
dm.setMatrix(J)

snes = PETSc.SNES().create(comm=COMM)
r = PETSc.Vec().createSeq(N)
x = PETSc.Vec().createSeq(N)
b = PETSc.Vec().createSeq(N)
snes.setFunction(f, r)
snes.setDM(dm)

# This enables coloring
snes.setUseFD(True)

x.setArray([1.1] * N)
b.set(0)