xref: /petsc/src/binding/petsc4py/src/petsc4py/PETSc/DM.pyx (revision b5f0bcd6e9e8ed97648738542f5163d94f7b1782)
1# --------------------------------------------------------------------
2
3class DMType(object):
4    """`DM` types."""
5    DA        = S_(DMDA_type)
6    COMPOSITE = S_(DMCOMPOSITE)
7    SLICED    = S_(DMSLICED)
8    SHELL     = S_(DMSHELL)
9    PLEX      = S_(DMPLEX)
10    REDUNDANT = S_(DMREDUNDANT)
11    PATCH     = S_(DMPATCH)
12    MOAB      = S_(DMMOAB)
13    NETWORK   = S_(DMNETWORK)
14    FOREST    = S_(DMFOREST)
15    P4EST     = S_(DMP4EST)
16    P8EST     = S_(DMP8EST)
17    SWARM     = S_(DMSWARM)
18    PRODUCT   = S_(DMPRODUCT)
19    STAG      = S_(DMSTAG)
20
21
22class DMBoundaryType(object):
23    """`DM` Boundary types."""
24    NONE     = DM_BOUNDARY_NONE
25    GHOSTED  = DM_BOUNDARY_GHOSTED
26    MIRROR   = DM_BOUNDARY_MIRROR
27    PERIODIC = DM_BOUNDARY_PERIODIC
28    TWIST    = DM_BOUNDARY_TWIST
29
30
31class DMPolytopeType(object):
32    """The `DM` cell types."""
33    POINT              = DM_POLYTOPE_POINT
34    SEGMENT            = DM_POLYTOPE_SEGMENT
35    POINT_PRISM_TENSOR = DM_POLYTOPE_POINT_PRISM_TENSOR
36    TRIANGLE           = DM_POLYTOPE_TRIANGLE
37    QUADRILATERAL      = DM_POLYTOPE_QUADRILATERAL
38    SEG_PRISM_TENSOR   = DM_POLYTOPE_SEG_PRISM_TENSOR
39    TETRAHEDRON        = DM_POLYTOPE_TETRAHEDRON
40    HEXAHEDRON         = DM_POLYTOPE_HEXAHEDRON
41    TRI_PRISM          = DM_POLYTOPE_TRI_PRISM
42    TRI_PRISM_TENSOR   = DM_POLYTOPE_TRI_PRISM_TENSOR
43    QUAD_PRISM_TENSOR  = DM_POLYTOPE_QUAD_PRISM_TENSOR
44    PYRAMID            = DM_POLYTOPE_PYRAMID
45    FV_GHOST           = DM_POLYTOPE_FV_GHOST
46    INTERIOR_GHOST     = DM_POLYTOPE_INTERIOR_GHOST
47    UNKNOWN            = DM_POLYTOPE_UNKNOWN
48    UNKNOWN_CELL       = DM_POLYTOPE_UNKNOWN_CELL
49    UNKNOWN_FACE       = DM_POLYTOPE_UNKNOWN_FACE
50
51
52class DMReorderDefaultFlag(object):
53    """The `DM` reordering default flags."""
54    NOTSET = DM_REORDER_DEFAULT_NOTSET
55    FALSE  = DM_REORDER_DEFAULT_FALSE
56    TRUE   = DM_REORDER_DEFAULT_TRUE
57
58# --------------------------------------------------------------------
59
60
61cdef class DM(Object):
62    """An object describing a computational grid or mesh."""
63
64    Type         = DMType
65    BoundaryType = DMBoundaryType
66    PolytopeType = DMPolytopeType
67
68    ReorderDefaultFlag = DMReorderDefaultFlag
69
70    #
71
72    def __cinit__(self):
73        self.obj = <PetscObject*> &self.dm
74        self.dm  = NULL
75
76    def view(self, Viewer viewer=None) -> None:
77        """View the `DM`.
78
79        Collective.
80
81        Parameters
82        ----------
83        viewer
84            The DM viewer.
85
86        See Also
87        --------
88        petsc.DMView
89
90        """
91        cdef PetscViewer vwr = NULL
92        if viewer is not None: vwr = viewer.vwr
93        CHKERR(DMView(self.dm, vwr))
94
95    def load(self, Viewer viewer) -> Self:
96        """Return a `DM` stored in binary.
97
98        Collective.
99
100        Parameters
101        ----------
102        viewer
103            Viewer used to store the `DM`,
104            like `Viewer.Type.BINARY` or `Viewer.Type.HDF5`.
105
106        Notes
107        -----
108        When using `Viewer.Type.HDF5` format, one can save multiple `DMPlex` meshes
109        in a single HDF5 files. This in turn requires one to name the `DMPlex`
110        object with `Object.setName` before saving it with `DM.view` and before
111        loading it with `DM.load` for identification of the mesh object.
112
113        See Also
114        --------
115        DM.view, DM.load, Object.setName, petsc.DMLoad
116
117        """
118        CHKERR(DMLoad(self.dm, viewer.vwr))
119        return self
120
121    def destroy(self) -> Self:
122        """Destroy the object.
123
124        Collective.
125
126        See Also
127        --------
128        petsc.DMDestroy
129
130        """
131        CHKERR(DMDestroy(&self.dm))
132        return self
133
134    def create(self, comm: Comm | None = None) -> Self:
135        """Return an empty `DM`.
136
137        Collective.
138
139        Parameters
140        ----------
141        comm
142            MPI communicator, defaults to `Sys.getDefaultComm`.
143
144        See Also
145        --------
146        petsc.DMCreate
147
148        """
149        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
150        cdef PetscDM newdm = NULL
151        CHKERR(DMCreate(ccomm, &newdm))
152        CHKERR(PetscCLEAR(self.obj)); self.dm = newdm
153        return self
154
155    def clone(self) -> DM:
156        """Return the cloned `DM` .
157
158        Collective.
159
160        See Also
161        --------
162        petsc.DMClone
163
164        """
165        cdef DM dm = type(self)()
166        CHKERR(DMClone(self.dm, &dm.dm))
167        return dm
168
169    def setType(self, dm_type: DM.Type | str) -> None:
170        """Build a `DM`.
171
172        Collective.
173
174        Parameters
175        ----------
176        dm_type
177            The type of `DM`.
178
179        Notes
180        -----
181        `DM` types are available in `DM.Type` class.
182
183        See Also
184        --------
185        DM.Type, petsc.DMSetType
186
187        """
188        cdef PetscDMType cval = NULL
189        dm_type = str2bytes(dm_type, &cval)
190        CHKERR(DMSetType(self.dm, cval))
191
192    def getType(self) -> str:
193        """Return the `DM` type name.
194
195        Not collective.
196
197        See Also
198        --------
199        petsc.DMGetType
200
201        """
202        cdef PetscDMType cval = NULL
203        CHKERR(DMGetType(self.dm, &cval))
204        return bytes2str(cval)
205
206    def getDimension(self) -> int:
207        """Return the topological dimension of the `DM`.
208
209        Not collective.
210
211        See Also
212        --------
213        petsc.DMGetDimension
214
215        """
216        cdef PetscInt dim = 0
217        CHKERR(DMGetDimension(self.dm, &dim))
218        return toInt(dim)
219
220    def setDimension(self, dim: int) -> None:
221        """Set the topological dimension of the `DM`.
222
223        Collective.
224
225        Parameters
226        ----------
227        dim
228            Topological dimension.
229
230        See Also
231        --------
232        petsc.DMSetDimension
233
234        """
235        cdef PetscInt cdim = asInt(dim)
236        CHKERR(DMSetDimension(self.dm, cdim))
237
238    def getCoordinateDim(self) -> int:
239        """Return the dimension of embedding space for coordinates values.
240
241        Not collective.
242
243        See Also
244        --------
245        petsc.DMGetCoordinateDim
246
247        """
248        cdef PetscInt dim = 0
249        CHKERR(DMGetCoordinateDim(self.dm, &dim))
250        return toInt(dim)
251
252    def setCoordinateDim(self, dim: int) -> None:
253        """Set the dimension of embedding space for coordinates values.
254
255        Not collective.
256
257        Parameters
258        ----------
259        dim
260            The embedding dimension.
261
262        See Also
263        --------
264        petsc.DMSetCoordinateDim
265
266        """
267        cdef PetscInt cdim = asInt(dim)
268        CHKERR(DMSetCoordinateDim(self.dm, cdim))
269
270    def setOptionsPrefix(self, prefix: str | None) -> None:
271        """Set the prefix used for searching for options in the database.
272
273        Logically collective.
274
275        See Also
276        --------
277        petsc_options, getOptionsPrefix, petsc.DMSetOptionsPrefix
278
279        """
280        cdef const char *cval = NULL
281        prefix = str2bytes(prefix, &cval)
282        CHKERR(DMSetOptionsPrefix(self.dm, cval))
283
284    def getOptionsPrefix(self) -> str:
285        """Return the prefix used for searching for options in the database.
286
287        Not collective.
288
289        See Also
290        --------
291        petsc_options, setOptionsPrefix, petsc.DMGetOptionsPrefix
292
293        """
294        cdef const char *cval = NULL
295        CHKERR(DMGetOptionsPrefix(self.dm, &cval))
296        return bytes2str(cval)
297
298    def appendOptionsPrefix(self, prefix: str | None) -> None:
299        """Append to the prefix used for searching for options in the database.
300
301        Logically collective.
302
303        See Also
304        --------
305        petsc_options, setOptionsPrefix, petsc.DMAppendOptionsPrefix
306
307        """
308        cdef const char *cval = NULL
309        prefix = str2bytes(prefix, &cval)
310        CHKERR(DMAppendOptionsPrefix(self.dm, cval))
311
312    def setFromOptions(self) -> None:
313        """Configure the object from the options database.
314
315        Collective.
316
317        See Also
318        --------
319        petsc_options, petsc.DMSetFromOptions
320
321        """
322        CHKERR(DMSetFromOptions(self.dm))
323
324    def setUp(self) -> Self:
325        """Return the data structure.
326
327        Collective.
328
329        See Also
330        --------
331        petsc.DMSetUp
332
333        """
334        CHKERR(DMSetUp(self.dm))
335        return self
336
337    # --- application context ---
338
339    def setAppCtx(self, appctx: Any) -> None:
340        """Set the application context."""
341        self.set_attr('__appctx__', appctx)
342
343    def getAppCtx(self) -> Any:
344        """Return the application context."""
345        return self.get_attr('__appctx__')
346
347    #
348
349    def getUseNatural(self) -> bool:
350        """Get the flag for constructing a global-to-natural map.
351
352        Not collective.
353
354        Returns
355        -------
356        useNatural : bool
357            Whether a global-to-natural map is created.
358
359        See Also
360        --------
361        petsc.DMGetUseNatural
362
363        """
364        cdef PetscBool uN  = PETSC_FALSE
365        CHKERR(DMGetUseNatural(self.dm, &uN))
366        return toBool(uN)
367
368    def setUseNatural(self, useNatural : bool) -> None:
369        """Set the flag for constructing a global-to-natural map.
370
371        Not collective.
372
373        Parameters
374        ----------
375        useNatural : bool
376            Whether a global-to-natural map is created.
377
378        See Also
379        --------
380        petsc.DMSetUseNatural
381
382        """
383        cdef PetscBool uN  = useNatural
384        CHKERR(DMSetUseNatural(self.dm, uN))
385
386    def setBasicAdjacency(self, useCone: bool, useClosure: bool) -> None:
387        """Set the flags for determining variable influence.
388
389        Not collective.
390
391        Parameters
392        ----------
393        useCone : bool
394            Whether adjacency uses cone information.
395        useClosure : bool
396            Whether adjacency is computed using full closure information.
397
398        See Also
399        --------
400        petsc.DMSetBasicAdjacency
401
402        """
403        cdef PetscBool uC  = useCone
404        cdef PetscBool uCl = useClosure
405        CHKERR(DMSetBasicAdjacency(self.dm, uC, uCl))
406
407    def getBasicAdjacency(self) -> tuple[bool, bool]:
408        """Return the flags for determining variable influence.
409
410        Not collective.
411
412        Returns
413        -------
414        useCone : bool
415            Whether adjacency uses cone information.
416        useClosure : bool
417            Whether adjacency is computed using full closure information.
418
419        See Also
420        --------
421        petsc.DMGetBasicAdjacency
422
423        """
424        cdef PetscBool uC  = PETSC_FALSE
425        cdef PetscBool uCl = PETSC_FALSE
426        CHKERR(DMGetBasicAdjacency(self.dm, &uC, &uCl))
427        return toBool(uC), toBool(uCl)
428
429    def setFieldAdjacency(self, field: int, useCone: bool, useClosure: bool) -> None:
430        """Set the flags for determining variable influence.
431
432        Not collective.
433
434        Parameters
435        ----------
436        field : int
437            The field number.
438        useCone : bool
439            Whether adjacency uses cone information.
440        useClosure : bool
441            Whether adjacency is computed using full closure information.
442
443        See Also
444        --------
445        petsc.DMSetAdjacency
446
447        """
448        cdef PetscInt  f   = asInt(field)
449        cdef PetscBool uC  = useCone
450        cdef PetscBool uCl = useClosure
451        CHKERR(DMSetAdjacency(self.dm, f, uC, uCl))
452
453    def getFieldAdjacency(self, field: int) -> tuple[bool, bool]:
454        """Return the flags for determining variable influence.
455
456        Not collective.
457
458        Parameters
459        ----------
460        field
461            The field number.
462
463        Returns
464        -------
465        useCone : bool
466            Whether adjacency uses cone information.
467        useClosure : bool
468            Whether adjacency is computed using full closure information.
469
470        See Also
471        --------
472        petsc.DMGetAdjacency
473
474        """
475        cdef PetscInt  f   = asInt(field)
476        cdef PetscBool uC  = PETSC_FALSE
477        cdef PetscBool uCl = PETSC_FALSE
478        CHKERR(DMGetAdjacency(self.dm, f, &uC, &uCl))
479        return toBool(uC), toBool(uCl)
480
481    #
482
483    def createSubDM(self, fields: Sequence[int]) -> tuple[IS, DM]:
484        """Return `IS` and `DM` encapsulating a subproblem.
485
486        Not collective.
487
488        Returns
489        -------
490        iset : IS
491            The global indices for all the degrees of freedom.
492        subdm : DM
493            The `DM` for the subproblem.
494
495        See Also
496        --------
497        petsc.DMCreateSubDM
498
499        """
500        cdef IS iset = IS()
501        cdef DM subdm = DM()
502        cdef PetscInt *ifields = NULL
503        cdef PetscInt numFields = 0
504        fields = iarray_i(fields, &numFields, &ifields)
505        CHKERR(DMCreateSubDM(self.dm, numFields, ifields, &iset.iset, &subdm.dm))
506        return iset, subdm
507
508    #
509
510    def setAuxiliaryVec(self, Vec aux, label: DMLabel | None, value=0, part=0) -> None:
511        """Set an auxiliary vector for a specific region.
512
513        Not collective.
514
515        Parameters
516        ----------
517        aux
518            The auxiliary vector.
519        label
520            The name of the `DMLabel`.
521        value
522            Indicate the region.
523        part
524            The equation part, or 0 is unused.
525
526        See Also
527        --------
528        petsc.DMGetLabel, petsc.DMSetAuxiliaryVec
529
530        """
531        cdef PetscInt cvalue = asInt(value)
532        cdef PetscInt cpart = asInt(part)
533        cdef const char *cval = NULL
534        cdef PetscDMLabel clbl = NULL
535        label = str2bytes(label, &cval)
536        if cval == NULL: cval = b"" # XXX Should be fixed upstream
537        CHKERR(DMGetLabel(self.dm, cval, &clbl))
538        CHKERR(DMSetAuxiliaryVec(self.dm, clbl, cvalue, cpart, aux.vec))
539
540    def getAuxiliaryVec(self, label: str | None = None, value: int | None = 0, part: int | None = 0) -> Vec:
541        """Return an auxiliary vector for region.
542
543        Not collective.
544
545        Parameters
546        ----------
547        label
548            The name of the `DMLabel`.
549        value
550            Indicate the region.
551        part
552            The equation part, or 0 is unused.
553
554        See Also
555        --------
556        DM.getLabel, petsc.DMGetAuxiliaryVec
557
558        """
559        cdef PetscInt cvalue = asInt(value)
560        cdef PetscInt cpart = asInt(part)
561        cdef const char *cval = NULL
562        cdef PetscDMLabel clbl = NULL
563        cdef Vec aux = Vec()
564        label = str2bytes(label, &cval)
565        if cval == NULL: cval = b"" # XXX Should be fixed upstream
566        CHKERR(DMGetLabel(self.dm, cval, &clbl))
567        CHKERR(DMGetAuxiliaryVec(self.dm, clbl, cvalue, cpart, &aux.vec))
568        return aux
569
570    def setNumFields(self, numFields: int) -> None:
571        """Set the number of fields in the `DM`.
572
573        Logically collective.
574
575        See Also
576        --------
577        petsc.DMSetNumFields
578
579        """
580        cdef PetscInt cnum = asInt(numFields)
581        CHKERR(DMSetNumFields(self.dm, cnum))
582
583    def getNumFields(self) -> int:
584        """Return the number of fields in the `DM`.
585
586        Not collective.
587
588        See Also
589        --------
590        petsc.DMGetNumFields
591
592        """
593        cdef PetscInt cnum = 0
594        CHKERR(DMGetNumFields(self.dm, &cnum))
595        return toInt(cnum)
596
597    def setField(self, index: int, Object field, label: str | None = None) -> None:
598        """Set the discretization object for a given `DM` field.
599
600        Logically collective.
601
602        Parameters
603        ----------
604        index
605            The field number.
606        field
607            The discretization object.
608        label
609            The name of the label indicating the support of the field,
610            or `None` for the entire mesh.
611
612        See Also
613        --------
614        petsc.DMSetField
615
616        """
617        cdef PetscInt     cidx = asInt(index)
618        cdef PetscObject  cobj = field.obj[0]
619        cdef PetscDMLabel clbl = NULL
620        assert label is None
621        CHKERR(DMSetField(self.dm, cidx, clbl, cobj))
622
623    def getField(self, index: int) -> tuple[Object, None]:
624        """Return the discretization object for a given `DM` field.
625
626        Not collective.
627
628        Parameters
629        ----------
630        index
631            The field number.
632
633        See Also
634        --------
635        petsc.DMGetField
636
637        """
638        cdef PetscInt     cidx = asInt(index)
639        cdef PetscObject  cobj = NULL
640        cdef PetscDMLabel clbl = NULL
641        CHKERR(DMGetField(self.dm, cidx, &clbl, &cobj))
642        assert clbl == NULL
643        cdef Object field = subtype_Object(cobj)()
644        field.obj[0] = cobj
645        CHKERR(PetscINCREF(field.obj))
646        return (field, None) # TODO REVIEW
647
648    def addField(self, Object field, label: str | None = None) -> None:
649        """Add a field to a `DM` object.
650
651        Logically collective.
652
653        Parameters
654        ----------
655        field
656            The discretization object.
657        label
658            The name of the label indicating the support of the field,
659            or `None` for the entire mesh.
660
661        See Also
662        --------
663        petsc.DMAddField
664
665        """
666        cdef PetscObject  cobj = field.obj[0]
667        cdef PetscDMLabel clbl = NULL
668        assert label is None
669        CHKERR(DMAddField(self.dm, clbl, cobj))
670
671    def clearFields(self) -> None:
672        """Remove all fields from the `DM`.
673
674        Logically collective.
675
676        See Also
677        --------
678        petsc.DMClearFields
679
680        """
681        CHKERR(DMClearFields(self.dm))
682
683    def copyFields(self, DM dm, minDegree = None, maxDegree = None) -> None:
684        """Copy the discretizations of this `DM` into another `DM`.
685
686        Collective.
687
688        Parameters
689        ----------
690        dm
691            The `DM` that the fields are copied into.
692        minDegree
693            The minimum polynommial degree for the discretization,
694            or `None` for no limit
695        maxDegree
696            The maximum polynommial degree for the discretization,
697            or `None` for no limit
698
699        See Also
700        --------
701        petsc.DMCopyFields
702
703        """
704        cdef PetscInt mindeg = PETSC_DETERMINE
705        if minDegree is None:
706            pass
707        else:
708            mindeg = asInt(minDegree)
709        cdef PetscInt maxdeg = PETSC_DETERMINE
710        if maxDegree is None:
711            pass
712        else:
713            maxdeg = asInt(maxDegree)
714        CHKERR(DMCopyFields(self.dm, mindeg, maxdeg, dm.dm))
715
716    def createDS(self) -> None:
717        """Create discrete systems.
718
719        Collective.
720
721        See Also
722        --------
723        petsc.DMCreateDS
724
725        """
726        CHKERR(DMCreateDS(self.dm))
727
728    def clearDS(self) -> None:
729        """Remove all discrete systems from the `DM`.
730
731        Logically collective.
732
733        See Also
734        --------
735        petsc.DMClearDS
736
737        """
738        CHKERR(DMClearDS(self.dm))
739
740    def getDS(self) -> DS:
741        """Return default `DS`.
742
743        Not collective.
744
745        See Also
746        --------
747        petsc.DMGetDS
748
749        """
750        cdef DS ds = DS()
751        CHKERR(DMGetDS(self.dm, &ds.ds))
752        CHKERR(PetscINCREF(ds.obj))
753        return ds
754
755    def copyDS(self, DM dm, minDegree = None, maxDegree = None) -> None:
756        """Copy the discrete systems for this `DM` into another `DM`.
757
758        Collective.
759
760        Parameters
761        ----------
762        dm
763            The `DM` that the discrete fields are copied into.
764        minDegree
765            The minimum polynommial degree for the discretization,
766            or `None` for no limit
767        maxDegree
768            The maximum polynommial degree for the discretization,
769            or `None` for no limit
770
771        See Also
772        --------
773        petsc.DMCopyDS
774
775        """
776        cdef PetscInt mindeg = PETSC_DETERMINE
777        if minDegree is None:
778            pass
779        else:
780            mindeg = asInt(minDegree)
781        cdef PetscInt maxdeg = PETSC_DETERMINE
782        if maxDegree is None:
783            pass
784        else:
785            maxdeg = asInt(maxDegree)
786        CHKERR(DMCopyDS(self.dm, mindeg, maxdeg, dm.dm))
787
788    def copyDisc(self, DM dm) -> None:
789        """Copy fields and discrete systems of a `DM` into another `DM`.
790
791        Collective.
792
793        Parameters
794        ----------
795        dm
796            The `DM` that the fields and discrete systems are copied into.
797
798        See Also
799        --------
800        petsc.DMCopyDisc
801
802        """
803        CHKERR(DMCopyDisc(self.dm, dm.dm))
804
805    #
806
807    def getBlockSize(self) -> int:
808        """Return the inherent block size associated with a `DM`.
809
810        Not collective.
811
812        See Also
813        --------
814        petsc.DMGetBlockSize
815
816        """
817        cdef PetscInt bs = 1
818        CHKERR(DMGetBlockSize(self.dm, &bs))
819        return toInt(bs)
820
821    def setVecType(self, vec_type: Vec.Type | str) -> None:
822        """Set the type of vector.
823
824        Logically collective.
825
826        See Also
827        --------
828        Vec.Type, petsc.DMSetVecType
829
830        """
831        cdef PetscVecType vtype = NULL
832        vec_type = str2bytes(vec_type, &vtype)
833        CHKERR(DMSetVecType(self.dm, vtype))
834
835    def createGlobalVec(self) -> Vec:
836        """Return a global vector.
837
838        Collective.
839
840        See Also
841        --------
842        petsc.DMCreateGlobalVector
843
844        """
845        cdef Vec vg = Vec()
846        CHKERR(DMCreateGlobalVector(self.dm, &vg.vec))
847        return vg
848
849    def createLocalVec(self) -> Vec:
850        """Return a local vector.
851
852        Not collective.
853
854        See Also
855        --------
856        petsc.DMCreateLocalVector
857
858        """
859        cdef Vec vl = Vec()
860        CHKERR(DMCreateLocalVector(self.dm, &vl.vec))
861        return vl
862
863    def getGlobalVec(self, name : str | None = None) -> Vec:
864        """Return a global vector.
865
866        Collective.
867
868        Parameters
869        ----------
870        name
871            The optional name to retrieve a persistent vector.
872
873        Notes
874        -----
875        When done with the vector, it must be restored using `restoreGlobalVec`.
876
877        See Also
878        --------
879        restoreGlobalVec, petsc.DMGetGlobalVector, petsc.DMGetNamedGlobalVector
880
881        """
882        cdef Vec vg = Vec()
883        cdef const char *cname = NULL
884        str2bytes(name, &cname)
885        if cname != NULL:
886            CHKERR(DMGetNamedGlobalVector(self.dm, cname, &vg.vec))
887        else:
888            CHKERR(DMGetGlobalVector(self.dm, &vg.vec))
889        CHKERR(PetscINCREF(vg.obj))
890        return vg
891
892    def restoreGlobalVec(self, Vec vg, name : str | None = None) -> None:
893        """Restore a global vector obtained with `getGlobalVec`.
894
895        Logically collective.
896
897        Parameters
898        ----------
899        vg
900            The global vector.
901        name
902            The name used to retrieve the persistent vector, if any.
903
904        See Also
905        --------
906        getGlobalVec, petsc.DMRestoreGlobalVector, petsc.DMRestoreNamedGlobalVector
907
908        """
909        cdef const char *cname = NULL
910        str2bytes(name, &cname)
911        CHKERR(PetscDECREF(vg.obj))
912        if cname != NULL:
913            CHKERR(DMRestoreNamedGlobalVector(self.dm, cname, &vg.vec))
914        else:
915            CHKERR(DMRestoreGlobalVector(self.dm, &vg.vec))
916
917    def getLocalVec(self, name : str | None = None) -> Vec:
918        """Return a local vector.
919
920        Not collective.
921
922        Parameters
923        ----------
924        name
925            The optional name to retrieve a persistent vector.
926
927        Notes
928        -----
929        When done with the vector, it must be restored using `restoreLocalVec`.
930
931        See Also
932        --------
933        restoreLocalVec, petsc.DMGetLocalVector, petsc.DMGetNamedLocalVector
934
935        """
936        cdef Vec vl = Vec()
937        cdef const char *cname = NULL
938        str2bytes(name, &cname)
939        if cname != NULL:
940            CHKERR(DMGetNamedLocalVector(self.dm, cname, &vl.vec))
941        else:
942            CHKERR(DMGetLocalVector(self.dm, &vl.vec))
943        CHKERR(PetscINCREF(vl.obj))
944        return vl
945
946    def restoreLocalVec(self, Vec vl, name : str | None = None) -> None:
947        """Restore a local vector obtained with `getLocalVec`.
948
949        Not collective.
950
951        Parameters
952        ----------
953        vl
954            The local vector.
955        name
956            The name used to retrieve the persistent vector, if any.
957
958        See Also
959        --------
960        getLocalVec, petsc.DMRestoreLocalVector, petsc.DMRestoreNamedLocalVector
961
962        """
963        cdef const char *cname = NULL
964        str2bytes(name, &cname)
965        CHKERR(PetscDECREF(vl.obj))
966        if cname != NULL:
967            CHKERR(DMRestoreNamedLocalVector(self.dm, cname, &vl.vec))
968        else:
969            CHKERR(DMRestoreLocalVector(self.dm, &vl.vec))
970
971    def globalToLocal(self, Vec vg, Vec vl, addv: InsertModeSpec | None = None) -> None:
972        """Update local vectors from global vector.
973
974        Neighborwise collective.
975
976        Parameters
977        ----------
978        vg
979            The global vector.
980        vl
981            The local vector.
982        addv
983            Insertion mode.
984
985        See Also
986        --------
987        petsc.DMGlobalToLocalBegin, petsc.DMGlobalToLocalEnd
988
989        """
990        cdef PetscInsertMode im = insertmode(addv)
991        CHKERR(DMGlobalToLocalBegin(self.dm, vg.vec, im, vl.vec))
992        CHKERR(DMGlobalToLocalEnd  (self.dm, vg.vec, im, vl.vec))
993
994    def localToGlobal(self, Vec vl, Vec vg, addv: InsertModeSpec | None = None) -> None:
995        """Update global vectors from local vector.
996
997        Neighborwise collective.
998
999        Parameters
1000        ----------
1001        vl
1002            The local vector.
1003        vg
1004            The global vector.
1005        addv
1006            Insertion mode.
1007
1008        See Also
1009        --------
1010        petsc.DMLocalToGlobalBegin, petsc.DMLocalToGlobalEnd
1011
1012        """
1013        cdef PetscInsertMode im = insertmode(addv)
1014        CHKERR(DMLocalToGlobalBegin(self.dm, vl.vec, im, vg.vec))
1015        CHKERR(DMLocalToGlobalEnd(self.dm, vl.vec, im, vg.vec))
1016
1017    def localToLocal(self, Vec vl, Vec vlg, addv: InsertModeSpec | None = None) -> None:
1018        """Map the values from a local vector to another local vector.
1019
1020        Neighborwise collective.
1021
1022        Parameters
1023        ----------
1024        vl
1025            The local vector.
1026        vlg
1027            The global vector.
1028        addv
1029            Insertion mode.
1030
1031        See Also
1032        --------
1033        petsc.DMLocalToLocalBegin, petsc.DMLocalToLocalEnd
1034
1035        """
1036        cdef PetscInsertMode im = insertmode(addv)
1037        CHKERR(DMLocalToLocalBegin(self.dm, vl.vec, im, vlg.vec))
1038        CHKERR(DMLocalToLocalEnd  (self.dm, vl.vec, im, vlg.vec))
1039
1040    def getLGMap(self) -> LGMap:
1041        """Return local mapping to global mapping.
1042
1043        Collective.
1044
1045        See Also
1046        --------
1047        petsc.DMGetLocalToGlobalMapping
1048
1049        """
1050        cdef LGMap lgm = LGMap()
1051        CHKERR(DMGetLocalToGlobalMapping(self.dm, &lgm.lgm))
1052        CHKERR(PetscINCREF(lgm.obj))
1053        return lgm
1054
1055    #
1056
1057    def getCoarseDM(self) -> DM:
1058        """Return the coarse `DM`.
1059
1060        Collective.
1061
1062        See Also
1063        --------
1064        petsc.DMGetCoarseDM
1065
1066        """
1067        cdef DM cdm = type(self)()
1068        CHKERR(DMGetCoarseDM(self.dm, &cdm.dm))
1069        CHKERR(PetscINCREF(cdm.obj))
1070        return cdm
1071
1072    def setCoarseDM(self, DM dm) -> None:
1073        """Set the coarse `DM`.
1074
1075        Collective.
1076
1077        See Also
1078        --------
1079        petsc.DMSetCoarseDM
1080
1081        """
1082        CHKERR(DMSetCoarseDM(self.dm, dm.dm))
1083        return
1084
1085    def getCoordinateDM(self) -> DM:
1086        """Return the coordinate `DM`.
1087
1088        Collective.
1089
1090        See Also
1091        --------
1092        petsc.DMGetCoordinateDM
1093
1094        """
1095        cdef DM cdm = type(self)()
1096        CHKERR(DMGetCoordinateDM(self.dm, &cdm.dm))
1097        CHKERR(PetscINCREF(cdm.obj))
1098        return cdm
1099
1100    def getCoordinateSection(self) -> Section:
1101        """Return coordinate values layout over the mesh.
1102
1103        Collective.
1104
1105        See Also
1106        --------
1107        petsc.DMGetCoordinateSection
1108
1109        """
1110        cdef Section sec = Section()
1111        CHKERR(DMGetCoordinateSection(self.dm, &sec.sec))
1112        CHKERR(PetscINCREF(sec.obj))
1113        return sec
1114
1115    def setCoordinates(self, Vec c) -> None:
1116        """Set a global vector that holds the coordinates.
1117
1118        Collective.
1119
1120        Parameters
1121        ----------
1122        c
1123            Coordinate Vector.
1124
1125        See Also
1126        --------
1127        petsc.DMSetCoordinates
1128
1129        """
1130        CHKERR(DMSetCoordinates(self.dm, c.vec))
1131
1132    def getCoordinates(self) -> Vec:
1133        """Return a global vector with the coordinates associated.
1134
1135        Collective.
1136
1137        See Also
1138        --------
1139        petsc.DMGetCoordinates
1140
1141        """
1142        cdef Vec c = Vec()
1143        CHKERR(DMGetCoordinates(self.dm, &c.vec))
1144        CHKERR(PetscINCREF(c.obj))
1145        return c
1146
1147    def setCoordinatesLocal(self, Vec c) -> None:
1148        """Set a local vector with the ghost point holding the coordinates.
1149
1150        Not collective.
1151
1152        Parameters
1153        ----------
1154        c
1155            Coordinate Vector.
1156
1157        See Also
1158        --------
1159        petsc.DMSetCoordinatesLocal
1160
1161        """
1162        CHKERR(DMSetCoordinatesLocal(self.dm, c.vec))
1163
1164    def getCoordinatesLocal(self) -> Vec:
1165        """Return a local vector with the coordinates associated.
1166
1167        Collective the first time it is called.
1168
1169        See Also
1170        --------
1171        petsc.DMGetCoordinatesLocal
1172
1173        """
1174        cdef Vec c = Vec()
1175        CHKERR(DMGetCoordinatesLocal(self.dm, &c.vec))
1176        CHKERR(PetscINCREF(c.obj))
1177        return c
1178
1179    def setCellCoordinateDM(self, DM dm) -> None:
1180        """Set the cell coordinate `DM`.
1181
1182        Collective.
1183
1184        Parameters
1185        ----------
1186        dm
1187            The cell coordinate `DM`.
1188
1189        See Also
1190        --------
1191        petsc.DMSetCellCoordinateDM
1192
1193        """
1194        CHKERR(DMSetCellCoordinateDM(self.dm, dm.dm))
1195
1196    def getCellCoordinateDM(self) -> DM:
1197        """Return the cell coordinate `DM`.
1198
1199        Collective.
1200
1201        See Also
1202        --------
1203        petsc.DMGetCellCoordinateDM
1204
1205        """
1206        cdef DM cdm = type(self)()
1207        CHKERR(DMGetCellCoordinateDM(self.dm, &cdm.dm))
1208        CHKERR(PetscINCREF(cdm.obj))
1209        return cdm
1210
1211    def setCellCoordinateSection(self, dim: int, Section sec) -> None:
1212        """Set the cell coordinate layout over the `DM`.
1213
1214        Collective.
1215
1216        Parameters
1217        ----------
1218        dim
1219            The embedding dimension, or `DETERMINE`.
1220        sec
1221            The cell coordinate `Section`.
1222
1223        See Also
1224        --------
1225        petsc.DMSetCellCoordinateSection
1226
1227        """
1228        cdef PetscInt cdim = asInt(dim)
1229        CHKERR(DMSetCellCoordinateSection(self.dm, cdim, sec.sec))
1230
1231    def getCellCoordinateSection(self) -> Section:
1232        """Return the cell coordinate layout over the `DM`.
1233
1234        Collective.
1235
1236        See Also
1237        --------
1238        petsc.DMGetCellCoordinateSection
1239
1240        """
1241        cdef Section sec = Section()
1242        CHKERR(DMGetCellCoordinateSection(self.dm, &sec.sec))
1243        CHKERR(PetscINCREF(sec.obj))
1244        return sec
1245
1246    def setCellCoordinates(self, Vec c) -> None:
1247        """Set a global vector with the cellwise coordinates.
1248
1249        Collective.
1250
1251        Parameters
1252        ----------
1253        c
1254            The global cell coordinate vector.
1255
1256        See Also
1257        --------
1258        petsc.DMSetCellCoordinates
1259
1260        """
1261        CHKERR(DMSetCellCoordinates(self.dm, c.vec))
1262
1263    def getCellCoordinates(self) -> Vec:
1264        """Return a global vector with the cellwise coordinates.
1265
1266        Collective.
1267
1268        See Also
1269        --------
1270        petsc.DMGetCellCoordinates
1271
1272        """
1273        cdef Vec c = Vec()
1274        CHKERR(DMGetCellCoordinates(self.dm, &c.vec))
1275        CHKERR(PetscINCREF(c.obj))
1276        return c
1277
1278    def setCellCoordinatesLocal(self, Vec c) -> None:
1279        """Set a local vector with the cellwise coordinates.
1280
1281        Not collective.
1282
1283        Parameters
1284        ----------
1285        c
1286            The local cell coordinate vector.
1287
1288        See Also
1289        --------
1290        petsc.DMSetCellCoordinatesLocal
1291
1292        """
1293        CHKERR(DMSetCellCoordinatesLocal(self.dm, c.vec))
1294
1295    def getCellCoordinatesLocal(self) -> Vec:
1296        """Return a local vector with the cellwise coordinates.
1297
1298        Collective.
1299
1300        See Also
1301        --------
1302        petsc.DMGetCellCoordinatesLocal
1303
1304        """
1305        cdef Vec c = Vec()
1306        CHKERR(DMGetCellCoordinatesLocal(self.dm, &c.vec))
1307        CHKERR(PetscINCREF(c.obj))
1308        return c
1309
1310    def setCoordinateDisc(self, FE disc, localized: bool, project: bool) -> Self:
1311        """Project coordinates to a different space.
1312
1313        Collective.
1314
1315        Parameters
1316        ----------
1317        disc
1318            The new coordinates discretization.
1319        localized
1320            Set a localized (DG) coordinate space
1321        project
1322            Project coordinates to new discretization
1323
1324        See Also
1325        --------
1326        petsc.DMSetCoordinateDisc
1327
1328        """
1329        cdef PetscBool clocalized = localized
1330        cdef PetscBool pr = project
1331        CHKERR(DMSetCoordinateDisc(self.dm, disc.fe, clocalized, pr))
1332        return self
1333
1334    def getCoordinatesLocalized(self) -> bool:
1335        """Check if the coordinates have been localized for cells.
1336
1337        Not collective.
1338
1339        See Also
1340        --------
1341        petsc.DMGetCoordinatesLocalized
1342
1343        """
1344        cdef PetscBool flag = PETSC_FALSE
1345        CHKERR(DMGetCoordinatesLocalized(self.dm, &flag))
1346        return toBool(flag)
1347
1348    def getBoundingBox(self) -> tuple[tuple[float, float], ...]:
1349        """Return the dimension of embedding space for coordinates values.
1350
1351        Not collective.
1352
1353        See Also
1354        --------
1355        petsc.DMGetBoundingBox
1356
1357        """
1358        cdef PetscInt dim=0
1359        CHKERR(DMGetCoordinateDim(self.dm, &dim))
1360        cdef PetscReal gmin[3], gmax[3]
1361        CHKERR(DMGetBoundingBox(self.dm, gmin, gmax))
1362        return tuple([(toReal(gmin[i]), toReal(gmax[i]))
1363                      for i from 0 <= i < dim])
1364
1365    def getLocalBoundingBox(self) -> tuple[tuple[float, float], ...]:
1366        """Return the bounding box for the piece of the `DM`.
1367
1368        Not collective.
1369
1370        See Also
1371        --------
1372        petsc.DMGetLocalBoundingBox
1373
1374        """
1375        cdef PetscInt dim=0
1376        CHKERR(DMGetCoordinateDim(self.dm, &dim))
1377        cdef PetscReal lmin[3], lmax[3]
1378        CHKERR(DMGetLocalBoundingBox(self.dm, lmin, lmax))
1379        return tuple([(toReal(lmin[i]), toReal(lmax[i]))
1380                      for i from 0 <= i < dim])
1381
1382    def localizeCoordinates(self) -> None:
1383        """Create local coordinates for cells having periodic faces.
1384
1385        Collective.
1386
1387        Notes
1388        -----
1389        Used if the mesh is periodic.
1390
1391        See Also
1392        --------
1393        petsc.DMLocalizeCoordinates
1394
1395        """
1396        CHKERR(DMLocalizeCoordinates(self.dm))
1397    #
1398
1399    def getPeriodicity(self) -> tuple[ArrayReal, ArrayReal, ArrayReal]:
1400        """Set the description of mesh periodicity.
1401
1402        Not collective.
1403
1404        Parameters
1405        ----------
1406        maxCell
1407            The longest allowable cell dimension in each direction.
1408        Lstart
1409            The start of each coordinate direction, usually [0, 0, 0]
1410        L
1411            The periodic length of each coordinate direction, or -1 for non-periodic
1412
1413        See Also
1414        --------
1415        petsc.DMSetPeriodicity
1416
1417        """
1418        cdef PetscInt dim = 0
1419        CHKERR(DMGetDimension(self.dm, &dim))
1420        cdef PetscReal *MAXCELL = NULL
1421        cdef ndarray mcell = oarray_r(empty_r(dim), NULL, &MAXCELL)
1422        cdef PetscReal *LSTART = NULL
1423        cdef ndarray lstart = oarray_r(empty_r(dim), NULL, &LSTART)
1424        cdef PetscReal *LPY = NULL
1425        cdef ndarray lpy = oarray_r(empty_r(dim), NULL, &LPY)
1426        cdef const PetscReal *maxCell = NULL
1427        cdef const PetscReal *Lstart = NULL
1428        cdef const PetscReal *L = NULL
1429        CHKERR(DMGetPeriodicity(self.dm, &maxCell, &Lstart, &L))
1430        CHKERR(PetscMemcpy(MAXCELL, maxCell, <size_t>dim*sizeof(PetscReal)))
1431        CHKERR(PetscMemcpy(LSTART, Lstart, <size_t>dim*sizeof(PetscReal)))
1432        CHKERR(PetscMemcpy(LPY, L, <size_t>dim*sizeof(PetscReal)))
1433        return (mcell, lstart, lpy)
1434
1435    def setPeriodicity(self, maxCell: Sequence[float], Lstart: Sequence[float], L: Sequence[float]) -> None:
1436        """Set the description of mesh periodicity.
1437
1438        Logically collective.
1439
1440        Parameters
1441        ----------
1442        maxCell
1443            The longest allowable cell dimension in each direction.
1444        Lstart
1445            The start of each coordinate direction, usually [0, 0, 0]
1446        L
1447            The periodic length of each coordinate direction, or -1 for non-periodic
1448
1449        See Also
1450        --------
1451        petsc.DMSetPeriodicity
1452
1453        """
1454        cdef PetscReal *MAXCELL = NULL
1455        cdef ndarray unuseda = oarray_r(maxCell, NULL, &MAXCELL)
1456        cdef PetscReal *LSTART = NULL
1457        cdef ndarray unusedb = oarray_r(Lstart, NULL, &LSTART)
1458        cdef PetscReal *LPY = NULL
1459        cdef ndarray unusedc = oarray_r(L, NULL, &LPY)
1460        CHKERR(DMSetPeriodicity(self.dm, MAXCELL, LSTART, LPY))
1461
1462    def setMatType(self, mat_type: Mat.Type | str) -> None:
1463        """Set matrix type to be used by `DM.createMat`.
1464
1465        Logically collective.
1466
1467        Parameters
1468        ----------
1469        mat_type
1470            The matrix type.
1471
1472        Notes
1473        -----
1474        The option ``-dm_mat_type`` is used to set the matrix type.
1475
1476        See Also
1477        --------
1478        petsc_options, petsc.DMSetMatType
1479
1480        """
1481        cdef PetscMatType mtype = NULL
1482        mat_type = str2bytes(mat_type, &mtype)
1483        CHKERR(DMSetMatType(self.dm, mtype))
1484
1485    def createMat(self) -> Mat:
1486        """Return an empty matrix.
1487
1488        Collective.
1489
1490        See Also
1491        --------
1492        petsc.DMCreateMatrix
1493
1494        """
1495        cdef Mat mat = Mat()
1496        CHKERR(DMCreateMatrix(self.dm, &mat.mat))
1497        return mat
1498
1499    def createMassMatrix(self, DM dmf) -> Mat:
1500        """Return the mass matrix between this `DM` and the given `DM`.
1501
1502        Collective.
1503
1504        Parameters
1505        ----------
1506        dmf
1507            The second `DM`.
1508
1509        See Also
1510        --------
1511        petsc.DMCreateMassMatrix
1512
1513        """
1514        cdef Mat mat = Mat()
1515        CHKERR(DMCreateMassMatrix(self.dm, dmf.dm, &mat.mat))
1516        return mat
1517
1518    def createInterpolation(self, DM dm) -> tuple[Mat, Vec]:
1519        """Return the interpolation matrix to a finer `DM`.
1520
1521        Collective.
1522
1523        Parameters
1524        ----------
1525        dm
1526            The second, finer `DM`.
1527
1528        See Also
1529        --------
1530        petsc.DMCreateInterpolation
1531
1532        """
1533        cdef Mat A = Mat()
1534        cdef Vec scale = Vec()
1535        CHKERR(DMCreateInterpolation(self.dm, dm.dm,
1536                                     &A.mat, &scale.vec))
1537        return (A, scale)
1538
1539    def createInjection(self, DM dm) -> Mat:
1540        """Return the injection matrix into a finer `DM`.
1541
1542        Collective.
1543
1544        Parameters
1545        ----------
1546        dm
1547            The second, finer `DM` object.
1548
1549        See Also
1550        --------
1551        petsc.DMCreateInjection
1552
1553        """
1554        cdef Mat inject = Mat()
1555        CHKERR(DMCreateInjection(self.dm, dm.dm, &inject.mat))
1556        return inject
1557
1558    def createRestriction(self, DM dm) -> Mat:
1559        """Return the restriction matrix between this `DM` and the given `DM`.
1560
1561        Collective.
1562
1563        Parameters
1564        ----------
1565        dm
1566            The second, finer `DM` object.
1567
1568        See Also
1569        --------
1570        petsc.DMCreateRestriction
1571
1572        """
1573        cdef Mat mat = Mat()
1574        CHKERR(DMCreateRestriction(self.dm, dm.dm, &mat.mat))
1575        return mat
1576
1577    def convert(self, dm_type: DM.Type | str) -> DM:
1578        """Return a `DM` converted to another `DM`.
1579
1580        Collective.
1581
1582        Parameters
1583        ----------
1584        dm_type
1585            The new `DM.Type`, use ``“same”`` for the same type.
1586
1587        See Also
1588        --------
1589        DM.Type, petsc.DMConvert
1590
1591        """
1592        cdef PetscDMType cval = NULL
1593        dm_type = str2bytes(dm_type, &cval)
1594        cdef PetscDM newdm = NULL
1595        CHKERR(DMConvert(self.dm, cval, &newdm))
1596        cdef DM dm = <DM>subtype_DM(newdm)()
1597        dm.dm = newdm
1598        return dm
1599
1600    def refine(self, comm: Comm | None = None) -> DM:
1601        """Return a refined `DM` object.
1602
1603        Collective.
1604
1605        Parameters
1606        ----------
1607        comm
1608            MPI communicator, defaults to `Sys.getDefaultComm`.
1609
1610        See Also
1611        --------
1612        petsc.DMRefine
1613
1614        """
1615        cdef MPI_Comm dmcomm = MPI_COMM_NULL
1616        CHKERR(PetscObjectGetComm(<PetscObject>self.dm, &dmcomm))
1617        dmcomm = def_Comm(comm, dmcomm)
1618        cdef PetscDM newdm = NULL
1619        CHKERR(DMRefine(self.dm, dmcomm, &newdm))
1620        cdef DM dm = subtype_DM(newdm)()
1621        dm.dm = newdm
1622        return dm
1623
1624    def coarsen(self, comm: Comm | None = None) -> DM:
1625        """Return a coarsened `DM` object.
1626
1627        Collective.
1628
1629        Parameters
1630        ----------
1631        comm
1632            MPI communicator, defaults to `Sys.getDefaultComm`.
1633
1634        See Also
1635        --------
1636        petsc.DMCoarsen
1637
1638        """
1639        cdef MPI_Comm dmcomm = MPI_COMM_NULL
1640        CHKERR(PetscObjectGetComm(<PetscObject>self.dm, &dmcomm))
1641        dmcomm = def_Comm(comm, dmcomm)
1642        cdef PetscDM newdm = NULL
1643        CHKERR(DMCoarsen(self.dm, dmcomm, &newdm))
1644        cdef DM dm = subtype_DM(newdm)()
1645        dm.dm = newdm
1646        return dm
1647
1648    def refineHierarchy(self, nlevels: int) -> list:
1649        """Refine this `DM` and return the refined `DM` hierarchy.
1650
1651        Collective.
1652
1653        Parameters
1654        ----------
1655        nlevels
1656            The number of levels of refinement.
1657
1658        See Also
1659        --------
1660        petsc.DMRefineHierarchy
1661
1662        """
1663        cdef PetscInt i, n = asInt(nlevels)
1664        cdef PetscDM *newdmf = NULL
1665        cdef object unused = oarray_p(empty_p(<PetscInt>n), NULL, <void**>&newdmf)
1666        CHKERR(DMRefineHierarchy(self.dm, n, newdmf))
1667        cdef DM dmf = None
1668        cdef list hierarchy = []
1669        for i from 0 <= i < n:
1670            dmf = subtype_DM(newdmf[i])()
1671            dmf.dm = newdmf[i]
1672            hierarchy.append(dmf)
1673        return hierarchy
1674
1675    def coarsenHierarchy(self, nlevels: int) -> list:
1676        """Coarsen this `DM` and return the coarsened `DM` hierarchy.
1677
1678        Collective.
1679
1680        Parameters
1681        ----------
1682        nlevels
1683            The number of levels of coarsening.
1684
1685        See Also
1686        --------
1687        petsc.DMCoarsenHierarchy
1688
1689        """
1690        cdef PetscInt i, n = asInt(nlevels)
1691        cdef PetscDM *newdmc = NULL
1692        cdef object unused = oarray_p(empty_p(<PetscInt>n), NULL, <void**>&newdmc)
1693        CHKERR(DMCoarsenHierarchy(self.dm, n, newdmc))
1694        cdef DM dmc = None
1695        cdef list hierarchy = []
1696        for i from 0 <= i < n:
1697            dmc = subtype_DM(newdmc[i])()
1698            dmc.dm = newdmc[i]
1699            hierarchy.append(dmc)
1700        return hierarchy
1701
1702    def getRefineLevel(self) -> int:
1703        """Return the refinement level.
1704
1705        Not collective.
1706
1707        See Also
1708        --------
1709        petsc.DMGetRefineLevel
1710
1711        """
1712        cdef PetscInt n = 0
1713        CHKERR(DMGetRefineLevel(self.dm, &n))
1714        return toInt(n)
1715
1716    def setRefineLevel(self, level: int) -> None:
1717        """Set the number of refinements.
1718
1719        Not collective.
1720
1721        Parameters
1722        ----------
1723        nlevels
1724            The number of refinement.
1725
1726        See Also
1727        --------
1728        petsc.DMSetRefineLevel
1729
1730        """
1731        cdef PetscInt clevel = asInt(level)
1732        CHKERR(DMSetRefineLevel(self.dm, clevel))
1733
1734    def getCoarsenLevel(self) -> int:
1735        """Return the number of coarsenings.
1736
1737        Not collective.
1738
1739        See Also
1740        --------
1741        petsc.DMGetCoarsenLevel
1742
1743        """
1744        cdef PetscInt n = 0
1745        CHKERR(DMGetCoarsenLevel(self.dm, &n))
1746        return toInt(n)
1747
1748    #
1749
1750    def adaptLabel(self, label: str) -> DM:
1751        """Adapt a `DM` based on a `DMLabel`.
1752
1753        Collective.
1754
1755        Parameters
1756        ----------
1757        label
1758            The name of the `DMLabel`.
1759
1760        See Also
1761        --------
1762        petsc.DMAdaptLabel
1763
1764        """
1765        cdef const char *cval = NULL
1766        cdef PetscDMLabel clbl = NULL
1767        label = str2bytes(label, &cval)
1768        CHKERR(DMGetLabel(self.dm, cval, &clbl))
1769        cdef DM newdm = DMPlex()
1770        CHKERR(DMAdaptLabel(self.dm, clbl, &newdm.dm))
1771        return newdm
1772
1773    def adaptMetric(
1774        self,
1775        Vec metric,
1776        bdLabel: str | None = None,
1777        rgLabel: str | None = None) -> DM:
1778        """Return a mesh adapted to the specified metric field.
1779
1780        Collective.
1781
1782        Parameters
1783        ----------
1784        metric
1785            The metric to which the mesh is adapted, defined vertex-wise.
1786        bdLabel
1787            Label for boundary tags.
1788        rgLabel
1789            Label for cell tag.
1790
1791        See Also
1792        --------
1793        petsc.DMAdaptMetric
1794
1795        """
1796        cdef const char *cval = NULL
1797        cdef PetscDMLabel cbdlbl = NULL
1798        cdef PetscDMLabel crglbl = NULL
1799        bdLabel = str2bytes(bdLabel, &cval)
1800        if cval == NULL: cval = b"" # XXX Should be fixed upstream
1801        CHKERR(DMGetLabel(self.dm, cval, &cbdlbl))
1802        rgLabel = str2bytes(rgLabel, &cval)
1803        if cval == NULL: cval = b"" # XXX Should be fixed upstream
1804        CHKERR(DMGetLabel(self.dm, cval, &crglbl))
1805        cdef DM newdm = DMPlex()
1806        CHKERR(DMAdaptMetric(self.dm, metric.vec, cbdlbl, crglbl, &newdm.dm))
1807        return newdm
1808
1809    def getLabel(self, name: str) -> DMLabel:
1810        """Return the label of a given name.
1811
1812        Not collective.
1813
1814        See Also
1815        --------
1816        petsc.DMGetLabel
1817
1818        """
1819        cdef const char *cname = NULL
1820        cdef DMLabel dmlabel = DMLabel()
1821        name = str2bytes(name, &cname)
1822        CHKERR(DMGetLabel(self.dm, cname, &dmlabel.dmlabel))
1823        CHKERR(PetscINCREF(dmlabel.obj))
1824        return dmlabel
1825
1826    #
1827
1828    def setLocalSection(self, Section sec) -> None:
1829        """Set the `Section` encoding the local data layout for the `DM`.
1830
1831        Collective.
1832
1833        See Also
1834        --------
1835        petsc.DMSetLocalSection
1836
1837        """
1838        CHKERR(DMSetLocalSection(self.dm, sec.sec))
1839
1840    def getLocalSection(self) -> Section:
1841        """Return the `Section` encoding the local data layout for the `DM`.
1842
1843        Not collective.
1844
1845        See Also
1846        --------
1847        petsc.DMGetGlobalSection
1848
1849        """
1850        cdef Section sec = Section()
1851        CHKERR(DMGetLocalSection(self.dm, &sec.sec))
1852        CHKERR(PetscINCREF(sec.obj))
1853        return sec
1854
1855    def setGlobalSection(self, Section sec) -> None:
1856        """Set the `Section` encoding the global data layout for the `DM`.
1857
1858        Collective.
1859
1860        See Also
1861        --------
1862        petsc.DMSetGlobalSection
1863
1864        """
1865        CHKERR(DMSetGlobalSection(self.dm, sec.sec))
1866
1867    def getGlobalSection(self) -> Section:
1868        """Return the `Section` encoding the global data layout for the `DM`.
1869
1870        Collective the first time it is called.
1871
1872        See Also
1873        --------
1874        petsc.DMGetGlobalSection
1875
1876        """
1877        cdef Section sec = Section()
1878        CHKERR(DMGetGlobalSection(self.dm, &sec.sec))
1879        CHKERR(PetscINCREF(sec.obj))
1880        return sec
1881
1882    setSection = setLocalSection
1883    getSection = getLocalSection
1884    setDefaultSection = setLocalSection
1885    getDefaultSection = getLocalSection
1886    setDefaultLocalSection = setLocalSection
1887    getDefaultLocalSection = getLocalSection
1888    setDefaultGlobalSection = setGlobalSection
1889    getDefaultGlobalSection = getGlobalSection
1890
1891    def createSectionSF(self, Section localsec, Section globalsec) -> None:
1892        """Create the `SF` encoding the parallel DOF overlap for the `DM`.
1893
1894        Collective.
1895
1896        Parameters
1897        ----------
1898        localsec
1899            Describe the local data layout.
1900        globalsec
1901            Describe the global data layout.
1902
1903        Notes
1904        -----
1905        Encoding based on the `Section` describing the data layout.
1906
1907        See Also
1908        --------
1909        DM.getSectionSF, petsc.DMCreateSectionSF
1910
1911        """
1912        CHKERR(DMCreateSectionSF(self.dm, localsec.sec, globalsec.sec))
1913
1914    def getSectionSF(self) -> SF:
1915        """Return the `Section` encoding the parallel DOF overlap.
1916
1917        Collective the first time it is called.
1918
1919        See Also
1920        --------
1921        petsc.DMGetSectionSF
1922
1923        """
1924        cdef SF sf = SF()
1925        CHKERR(DMGetSectionSF(self.dm, &sf.sf))
1926        CHKERR(PetscINCREF(sf.obj))
1927        return sf
1928
1929    def setSectionSF(self, SF sf) -> None:
1930        """Set the `Section` encoding the parallel DOF overlap for the `DM`.
1931
1932        Logically collective.
1933
1934        See Also
1935        --------
1936        petsc.DMSetSectionSF
1937
1938        """
1939        CHKERR(DMSetSectionSF(self.dm, sf.sf))
1940
1941    createDefaultSF = createSectionSF
1942    getDefaultSF = getSectionSF
1943    setDefaultSF = setSectionSF
1944
1945    def getPointSF(self) -> SF:
1946        """Return the `SF` encoding the parallel DOF overlap for the `DM`.
1947
1948        Not collective.
1949
1950        See Also
1951        --------
1952        petsc.DMGetPointSF
1953
1954        """
1955        cdef SF sf = SF()
1956        CHKERR(DMGetPointSF(self.dm, &sf.sf))
1957        CHKERR(PetscINCREF(sf.obj))
1958        return sf
1959
1960    def setPointSF(self, SF sf) -> None:
1961        """Set the `SF` encoding the parallel DOF overlap for the `DM`.
1962
1963        Logically collective.
1964
1965        See Also
1966        --------
1967        petsc.DMSetPointSF
1968
1969        """
1970        CHKERR(DMSetPointSF(self.dm, sf.sf))
1971
1972    def getNumLabels(self) -> int:
1973        """Return the number of labels defined by on the `DM`.
1974
1975        Not collective.
1976
1977        See Also
1978        --------
1979        petsc.DMGetNumLabels
1980
1981        """
1982        cdef PetscInt nLabels = 0
1983        CHKERR(DMGetNumLabels(self.dm, &nLabels))
1984        return toInt(nLabels)
1985
1986    def getLabelName(self, index: int) -> str:
1987        """Return the name of nth label.
1988
1989        Not collective.
1990
1991        Parameters
1992        ----------
1993        index
1994            The label number.
1995
1996        See Also
1997        --------
1998        petsc.DMGetLabelName
1999
2000        """
2001        cdef PetscInt cindex = asInt(index)
2002        cdef const char *cname = NULL
2003        CHKERR(DMGetLabelName(self.dm, cindex, &cname))
2004        return bytes2str(cname)
2005
2006    def hasLabel(self, name: str) -> bool:
2007        """Determine whether the `DM` has a label.
2008
2009        Not collective.
2010
2011        Parameters
2012        ----------
2013        name
2014            The label name.
2015
2016        See Also
2017        --------
2018        petsc.DMHasLabel
2019
2020        """
2021        cdef PetscBool flag = PETSC_FALSE
2022        cdef const char *cname = NULL
2023        name = str2bytes(name, &cname)
2024        CHKERR(DMHasLabel(self.dm, cname, &flag))
2025        return toBool(flag)
2026
2027    def createLabel(self, name: str) -> None:
2028        """Create a label of the given name if it does not already exist.
2029
2030        Not collective.
2031
2032        Parameters
2033        ----------
2034        name
2035            The label name.
2036
2037        See Also
2038        --------
2039        petsc.DMCreateLabel
2040
2041        """
2042        cdef const char *cname = NULL
2043        name = str2bytes(name, &cname)
2044        CHKERR(DMCreateLabel(self.dm, cname))
2045
2046    def removeLabel(self, name: str) -> None:
2047        """Remove and destroy the label by name.
2048
2049        Not collective.
2050
2051        Parameters
2052        ----------
2053        name
2054            The label name.
2055
2056        See Also
2057        --------
2058        petsc.DMRemoveLabel
2059
2060        """
2061        cdef const char *cname = NULL
2062        cdef PetscDMLabel clbl = NULL
2063        name = str2bytes(name, &cname)
2064        CHKERR(DMRemoveLabel(self.dm, cname, &clbl))
2065        # TODO: Once DMLabel is wrapped, this should return the label, like the C function.
2066        CHKERR(DMLabelDestroy(&clbl))
2067
2068    def getLabelValue(self, name: str, point: int) -> int:
2069        """Return the value in `DMLabel` for the given point.
2070
2071        Not collective.
2072
2073        Parameters
2074        ----------
2075        name
2076            The label name.
2077        point
2078            The mesh point
2079
2080        See Also
2081        --------
2082        petsc.DMGetLabelValue
2083
2084        """
2085        cdef PetscInt cpoint = asInt(point), value = 0
2086        cdef const char *cname = NULL
2087        name = str2bytes(name, &cname)
2088        CHKERR(DMGetLabelValue(self.dm, cname, cpoint, &value))
2089        return toInt(value)
2090
2091    def setLabelValue(self, name: str, point: int, value: int) -> None:
2092        """Set a point to a `DMLabel` with a given value.
2093
2094        Not collective.
2095
2096        Parameters
2097        ----------
2098        name
2099            The label name.
2100        point
2101            The mesh point.
2102        value
2103            The label value for the point.
2104
2105        See Also
2106        --------
2107        petsc.DMSetLabelValue
2108
2109        """
2110        cdef PetscInt cpoint = asInt(point), cvalue = asInt(value)
2111        cdef const char *cname = NULL
2112        name = str2bytes(name, &cname)
2113        CHKERR(DMSetLabelValue(self.dm, cname, cpoint, cvalue))
2114
2115    def clearLabelValue(self, name: str, point: int, value: int) -> None:
2116        """Remove a point from a `DMLabel` with given value.
2117
2118        Not collective.
2119
2120        Parameters
2121        ----------
2122        name
2123            The label name.
2124        point
2125            The mesh point.
2126        value
2127            The label value for the point.
2128
2129        See Also
2130        --------
2131        petsc.DMClearLabelValue
2132
2133        """
2134        cdef PetscInt cpoint = asInt(point), cvalue = asInt(value)
2135        cdef const char *cname = NULL
2136        name = str2bytes(name, &cname)
2137        CHKERR(DMClearLabelValue(self.dm, cname, cpoint, cvalue))
2138
2139    def getLabelSize(self, name: str) -> int:
2140        """Return the number of values that the `DMLabel` takes.
2141
2142        Not collective.
2143
2144        Parameters
2145        ----------
2146        name
2147            The label name.
2148
2149        See Also
2150        --------
2151        petsc.DMLabelGetNumValues, petsc.DMGetLabelSize
2152
2153        """
2154        cdef PetscInt size = 0
2155        cdef const char *cname = NULL
2156        name = str2bytes(name, &cname)
2157        CHKERR(DMGetLabelSize(self.dm, cname, &size))
2158        return toInt(size)
2159
2160    def getLabelIdIS(self, name: str) -> IS:
2161        """Return an `IS` of all values that the `DMLabel` takes.
2162
2163        Not collective.
2164
2165        Parameters
2166        ----------
2167        name
2168            The label name.
2169
2170        See Also
2171        --------
2172        petsc.DMLabelGetValueIS, petsc.DMGetLabelIdIS
2173
2174        """
2175        cdef const char *cname = NULL
2176        name = str2bytes(name, &cname)
2177        cdef IS lis = IS()
2178        CHKERR(DMGetLabelIdIS(self.dm, cname, &lis.iset))
2179        return lis
2180
2181    def getStratumSize(self, name: str, value: int) -> int:
2182        """Return the number of points in a label stratum.
2183
2184        Not collective.
2185
2186        Parameters
2187        ----------
2188        name
2189            The label name.
2190        value
2191            The stratum value.
2192
2193        See Also
2194        --------
2195        petsc.DMGetStratumSize
2196
2197        """
2198        cdef PetscInt size = 0
2199        cdef PetscInt cvalue = asInt(value)
2200        cdef const char *cname = NULL
2201        name = str2bytes(name, &cname)
2202        CHKERR(DMGetStratumSize(self.dm, cname, cvalue, &size))
2203        return toInt(size)
2204
2205    def getStratumIS(self, name: str, value: int) -> IS:
2206        """Return the points in a label stratum.
2207
2208        Not collective.
2209
2210        Parameters
2211        ----------
2212        name
2213            The label name.
2214        value
2215            The stratum value.
2216
2217        See Also
2218        --------
2219        petsc.DMGetStratumIS
2220
2221        """
2222        cdef PetscInt cvalue = asInt(value)
2223        cdef const char *cname = NULL
2224        name = str2bytes(name, &cname)
2225        cdef IS sis = IS()
2226        CHKERR(DMGetStratumIS(self.dm, cname, cvalue, &sis.iset))
2227        return sis
2228
2229    def clearLabelStratum(self, name: str, value: int) -> None:
2230        """Remove all points from a stratum.
2231
2232        Not collective.
2233
2234        Parameters
2235        ----------
2236        name
2237            The label name.
2238        value
2239            The stratum value.
2240
2241        See Also
2242        --------
2243        petsc.DMClearLabelStratum
2244
2245        """
2246        cdef PetscInt cvalue = asInt(value)
2247        cdef const char *cname = NULL
2248        name = str2bytes(name, &cname)
2249        CHKERR(DMClearLabelStratum(self.dm, cname, cvalue))
2250
2251    def setLabelOutput(self, name: str, output: bool) -> None:
2252        """Set if a given label should be saved to a view.
2253
2254        Not collective.
2255
2256        Parameters
2257        ----------
2258        name
2259            The label name.
2260        output
2261            If `True`, the label is saved to the viewer.
2262
2263        See Also
2264        --------
2265        petsc.DMSetLabelOutput
2266
2267        """
2268        cdef const char *cname = NULL
2269        name = str2bytes(name, &cname)
2270        cdef PetscBool coutput = output
2271        CHKERR(DMSetLabelOutput(self.dm, cname, coutput))
2272
2273    def getLabelOutput(self, name: str) -> bool:
2274        """Return the output flag for a given label.
2275
2276        Not collective.
2277
2278        Parameters
2279        ----------
2280        name
2281            The label name.
2282
2283        See Also
2284        --------
2285        petsc.DMGetLabelOutput
2286
2287        """
2288        cdef const char *cname = NULL
2289        name = str2bytes(name, &cname)
2290        cdef PetscBool coutput = PETSC_FALSE
2291        CHKERR(DMGetLabelOutput(self.dm, cname, &coutput))
2292        return coutput
2293
2294    # backward compatibility
2295    createGlobalVector = createGlobalVec
2296    createLocalVector = createLocalVec
2297    getMatrix = createMatrix = createMat
2298
2299    def setKSPComputeOperators(
2300        self, operators,
2301        args: tuple[Any, ...] | None = None,
2302        kargs: dict[str, Any] | None = None) -> None:
2303        """Matrix associated with the linear system.
2304
2305        Collective.
2306
2307        Parameters
2308        ----------
2309        operator
2310            Callback function to compute the operators.
2311        args
2312            Positional arguments for the callback.
2313        kargs
2314            Keyword arguments for the callback.
2315
2316        See Also
2317        --------
2318        petsc.DMKSPSetComputeOperators
2319
2320        """
2321        if args  is None: args  = ()
2322        if kargs is None: kargs = {}
2323        context = (operators, args, kargs)
2324        self.set_attr('__operators__', context)
2325        CHKERR(DMKSPSetComputeOperators(self.dm, KSP_ComputeOps, <void*>context))
2326
2327    def createFieldDecomposition(self) -> tuple[list, list, list]:
2328        """Return field splitting information.
2329
2330        Collective.
2331
2332        See Also
2333        --------
2334        petsc.DMCreateFieldDecomposition
2335
2336        """
2337        cdef PetscInt clen = 0
2338        cdef PetscIS *cis = NULL
2339        cdef PetscDM *cdm = NULL
2340        cdef char** cnamelist = NULL
2341
2342        CHKERR(DMCreateFieldDecomposition(self.dm, &clen, &cnamelist, &cis, &cdm))
2343
2344        cdef list isets = [ref_IS(cis[i]) for i from 0 <= i < clen]
2345        cdef list dms   = []
2346        cdef list names = []
2347        cdef DM dm = None
2348
2349        for i from 0 <= i < clen:
2350            if cdm != NULL:
2351                dm = subtype_DM(cdm[i])()
2352                dm.dm = cdm[i]
2353                CHKERR(PetscINCREF(dm.obj))
2354                dms.append(dm)
2355                CHKERR(DMDestroy(&cdm[i]))
2356            else:
2357                dms.append(None)
2358
2359            if cnamelist != NULL:
2360                name = bytes2str(cnamelist[i])
2361                names.append(name)
2362                CHKERR(PetscFree(cnamelist[i]))
2363            else:
2364                names.append(None)
2365
2366            CHKERR(ISDestroy(&cis[i]))
2367
2368        CHKERR(PetscFree(cis))
2369        CHKERR(PetscFree(cdm))
2370        CHKERR(PetscFree(cnamelist))
2371
2372        return (names, isets, dms)
2373
2374    def setSNESFunction(
2375        self, function: SNESFunction,
2376        args: tuple[Any, ...] | None = None,
2377        kargs: dict[str, Any] | None = None) -> None:
2378
2379        """Set `SNES` residual evaluation function.
2380
2381        Not collective.
2382
2383        Parameters
2384        ----------
2385        function
2386            The callback.
2387        args
2388            Positional arguments for the callback.
2389        kargs
2390            Keyword arguments for the callback.
2391
2392        See Also
2393        --------
2394        SNES.setFunction, petsc.DMSNESSetFunction
2395
2396        """
2397        if function is not None:
2398            if args  is None: args  = ()
2399            if kargs is None: kargs = {}
2400            context = (function, args, kargs)
2401            self.set_attr('__function__', context)
2402            CHKERR(DMSNESSetFunction(self.dm, SNES_Function, <void*>context))
2403        else:
2404            CHKERR(DMSNESSetFunction(self.dm, NULL, NULL))
2405
2406    def setSNESJacobian(
2407        self, jacobian: SNESJacobianFunction,
2408        args: tuple[Any, ...] | None = None,
2409        kargs: dict[str, Any] | None = None) -> None:
2410        """Set the `SNES` Jacobian evaluation function.
2411
2412        Not collective.
2413
2414        Parameters
2415        ----------
2416        jacobian
2417            The Jacobian callback.
2418        args
2419            Positional arguments for the callback.
2420        kargs
2421            Keyword arguments for the callback.
2422
2423        See Also
2424        --------
2425        SNES.setJacobian, petsc.DMSNESSetJacobian
2426
2427        """
2428        if jacobian is not None:
2429            if args  is None: args  = ()
2430            if kargs is None: kargs = {}
2431            context = (jacobian, args, kargs)
2432            self.set_attr('__jacobian__', context)
2433            CHKERR(DMSNESSetJacobian(self.dm, SNES_Jacobian, <void*>context))
2434        else:
2435            CHKERR(DMSNESSetJacobian(self.dm, NULL, NULL))
2436
2437    def addCoarsenHook(
2438        self,
2439        coarsenhook: DMCoarsenHookFunction,
2440        restricthook: DMRestrictHookFunction,
2441        args: tuple[Any, ...] | None = None,
2442        kargs: dict[str, Any] | None = None) -> None:
2443        """Add a callback to be executed when restricting to a coarser grid.
2444
2445        Logically collective.
2446
2447        Parameters
2448        ----------
2449        coarsenhook
2450            The coarsen hook function.
2451        restricthook
2452            The restrict hook function.
2453        args
2454            Positional arguments for the hooks.
2455        kargs
2456            Keyword arguments for the hooks.
2457
2458        See Also
2459        --------
2460        petsc.DMCoarsenHookAdd
2461
2462        """
2463        if args  is None: args  = ()
2464        if kargs is None: kargs = {}
2465
2466        if coarsenhook is not None:
2467            coarsencontext = (coarsenhook, args, kargs)
2468
2469            coarsenhooks = self.get_attr('__coarsenhooks__')
2470            if coarsenhooks is None:
2471                coarsenhooks = [coarsencontext]
2472                CHKERR(DMCoarsenHookAdd(self.dm, DM_PyCoarsenHook, NULL, <void*>NULL))
2473            else:
2474                coarsenhooks.append(coarsencontext)
2475            self.set_attr('__coarsenhooks__', coarsenhooks)
2476
2477        if restricthook is not None:
2478            restrictcontext = (restricthook, args, kargs)
2479
2480            restricthooks = self.get_attr('__restricthooks__')
2481            if restricthooks is None:
2482                restricthooks = [restrictcontext]
2483                CHKERR(DMCoarsenHookAdd(self.dm, NULL, DM_PyRestrictHook, <void*>NULL))
2484            else:
2485                restricthooks.append(restrictcontext)
2486            self.set_attr('__restricthooks__', restricthooks)
2487
2488    # --- application context ---
2489
2490    property appctx:
2491        """Application context."""
2492        def __get__(self) -> object:
2493            return self.getAppCtx()
2494
2495        def __set__(self, value):
2496            self.setAppCtx(value)
2497
2498    # --- discretization space ---
2499
2500    property ds:
2501        """Discrete space."""
2502        def __get__(self) -> DS:
2503            return self.getDS()
2504
2505        def __set__(self, value):
2506            self.setDS(value)
2507
2508# --------------------------------------------------------------------
2509
2510del DMType
2511del DMBoundaryType
2512del DMPolytopeType
2513del DMReorderDefaultFlag
2514
2515# --------------------------------------------------------------------
2516