xref: /petsc/src/binding/petsc4py/src/petsc4py/PETSc/DMLabel.pyx (revision 552edb6364df478b294b3111f33a8f37ca096b20)
1
2cdef class DMLabel(Object):
3    """An object representing a subset of mesh entities from a `DM`."""
4    def __cinit__(self):
5        self.obj = <PetscObject*> &self.dmlabel
6        self.dmlabel  = NULL
7
8    def destroy(self) -> Self:
9        """Destroy the label.
10
11        Collective.
12
13        See Also
14        --------
15        petsc.DMLabelDestroy
16
17        """
18        CHKERR(DMLabelDestroy(&self.dmlabel))
19        return self
20
21    def view(self, Viewer viewer=None) -> None:
22        """View the label.
23
24        Collective.
25
26        Parameters
27        ----------
28        viewer
29            A `Viewer` to display the graph.
30
31        See Also
32        --------
33        petsc.DMLabelView
34
35        """
36        cdef PetscViewer vwr = NULL
37        if viewer is not None: vwr = viewer.vwr
38        CHKERR(DMLabelView(self.dmlabel, vwr))
39
40    def create(self, name: str, comm: Comm | None = None) -> Self:
41        """Create a `DMLabel` object, which is a multimap.
42
43        Collective.
44
45        Parameters
46        ----------
47        name
48            The label name.
49        comm
50            MPI communicator, defaults to `COMM_SELF`.
51
52        See Also
53        --------
54        petsc.DMLabelCreate
55
56        """
57        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_SELF)
58        cdef PetscDMLabel newdmlabel = NULL
59        cdef const char *cname = NULL
60        name = str2bytes(name, &cname)
61        CHKERR(DMLabelCreate(ccomm, cname, &newdmlabel))
62        CHKERR(PetscCLEAR(self.obj)); self.dmlabel = newdmlabel
63        return self
64
65    def duplicate(self) -> DMLabel:
66        """Duplicate the `DMLabel`.
67
68        Collective.
69
70        See Also
71        --------
72        petsc.DMLabelDuplicate
73
74        """
75        cdef DMLabel new = DMLabel()
76        CHKERR(DMLabelDuplicate(self.dmlabel, &new.dmlabel))
77        return new
78
79    def reset(self) -> None:
80        """Destroy internal data structures in the `DMLabel`.
81
82        Not collective.
83
84        See Also
85        --------
86        petsc.DMLabelReset
87
88        """
89        CHKERR(DMLabelReset(self.dmlabel))
90
91    def insertIS(self, IS iset, value: int) -> Self:
92        """Set all points in the `IS` to a value.
93
94        Not collective.
95
96        Parameters
97        ----------
98        iset
99            The point IS.
100        value
101            The point value.
102
103        See Also
104        --------
105        petsc.DMLabelInsertIS
106
107        """
108        cdef PetscInt cvalue = asInt(value)
109        CHKERR(DMLabelInsertIS(self.dmlabel, iset.iset, cvalue))
110        return self
111
112    def setValue(self, point: int, value: int) -> None:
113        """Set the value a label assigns to a point.
114
115        Not collective.
116
117        If the value is the same as the label's default value (which is
118        initially ``-1``, and can be changed with `setDefaultValue`), this
119        function will do nothing.
120
121        Parameters
122        ----------
123        point
124            The point.
125        value
126            The point value.
127
128        See Also
129        --------
130        getValue, setDefaultValue, petsc.DMLabelSetValue
131
132        """
133        cdef PetscInt cpoint = asInt(point)
134        cdef PetscInt cvalue = asInt(value)
135        CHKERR(DMLabelSetValue(self.dmlabel, cpoint, cvalue))
136
137    def getValue(self, point: int) -> int:
138        """Return the value a label assigns to a point.
139
140        Not collective.
141
142        If no value was assigned, a default value will be returned
143        The default value, initially ``-1``, can be changed with
144        `setDefaultValue`.
145
146        Parameters
147        ----------
148        point
149            The point.
150
151        See Also
152        --------
153        setValue, setDefaultValue, petsc.DMLabelGetValue
154
155        """
156        cdef PetscInt cpoint = asInt(point)
157        cdef PetscInt cvalue = 0
158        CHKERR(DMLabelGetValue(self.dmlabel, cpoint, &cvalue))
159        return toInt(cvalue)
160
161    def getDefaultValue(self) -> int:
162        """Return the default value returned by `getValue`.
163
164        Not collective.
165
166        The default value is returned if a point has not been explicitly given
167        a value. When a label is created, it is initialized to ``-1``.
168
169        See Also
170        --------
171        setDefaultValue, petsc.DMLabelGetDefaultValue
172
173        """
174        cdef PetscInt cvalue = 0
175        CHKERR(DMLabelGetDefaultValue(self.dmlabel, &cvalue))
176        return toInt(cvalue)
177
178    def setDefaultValue(self, value: int) -> None:
179        """Set the default value returned by `getValue`.
180
181        Not collective.
182
183        The value is used if a point has not been explicitly given a value.
184        When a label is created, the default value is initialized to ``-1``.
185
186        Parameters
187        ----------
188        value
189            The default value.
190
191        See Also
192        --------
193        getDefaultValue, petsc.DMLabelSetDefaultValue
194
195        """
196        cdef PetscInt cvalue = asInt(value)
197        CHKERR(DMLabelSetDefaultValue(self.dmlabel, cvalue))
198
199    def clearValue(self, point: int, value: int) -> None:
200        """Clear the value a label assigns to a point.
201
202        Not collective.
203
204        Parameters
205        ----------
206        point
207            The point.
208        value
209            The point value.
210
211        See Also
212        --------
213        petsc.DMLabelClearValue
214
215        """
216        cdef PetscInt cpoint = asInt(point)
217        cdef PetscInt cvalue = asInt(value)
218        CHKERR(DMLabelClearValue(self.dmlabel, cpoint, cvalue))
219
220    def addStratum(self, value: int) -> None:
221        """Add a new stratum value in a `DMLabel`.
222
223        Not collective.
224
225        Parameters
226        ----------
227        value
228            The stratum value.
229
230        See Also
231        --------
232        addStrata, addStrataIS, petsc.DMLabelAddStratum
233
234        """
235        cdef PetscInt cvalue = asInt(value)
236        CHKERR(DMLabelAddStratum(self.dmlabel, cvalue))
237
238    def addStrata(self, strata: Sequence[int]) -> None:
239        """Add new stratum values in a `DMLabel`.
240
241        Not collective.
242
243        Parameters
244        ----------
245        strata
246            The stratum values.
247
248        See Also
249        --------
250        addStrataIS, addStratum, petsc.DMLabelAddStrata
251
252        """
253        cdef PetscInt *istrata = NULL
254        cdef PetscInt numStrata = 0
255        strata = iarray_i(strata, &numStrata, &istrata)
256        CHKERR(DMLabelAddStrata(self.dmlabel, numStrata, istrata))
257
258    def addStrataIS(self, IS iset) -> None:
259        """Add new stratum values in a `DMLabel`.
260
261        Not collective.
262
263        Parameters
264        ----------
265        iset
266            Index set with stratum values.
267
268        See Also
269        --------
270        addStrata, addStratum, petsc.DMLabelAddStrataIS
271
272        """
273        CHKERR(DMLabelAddStrataIS(self.dmlabel, iset.iset))
274
275    def getNumValues(self) -> int:
276        """Return the number of values that the `DMLabel` takes.
277
278        Not collective.
279
280        See Also
281        --------
282        petsc.DMLabelGetNumValues
283
284        """
285        cdef PetscInt numValues = 0
286        CHKERR(DMLabelGetNumValues(self.dmlabel, &numValues))
287        return toInt(numValues)
288
289    def getValueIS(self) -> IS:
290        """Return an `IS` of all values that the `DMLabel` takes.
291
292        Not collective.
293
294        See Also
295        --------
296        petsc.DMLabelGetValueIS
297
298        """
299        cdef IS iset = IS()
300        CHKERR(DMLabelGetValueIS(self.dmlabel, &iset.iset))
301        return iset
302
303    def stratumHasPoint(self, value: int, point: int) -> bool:
304        """Return whether the stratum contains a point.
305
306        Not collective.
307
308        Parameters
309        ----------
310        value
311            The stratum value.
312        point
313            The point.
314
315        See Also
316        --------
317        petsc.DMLabelStratumHasPoint
318
319        """
320        cdef PetscInt cpoint = asInt(point)
321        cdef PetscInt cvalue = asInt(value)
322        cdef PetscBool ccontains = PETSC_FALSE
323        CHKERR(DMLabelStratumHasPoint(self.dmlabel, cvalue, cpoint, &ccontains))
324        return toBool(ccontains)
325
326    def hasStratum(self, value: int) -> bool:
327        """Determine whether points exist with the given value.
328
329        Not collective.
330
331        Parameters
332        ----------
333        value
334            The stratum value.
335
336        See Also
337        --------
338        petsc.DMLabelHasStratum
339
340        """
341        cdef PetscInt cvalue = asInt(value)
342        cdef PetscBool cexists = PETSC_FALSE
343        CHKERR(DMLabelHasStratum(self.dmlabel, cvalue, &cexists))
344        return toBool(cexists)
345
346    def getStratumSize(self, stratum: int) -> int:
347        """Return the size of a stratum.
348
349        Not collective.
350
351        Parameters
352        ----------
353        stratum
354            The stratum value.
355
356        See Also
357        --------
358        petsc.DMLabelGetStratumSize
359
360        """
361        cdef PetscInt cstratum = asInt(stratum)
362        cdef PetscInt csize = 0
363        CHKERR(DMLabelGetStratumSize(self.dmlabel, cstratum, &csize))
364        return toInt(csize)
365
366    def getStratumIS(self, stratum: int) -> IS:
367        """Return an `IS` with the stratum points.
368
369        Not collective.
370
371        Parameters
372        ----------
373        stratum
374            The stratum value.
375
376        See Also
377        --------
378        setStratumIS, petsc.DMLabelGetStratumIS
379
380        """
381        cdef PetscInt cstratum = asInt(stratum)
382        cdef IS iset = IS()
383        CHKERR(DMLabelGetStratumIS(self.dmlabel, cstratum, &iset.iset))
384        return iset
385
386    def setStratumIS(self, stratum: int, IS iset) -> None:
387        """Set the stratum points using an `IS`.
388
389        Not collective.
390
391        Parameters
392        ----------
393        stratum
394            The stratum value.
395        iset
396            The stratum points.
397
398        See Also
399        --------
400        getStratumIS, petsc.DMLabelSetStratumIS
401
402        """
403        cdef PetscInt cstratum = asInt(stratum)
404        CHKERR(DMLabelSetStratumIS(self.dmlabel, cstratum, iset.iset))
405
406    def clearStratum(self, stratum: int) -> None:
407        """Remove a stratum.
408
409        Not collective.
410
411        Parameters
412        ----------
413        stratum
414            The stratum value.
415
416        See Also
417        --------
418        petsc.DMLabelClearStratum
419
420        """
421        cdef PetscInt cstratum = asInt(stratum)
422        CHKERR(DMLabelClearStratum(self.dmlabel, cstratum))
423
424    def computeIndex(self) -> None:
425        """Create an index structure for membership determination.
426
427        Not collective.
428
429        Automatically determines the bounds.
430
431        See Also
432        --------
433        petsc.DMLabelComputeIndex
434
435        """
436        CHKERR(DMLabelComputeIndex(self.dmlabel))
437
438    def createIndex(self, pStart: int, pEnd: int) -> None:
439        """Create an index structure for membership determination.
440
441        Not collective.
442
443        Parameters
444        ----------
445        pStart
446            The smallest point.
447        pEnd
448            The largest point + 1.
449
450        See Also
451        --------
452        destroyIndex, petsc.DMLabelCreateIndex
453
454        """
455        cdef PetscInt cpstart = asInt(pStart), cpend = asInt(pEnd)
456        CHKERR(DMLabelCreateIndex(self.dmlabel, cpstart, cpend))
457
458    def destroyIndex(self) -> None:
459        """Destroy the index structure.
460
461        Not collective.
462
463        See Also
464        --------
465        createIndex, petsc.DMLabelDestroyIndex
466
467        """
468        CHKERR(DMLabelDestroyIndex(self.dmlabel))
469
470    def hasValue(self, value: int) -> bool:
471        """Determine whether a label assigns the value to any point.
472
473        Not collective.
474
475        Parameters
476        ----------
477        value
478            The value.
479
480        See Also
481        --------
482        hasPoint, petsc.DMLabelHasValue
483
484        """
485        cdef PetscInt cvalue = asInt(value)
486        cdef PetscBool cexists = PETSC_FALSE
487        CHKERR(DMLabelHasValue(self.dmlabel, cvalue, &cexists))
488        return toBool(cexists)
489
490    def hasPoint(self, point: int) -> bool:
491        """Determine whether the label contains a point.
492
493        Not collective.
494
495        The user must call `createIndex` before this function.
496
497        Parameters
498        ----------
499        point
500            The point.
501
502        See Also
503        --------
504        hasValue, petsc.DMLabelHasPoint
505
506        """
507        cdef PetscInt cpoint = asInt(point)
508        cdef PetscBool cexists = PETSC_FALSE
509        CHKERR(DMLabelHasPoint(self.dmlabel, cpoint, &cexists))
510        return toBool(cexists)
511
512    def getBounds(self) -> tuple[int, int]:
513        """Return the smallest and largest point in the label.
514
515        Not collective.
516
517        The returned values are the smallest point and the largest point + 1.
518
519        See Also
520        --------
521        petsc.DMLabelGetBounds
522
523        """
524        cdef PetscInt cpstart = 0, cpend = 0
525        CHKERR(DMLabelGetBounds(self.dmlabel, &cpstart, &cpend))
526        return toInt(cpstart), toInt(cpend)
527
528    def filter(self, start: int, end: int) -> None:
529        """Remove all points outside of [start, end).
530
531        Not collective.
532
533        Parameters
534        ----------
535        start
536            The first point kept.
537        end
538            One more than the last point kept.
539
540        See Also
541        --------
542        petsc.DMLabelFilter
543
544        """
545        cdef PetscInt cstart = asInt(start), cend = asInt(end)
546        CHKERR(DMLabelFilter(self.dmlabel, cstart, cend))
547
548    def permute(self, IS permutation) -> DMLabel:
549        """Create a new label with permuted points.
550
551        Not collective.
552
553        Parameters
554        ----------
555        permutation
556            The point permutation.
557
558        See Also
559        --------
560        petsc.DMLabelPermute
561
562        """
563        cdef DMLabel new = DMLabel()
564        CHKERR(DMLabelPermute(self.dmlabel, permutation.iset, &new.dmlabel))
565        return new
566
567    def distribute(self, SF sf) -> DMLabel:
568        """Create a new label pushed forward over the `SF`.
569
570        Collective.
571
572        Parameters
573        ----------
574        sf
575            The map from old to new distribution.
576
577        See Also
578        --------
579        gather, petsc.DMLabelDistribute
580
581        """
582        cdef DMLabel new = DMLabel()
583        CHKERR(DMLabelDistribute(self.dmlabel, sf.sf, &new.dmlabel))
584        return new
585
586    def gather(self, SF sf) -> DMLabel:
587        """Gather all label values from leaves into roots.
588
589        Collective.
590
591        This is the inverse operation to `distribute`.
592
593        Parameters
594        ----------
595        sf
596            The `SF` communication map.
597
598        See Also
599        --------
600        distribute, petsc.DMLabelGather
601
602        """
603        cdef DMLabel new = DMLabel()
604        CHKERR(DMLabelGather(self.dmlabel, sf.sf, &new.dmlabel))
605        return new
606
607    def convertToSection(self) -> tuple[Section, IS]:
608        """Return a `Section` and `IS` that encode the label.
609
610        Not collective.
611
612        See Also
613        --------
614        petsc.DMLabelConvertToSection
615
616        """
617        cdef Section section = Section()
618        cdef IS iset = IS()
619        CHKERR(DMLabelConvertToSection(self.dmlabel, &section.sec, &iset.iset))
620        return section, iset
621
622    def getNonEmptyStratumValuesIS(self) -> IS:
623        """Return an `IS` of all values that the `DMLabel` takes.
624
625        Not collective.
626
627        See Also
628        --------
629        petsc.DMLabelGetNonEmptyStratumValuesIS
630
631        """
632        cdef IS iset = IS()
633        CHKERR(DMLabelGetNonEmptyStratumValuesIS(self.dmlabel, &iset.iset))
634        return iset
635