xref: /petsc/src/binding/petsc4py/src/petsc4py/PETSc/FE.pyx (revision b5f0bcd6e9e8ed97648738542f5163d94f7b1782)
1# --------------------------------------------------------------------
2
3class FEType(object):
4    """The finite element types."""
5    BASIC     = S_(PETSCFEBASIC)
6    OPENCL    = S_(PETSCFEOPENCL)
7    COMPOSITE = S_(PETSCFECOMPOSITE)
8
9# --------------------------------------------------------------------
10
11
12cdef class FE(Object):
13    """A PETSc object that manages a finite element space."""
14
15    Type = FEType
16
17    def __cinit__(self):
18        self.obj = <PetscObject*> &self.fe
19        self.fe = NULL
20
21    def view(self, Viewer viewer=None) -> None:
22        """View a `FE` object.
23
24        Collective.
25
26        Parameters
27        ----------
28        viewer
29            A `Viewer` to display the graph.
30
31        See Also
32        --------
33        petsc.PetscFEView
34
35        """
36        cdef PetscViewer vwr = NULL
37        if viewer is not None: vwr = viewer.vwr
38        CHKERR(PetscFEView(self.fe, vwr))
39
40    def destroy(self) -> Self:
41        """Destroy the `FE` object.
42
43        Collective.
44
45        See Also
46        --------
47        petsc.PetscFEDestroy
48
49        """
50        CHKERR(PetscFEDestroy(&self.fe))
51        return self
52
53    def create(self, comm: Comm | None = None) -> Self:
54        """Create an empty `FE` object.
55
56        Collective.
57
58        The type can then be set with `setType`.
59
60        Parameters
61        ----------
62        comm
63            MPI communicator, defaults to `Sys.getDefaultComm`.
64
65        See Also
66        --------
67        setType, petsc.PetscFECreate
68
69        """
70        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
71        cdef PetscFE newfe = NULL
72        CHKERR(PetscFECreate(ccomm, &newfe))
73        CHKERR(PetscCLEAR(self.obj)); self.fe = newfe
74        return self
75
76    def createDefault(
77        self,
78        dim: int,
79        nc: int,
80        isSimplex: bool,
81        qorder: int = DETERMINE,
82        prefix: str | None = None,
83        comm: Comm | None = None) -> Self:
84        """Create a `FE` for basic FEM computation.
85
86        Collective.
87
88        Parameters
89        ----------
90        dim
91            The spatial dimension.
92        nc
93            The number of components.
94        isSimplex
95            Flag for simplex reference cell, otherwise it's a tensor product.
96        qorder
97            The quadrature order or `DETERMINE` to use `Space` polynomial
98            degree.
99        prefix
100            The options prefix, or `None`.
101        comm
102            MPI communicator, defaults to `Sys.getDefaultComm`.
103
104        See Also
105        --------
106        petsc.PetscFECreateDefault
107
108        """
109        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
110        cdef PetscFE newfe = NULL
111        cdef PetscInt cdim = asInt(dim)
112        cdef PetscInt cnc = asInt(nc)
113        cdef PetscInt cqorder = asInt(qorder)
114        cdef PetscBool cisSimplex = asBool(isSimplex)
115        cdef const char *cprefix = NULL
116        if prefix:
117            prefix = str2bytes(prefix, &cprefix)
118        CHKERR(PetscFECreateDefault(ccomm, cdim, cnc, cisSimplex, cprefix, cqorder, &newfe))
119        CHKERR(PetscCLEAR(self.obj)); self.fe = newfe
120        return self
121
122    def createByCell(
123        self,
124        dim: int,
125        nc: int,
126        ctype: DM.PolytopeType,
127        qorder: int = DETERMINE,
128        prefix: str | None = None,
129        comm: Comm | None = None) -> Self:
130        """Create a `FE` for basic FEM computation.
131
132        Collective.
133
134        Parameters
135        ----------
136        dim
137            The spatial dimension.
138        nc
139            The number of components.
140        ctype
141            The cell type.
142        qorder
143            The quadrature order or `DETERMINE` to use `Space` polynomial
144            degree.
145        prefix
146            The options prefix, or `None`.
147        comm
148            MPI communicator, defaults to `Sys.getDefaultComm`.
149
150        See Also
151        --------
152        petsc.PetscFECreateByCell
153
154        """
155        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_SELF)
156        cdef PetscFE newfe = NULL
157        cdef PetscInt cdim = asInt(dim)
158        cdef PetscInt cnc = asInt(nc)
159        cdef PetscInt cqorder = asInt(qorder)
160        cdef PetscDMPolytopeType cCellType = ctype
161        cdef const char *cprefix = NULL
162        if prefix:
163            prefix = str2bytes(prefix, &cprefix)
164        CHKERR(PetscFECreateDefault(ccomm, cdim, cnc, cCellType, cprefix, cqorder, &newfe))
165        CHKERR(PetscCLEAR(self.obj)); self.fe = newfe
166        return self
167
168    def createLagrange(
169        self,
170        dim: int,
171        nc: int,
172        isSimplex: bool,
173        k: int,
174        qorder: int = DETERMINE,
175        comm: Comm | None = None) -> Self:
176        """Create a `FE` for the basic Lagrange space of degree k.
177
178        Collective.
179
180        Parameters
181        ----------
182        dim
183            The spatial dimension.
184        nc
185            The number of components.
186        isSimplex
187            Flag for simplex reference cell, otherwise it's a tensor product.
188        k
189            The degree of the space.
190        qorder
191            The quadrature order or `DETERMINE` to use `Space` polynomial
192            degree.
193        comm
194            MPI communicator, defaults to `Sys.getDefaultComm`.
195
196        See Also
197        --------
198        petsc.PetscFECreateLagrange
199
200        """
201        cdef MPI_Comm ccomm = def_Comm(comm, PETSC_COMM_DEFAULT)
202        cdef PetscFE newfe = NULL
203        cdef PetscInt cdim = asInt(dim)
204        cdef PetscInt cnc = asInt(nc)
205        cdef PetscInt ck = asInt(k)
206        cdef PetscInt cqorder = asInt(qorder)
207        cdef PetscBool cisSimplex = asBool(isSimplex)
208        CHKERR(PetscFECreateLagrange(ccomm, cdim, cnc, cisSimplex, ck, cqorder, &newfe))
209        CHKERR(PetscCLEAR(self.obj)); self.fe = newfe
210        return self
211
212    def getQuadrature(self) -> Quad:
213        """Return the `Quad` used to calculate inner products.
214
215        Not collective.
216
217        See Also
218        --------
219        setQuadrature, petsc.PetscFEGetQuadrature
220
221        """
222        cdef Quad quad = Quad()
223        CHKERR(PetscFEGetQuadrature(self.fe, &quad.quad))
224        CHKERR(PetscINCREF(quad.obj))
225        return quad
226
227    def getDimension(self) -> int:
228        """Return the dimension of the finite element space on a cell.
229
230        Not collective.
231
232        See Also
233        --------
234        petsc.PetscFEGetDimension
235
236        """
237        cdef PetscInt cdim = 0
238        CHKERR(PetscFEGetDimension(self.fe, &cdim))
239        return toInt(cdim)
240
241    def getSpatialDimension(self) -> int:
242        """Return the spatial dimension of the element.
243
244        Not collective.
245
246        See Also
247        --------
248        petsc.PetscFEGetSpatialDimension
249
250        """
251        cdef PetscInt csdim = 0
252        CHKERR(PetscFEGetSpatialDimension(self.fe, &csdim))
253        return toInt(csdim)
254
255    def getNumComponents(self) -> int:
256        """Return the number of components in the element.
257
258        Not collective.
259
260        See Also
261        --------
262        setNumComponents, petsc.PetscFEGetNumComponents
263
264        """
265        cdef PetscInt comp = 0
266        CHKERR(PetscFEGetNumComponents(self.fe, &comp))
267        return toInt(comp)
268
269    def setNumComponents(self, comp: int) -> None:
270        """Set the number of field components in the element.
271
272        Not collective.
273
274        Parameters
275        ----------
276        comp
277            The number of field components.
278
279        See Also
280        --------
281        getNumComponents, petsc.PetscFESetNumComponents
282
283        """
284        cdef PetscInt ccomp = asInt(comp)
285        CHKERR(PetscFESetNumComponents(self.fe, ccomp))
286
287    def getNumDof(self) -> ndarray:
288        """Return the number of DOFs.
289
290        Not collective.
291
292        Return the number of DOFs (dual basis vectors) associated with mesh
293        points on the reference cell of a given dimension.
294
295        See Also
296        --------
297        petsc.PetscFEGetNumDof
298
299        """
300        cdef const PetscInt *numDof = NULL
301        cdef PetscInt cdim = 0
302        CHKERR(PetscFEGetDimension(self.fe, &cdim))
303        CHKERR(PetscFEGetNumDof(self.fe, &numDof))
304        return array_i(cdim, numDof)
305
306    def getTileSizes(self) -> tuple[int, int, int, int]:
307        """Return the tile sizes for evaluation.
308
309        Not collective.
310
311        Returns
312        -------
313        blockSize : int
314            The number of elements in a block.
315        numBlocks : int
316            The number of blocks in a batch.
317        batchSize : int
318            The number of elements in a batch.
319        numBatches : int
320            The number of batches in a chunk.
321
322        See Also
323        --------
324        setTileSizes, petsc.PetscFEGetTileSizes
325
326        """
327        cdef PetscInt blockSize = 0, numBlocks = 0
328        cdef PetscInt batchSize = 0, numBatches = 0
329        CHKERR(PetscFEGetTileSizes(self.fe, &blockSize, &numBlocks, &batchSize, &numBatches))
330        return toInt(blockSize), toInt(numBlocks), toInt(batchSize), toInt(numBatches)
331
332    def setTileSizes(
333        self,
334        blockSize: int,
335        numBlocks: int,
336        batchSize: int,
337        numBatches: int) -> None:
338        """Set the tile sizes for evaluation.
339
340        Not collective.
341
342        Parameters
343        ----------
344        blockSize
345            The number of elements in a block.
346        numBlocks
347            The number of blocks in a batch.
348        batchSize
349            The number of elements in a batch.
350        numBatches
351            The number of batches in a chunk.
352
353        See Also
354        --------
355        getTileSizes, petsc.PetscFESetTileSizes
356
357        """
358        cdef PetscInt cblockSize = asInt(blockSize), cnumBlocks = asInt(numBlocks)
359        cdef PetscInt cbatchSize = asInt(batchSize), cnumBatches = asInt(numBatches)
360        CHKERR(PetscFESetTileSizes(self.fe, cblockSize, cnumBlocks, cbatchSize, cnumBatches))
361
362    def getFaceQuadrature(self) -> Quad:
363        """Return the `Quad` used to calculate inner products on faces.
364
365        Not collective.
366
367        See Also
368        --------
369        setFaceQuadrature, petsc.PetscFEGetFaceQuadrature
370
371        """
372        cdef Quad quad = Quad()
373        CHKERR(PetscFEGetFaceQuadrature(self.fe, &quad.quad))
374        CHKERR(PetscINCREF(quad.obj))
375        return quad
376
377    def setQuadrature(self, Quad quad) -> Self:
378        """Set the `Quad` used to calculate inner products.
379
380        Not collective.
381
382        Parameters
383        ----------
384        quad
385            The `Quad` object.
386
387        See Also
388        --------
389        getQuadrature, petsc.PetscFESetQuadrature
390
391        """
392        CHKERR(PetscFESetQuadrature(self.fe, quad.quad))
393        return self
394
395    def setFaceQuadrature(self, Quad quad) -> Quad:
396        """Set the `Quad` used to calculate inner products on faces.
397
398        Not collective.
399
400        Parameters
401        ----------
402        quad
403            The `Quad` object.
404
405        See Also
406        --------
407        getFaceQuadrature, petsc.PetscFESetFaceQuadrature
408
409        """
410        CHKERR(PetscFESetFaceQuadrature(self.fe, quad.quad))
411        return self
412
413    def setType(self, fe_type: Type | str) -> Self:
414        """Build a particular `FE`.
415
416        Collective.
417
418        Parameters
419        ----------
420        fe_type
421            The kind of FEM space.
422
423        See Also
424        --------
425        petsc.PetscFESetType
426
427        """
428        cdef PetscFEType cval = NULL
429        fe_type = str2bytes(fe_type, &cval)
430        CHKERR(PetscFESetType(self.fe, cval))
431        return self
432
433    def getBasisSpace(self) -> Space:
434        """Return the `Space` used for the approximation of the `FE` solution.
435
436        Not collective.
437
438        See Also
439        --------
440        setBasisSpace, petsc.PetscFEGetBasisSpace
441
442        """
443        cdef Space sp = Space()
444        CHKERR(PetscFEGetBasisSpace(self.fe, &sp.space))
445        CHKERR(PetscINCREF(sp.obj))
446        return sp
447
448    def setBasisSpace(self, Space sp) -> None:
449        """Set the `Space` used for the approximation of the solution.
450
451        Not collective.
452
453        Parameters
454        ----------
455        sp
456            The `Space` object.
457
458        See Also
459        --------
460        getBasisSpace, petsc.PetscFESetBasisSpace
461
462        """
463        CHKERR(PetscFESetBasisSpace(self.fe, sp.space))
464
465    def setFromOptions(self) -> None:
466        """Set parameters in a `FE` from the options database.
467
468        Collective.
469
470        See Also
471        --------
472        petsc_options, petsc.PetscFESetFromOptions
473
474        """
475        CHKERR(PetscFESetFromOptions(self.fe))
476
477    def setUp(self) -> None:
478        """Construct data structures for the `FE` after the `Type` has been set.
479
480        Collective.
481
482        See Also
483        --------
484        petsc.PetscFESetUp
485
486        """
487        CHKERR(PetscFESetUp(self.fe))
488
489    def getDualSpace(self) -> DualSpace:
490        """Return the `DualSpace` used to define the inner product for the `FE`.
491
492        Not collective.
493
494        See Also
495        --------
496        setDualSpace, DualSpace, petsc.PetscFEGetDualSpace
497
498        """
499        cdef DualSpace dspace = DualSpace()
500        CHKERR(PetscFEGetDualSpace(self.fe, &dspace.dualspace))
501        CHKERR(PetscINCREF(dspace.obj))
502        return dspace
503
504    def setDualSpace(self, DualSpace dspace) -> None:
505        """Set the `DualSpace` used to define the inner product.
506
507        Not collective.
508
509        Parameters
510        ----------
511        dspace
512            The `DualSpace` object.
513
514        See Also
515        --------
516        getDualSpace, DualSpace, petsc.PetscFESetDualSpace
517
518        """
519        CHKERR(PetscFESetDualSpace(self.fe, dspace.dualspace))
520
521# --------------------------------------------------------------------
522
523del FEType
524
525# --------------------------------------------------------------------
526