xref: /petsc/src/binding/petsc4py/src/petsc4py/PETSc/petscdm.pxi (revision b5f0bcd6e9e8ed97648738542f5163d94f7b1782)
1# --------------------------------------------------------------------
2
3cdef extern from * nogil:
4
5    ctypedef const char* PetscDMType "DMType"
6    PetscDMType DMDA_type "DMDA"
7    PetscDMType DMCOMPOSITE
8    PetscDMType DMSLICED
9    PetscDMType DMSHELL
10    PetscDMType DMPLEX
11    PetscDMType DMREDUNDANT
12    PetscDMType DMPATCH
13    PetscDMType DMMOAB
14    PetscDMType DMNETWORK
15    PetscDMType DMFOREST
16    PetscDMType DMP4EST
17    PetscDMType DMP8EST
18    PetscDMType DMSWARM
19    PetscDMType DMPRODUCT
20    PetscDMType DMSTAG
21
22    ctypedef enum PetscDMBoundaryType"DMBoundaryType":
23        DM_BOUNDARY_NONE
24        DM_BOUNDARY_GHOSTED
25        DM_BOUNDARY_MIRROR
26        DM_BOUNDARY_PERIODIC
27        DM_BOUNDARY_TWIST
28
29    ctypedef enum PetscDMPolytopeType "DMPolytopeType":
30        DM_POLYTOPE_POINT
31        DM_POLYTOPE_SEGMENT
32        DM_POLYTOPE_POINT_PRISM_TENSOR
33        DM_POLYTOPE_TRIANGLE
34        DM_POLYTOPE_QUADRILATERAL
35        DM_POLYTOPE_SEG_PRISM_TENSOR
36        DM_POLYTOPE_TETRAHEDRON
37        DM_POLYTOPE_HEXAHEDRON
38        DM_POLYTOPE_TRI_PRISM
39        DM_POLYTOPE_TRI_PRISM_TENSOR
40        DM_POLYTOPE_QUAD_PRISM_TENSOR
41        DM_POLYTOPE_PYRAMID
42        DM_POLYTOPE_FV_GHOST
43        DM_POLYTOPE_INTERIOR_GHOST
44        DM_POLYTOPE_UNKNOWN
45        DM_POLYTOPE_UNKNOWN_CELL
46        DM_POLYTOPE_UNKNOWN_FACE
47        DM_NUM_POLYTOPES
48
49    ctypedef enum PetscDMReorderDefaultFlag "DMReorderDefaultFlag":
50        DM_REORDER_DEFAULT_NOTSET
51        DM_REORDER_DEFAULT_FALSE
52        DM_REORDER_DEFAULT_TRUE
53
54    ctypedef PetscErrorCode (*PetscDMCoarsenHook)(PetscDM,
55                                                  PetscDM,
56                                                  void*) except PETSC_ERR_PYTHON
57    ctypedef PetscErrorCode (*PetscDMRestrictHook)(PetscDM,
58                                                   PetscMat,
59                                                   PetscVec,
60                                                   PetscMat,
61                                                   PetscDM,
62                                                   void*) except PETSC_ERR_PYTHON
63
64    PetscErrorCode DMCreate(MPI_Comm, PetscDM*)
65    PetscErrorCode DMClone(PetscDM, PetscDM*)
66    PetscErrorCode DMDestroy(PetscDM*)
67    PetscErrorCode DMView(PetscDM, PetscViewer)
68    PetscErrorCode DMLoad(PetscDM, PetscViewer)
69    PetscErrorCode DMSetType(PetscDM, PetscDMType)
70    PetscErrorCode DMGetType(PetscDM, PetscDMType*)
71    PetscErrorCode DMGetDimension(PetscDM, PetscInt*)
72    PetscErrorCode DMSetDimension(PetscDM, PetscInt)
73    PetscErrorCode DMSetOptionsPrefix(PetscDM, char[])
74    PetscErrorCode DMGetOptionsPrefix(PetscDM, char*[])
75    PetscErrorCode DMAppendOptionsPrefix(PetscDM, char[])
76    PetscErrorCode DMSetFromOptions(PetscDM)
77    PetscErrorCode DMSetUp(PetscDM)
78
79    PetscErrorCode DMGetAdjacency(PetscDM, PetscInt, PetscBool*, PetscBool*)
80    PetscErrorCode DMSetAdjacency(PetscDM, PetscInt, PetscBool, PetscBool)
81    PetscErrorCode DMGetBasicAdjacency(PetscDM, PetscBool*, PetscBool*)
82    PetscErrorCode DMSetBasicAdjacency(PetscDM, PetscBool, PetscBool)
83
84    PetscErrorCode DMSetNumFields(PetscDM, PetscInt)
85    PetscErrorCode DMGetNumFields(PetscDM, PetscInt*)
86    PetscErrorCode DMSetField(PetscDM, PetscInt, PetscDMLabel, PetscObject)
87    PetscErrorCode DMAddField(PetscDM, PetscDMLabel, PetscObject)
88    PetscErrorCode DMGetField(PetscDM, PetscInt, PetscDMLabel*, PetscObject*)
89    PetscErrorCode DMClearFields(PetscDM)
90    PetscErrorCode DMCopyFields(PetscDM, PetscInt, PetscInt, PetscDM)
91    PetscErrorCode DMCreateDS(PetscDM)
92    PetscErrorCode DMClearDS(PetscDM)
93    PetscErrorCode DMGetDS(PetscDM, PetscDS*)
94    PetscErrorCode DMCopyDS(PetscDM, PetscInt, PetscInt, PetscDM)
95    PetscErrorCode DMCopyDisc(PetscDM, PetscDM)
96
97    PetscErrorCode DMGetBlockSize(PetscDM, PetscInt*)
98    PetscErrorCode DMSetVecType(PetscDM, PetscVecType)
99    PetscErrorCode DMCreateLocalVector(PetscDM, PetscVec*)
100    PetscErrorCode DMCreateGlobalVector(PetscDM, PetscVec*)
101    PetscErrorCode DMGetLocalVector(PetscDM, PetscVec*)
102    PetscErrorCode DMRestoreLocalVector(PetscDM, PetscVec*)
103    PetscErrorCode DMGetGlobalVector(PetscDM, PetscVec*)
104    PetscErrorCode DMRestoreGlobalVector(PetscDM, PetscVec*)
105    PetscErrorCode DMGetNamedLocalVector(PetscDM, const char[], PetscVec*)
106    PetscErrorCode DMRestoreNamedLocalVector(PetscDM, const char[], PetscVec*)
107    PetscErrorCode DMGetNamedGlobalVector(PetscDM, const char[], PetscVec*)
108    PetscErrorCode DMRestoreNamedGlobalVector(PetscDM, const char[], PetscVec*)
109    PetscErrorCode DMSetMatType(PetscDM, PetscMatType)
110    PetscErrorCode DMCreateMatrix(PetscDM, PetscMat*)
111    PetscErrorCode DMCreateMassMatrix(PetscDM, PetscDM, PetscMat*)
112
113    PetscErrorCode DMGetCoordinateDM(PetscDM, PetscDM*)
114    PetscErrorCode DMGetCoordinateSection(PetscDM, PetscSection*)
115    PetscErrorCode DMSetCoordinates(PetscDM, PetscVec)
116    PetscErrorCode DMGetCoordinates(PetscDM, PetscVec*)
117    PetscErrorCode DMSetCoordinatesLocal(PetscDM, PetscVec)
118    PetscErrorCode DMGetCoordinatesLocal(PetscDM, PetscVec*)
119    PetscErrorCode DMGetCoordinateDim(PetscDM, PetscInt*)
120    PetscErrorCode DMSetCoordinateDim(PetscDM, PetscInt)
121    PetscErrorCode DMLocalizeCoordinates(PetscDM)
122    PetscErrorCode DMSetCoordinateDisc(PetscDM, PetscFE, PetscBool, PetscBool)
123    PetscErrorCode DMSetCellCoordinateDM(PetscDM, PetscDM)
124    PetscErrorCode DMGetCellCoordinateDM(PetscDM, PetscDM*)
125    PetscErrorCode DMSetCellCoordinateSection(PetscDM, PetscInt, PetscSection)
126    PetscErrorCode DMGetCellCoordinateSection(PetscDM, PetscSection*)
127    PetscErrorCode DMSetCellCoordinates(PetscDM, PetscVec)
128    PetscErrorCode DMGetCellCoordinates(PetscDM, PetscVec*)
129    PetscErrorCode DMSetCellCoordinatesLocal(PetscDM, PetscVec)
130    PetscErrorCode DMGetCellCoordinatesLocal(PetscDM, PetscVec*)
131    PetscErrorCode DMGetCoordinatesLocalized(PetscDM, PetscBool*)
132    PetscErrorCode DMGetPeriodicity(PetscDM, const PetscReal *[], const PetscReal *[], const PetscReal *[])
133    PetscErrorCode DMSetPeriodicity(PetscDM, const PetscReal[], const PetscReal[], const PetscReal[])
134
135    PetscErrorCode DMCreateInterpolation(PetscDM, PetscDM, PetscMat*, PetscVec*)
136    PetscErrorCode DMCreateInjection(PetscDM, PetscDM, PetscMat*)
137    PetscErrorCode DMCreateRestriction(PetscDM, PetscDM, PetscMat*)
138
139    PetscErrorCode DMConvert(PetscDM, PetscDMType, PetscDM*)
140    PetscErrorCode DMRefine(PetscDM, MPI_Comm, PetscDM*)
141    PetscErrorCode DMCoarsen(PetscDM, MPI_Comm, PetscDM*)
142    PetscErrorCode DMRefineHierarchy(PetscDM, PetscInt, PetscDM[])
143    PetscErrorCode DMCoarsenHierarchy(PetscDM, PetscInt, PetscDM[])
144    PetscErrorCode DMGetRefineLevel(PetscDM, PetscInt*)
145    PetscErrorCode DMSetRefineLevel(PetscDM, PetscInt)
146    PetscErrorCode DMGetCoarsenLevel(PetscDM, PetscInt*)
147    PetscErrorCode DMGetCoarseDM(PetscDM, PetscDM*)
148    PetscErrorCode DMSetCoarseDM(PetscDM, PetscDM)
149
150    PetscErrorCode DMAdaptLabel(PetscDM, PetscDMLabel, PetscDM*)
151    PetscErrorCode DMAdaptMetric(PetscDM, PetscVec, PetscDMLabel, PetscDMLabel, PetscDM*)
152
153    PetscErrorCode DMGlobalToLocalBegin(PetscDM, PetscVec, PetscInsertMode, PetscVec)
154    PetscErrorCode DMGlobalToLocalEnd(PetscDM, PetscVec, PetscInsertMode, PetscVec)
155    PetscErrorCode DMLocalToGlobalBegin(PetscDM, PetscVec, PetscInsertMode, PetscVec)
156    PetscErrorCode DMLocalToGlobalEnd(PetscDM, PetscVec, PetscInsertMode, PetscVec)
157    PetscErrorCode DMLocalToLocalBegin(PetscDM, PetscVec, PetscInsertMode, PetscVec)
158    PetscErrorCode DMLocalToLocalEnd(PetscDM, PetscVec, PetscInsertMode, PetscVec)
159
160    PetscErrorCode DMGetLocalToGlobalMapping(PetscDM, PetscLGMap*)
161
162    PetscErrorCode DMSetLocalSection(PetscDM, PetscSection)
163    PetscErrorCode DMGetLocalSection(PetscDM, PetscSection*)
164    PetscErrorCode DMSetGlobalSection(PetscDM, PetscSection)
165    PetscErrorCode DMGetGlobalSection(PetscDM, PetscSection*)
166    PetscErrorCode DMCreateSectionSF(PetscDM, PetscSection, PetscSection)
167    PetscErrorCode DMGetSectionSF(PetscDM, PetscSF*)
168    PetscErrorCode DMSetSectionSF(PetscDM, PetscSF)
169    PetscErrorCode DMGetPointSF(PetscDM, PetscSF*)
170    PetscErrorCode DMSetPointSF(PetscDM, PetscSF)
171    PetscErrorCode DMGetUseNatural(PetscDM, PetscBool *)
172    PetscErrorCode DMSetUseNatural(PetscDM, PetscBool)
173
174    PetscErrorCode DMCreateSubDM(PetscDM, PetscInt, const PetscInt[], PetscIS*, PetscDM*)
175    PetscErrorCode DMSetAuxiliaryVec(PetscDM, PetscDMLabel, PetscInt, PetscInt, PetscVec)
176    PetscErrorCode DMGetAuxiliaryVec(PetscDM, PetscDMLabel, PetscInt, PetscInt, PetscVec*)
177
178    PetscErrorCode DMCreateLabel(PetscDM, const char[])
179    PetscErrorCode DMGetLabelValue(PetscDM, const char[], PetscInt, PetscInt*)
180    PetscErrorCode DMSetLabelValue(PetscDM, const char[], PetscInt, PetscInt)
181    PetscErrorCode DMHasLabel(PetscDM, const char[], PetscBool*)
182    PetscErrorCode DMClearLabelValue(PetscDM, const char[], PetscInt, PetscInt)
183    PetscErrorCode DMGetLabelSize(PetscDM, const char[], PetscInt*)
184    PetscErrorCode DMGetLabelIdIS(PetscDM, const char[], PetscIS*)
185    PetscErrorCode DMGetStratumSize(PetscDM, const char[], PetscInt, PetscInt*)
186    PetscErrorCode DMGetStratumIS(PetscDM, const char[], PetscInt, PetscIS*)
187    PetscErrorCode DMClearLabelStratum(PetscDM, const char[], PetscInt)
188    PetscErrorCode DMSetLabelOutput(PetscDM, const char[], PetscBool)
189    PetscErrorCode DMGetLabelOutput(PetscDM, const char[], PetscBool*)
190    PetscErrorCode DMGetNumLabels(PetscDM, PetscInt*)
191    PetscErrorCode DMGetLabelName(PetscDM, PetscInt, const char**)
192    PetscErrorCode DMHasLabel(PetscDM, const char[], PetscBool*)
193    PetscErrorCode DMGetLabel(PetscDM, const char*, PetscDMLabel*)
194    PetscErrorCode DMAddLabel(PetscDM, PetscDMLabel)
195    PetscErrorCode DMRemoveLabel(PetscDM, const char[], PetscDMLabel*)
196    PetscErrorCode DMLabelDestroy(PetscDMLabel *)
197
198    PetscErrorCode DMShellSetGlobalVector(PetscDM, PetscVec)
199    PetscErrorCode DMShellSetLocalVector(PetscDM, PetscVec)
200
201    PetscErrorCode DMKSPSetComputeOperators(PetscDM, PetscKSPComputeOpsFunction, void*)
202
203    PetscErrorCode DMCreateFieldDecomposition(PetscDM, PetscInt*, char***, PetscIS**, PetscDM**)
204
205    PetscErrorCode DMSNESSetFunction(PetscDM, PetscSNESFunctionFunction, void*)
206    PetscErrorCode DMSNESSetJacobian(PetscDM, PetscSNESJacobianFunction, void*)
207
208    PetscErrorCode DMCoarsenHookAdd(PetscDM, PetscDMCoarsenHook, PetscDMRestrictHook, void*)
209
210# --------------------------------------------------------------------
211
212cdef inline PetscDMBoundaryType asBoundaryType(object boundary) \
213    except <PetscDMBoundaryType>(-1):
214    if boundary is None:
215        return DM_BOUNDARY_NONE
216    if boundary is False:
217        return DM_BOUNDARY_NONE
218    if boundary is True:
219        return DM_BOUNDARY_PERIODIC
220    if isinstance(boundary, str):
221        if boundary == 'none':
222            return DM_BOUNDARY_NONE
223        elif boundary == 'ghosted':
224            return DM_BOUNDARY_GHOSTED
225        elif boundary == 'mirror':
226            return DM_BOUNDARY_MIRROR
227        elif boundary == 'periodic':
228            return DM_BOUNDARY_PERIODIC
229        elif boundary == 'twist':
230            return DM_BOUNDARY_TWIST
231        else:
232            raise ValueError("unknown boundary type: %s" % boundary)
233    return boundary
234
235cdef inline PetscInt asBoundary(object boundary,
236                                PetscDMBoundaryType *_x,
237                                PetscDMBoundaryType *_y,
238                                PetscDMBoundaryType *_z) except -1:
239    cdef PetscInt dim = 0
240    cdef object x=None, y=None, z=None
241    # Use `type(0)` instead of `int` to workaround
242    # Cython 3.1 failing to interpret `int` as a type
243    cdef type pyint = type(0)
244    if boundary is None or \
245       isinstance(boundary, str) or \
246       isinstance(boundary, pyint):
247        _x[0] = _y[0] = _z[0] = asBoundaryType(boundary)
248    else:
249        _x[0] = _y[0] = _z[0] = DM_BOUNDARY_NONE
250        boundary = tuple(boundary)
251        dim = <PetscInt>len(boundary)
252        if   dim == 0: pass
253        elif dim == 1: (x,) = boundary
254        elif dim == 2: (x, y) = boundary
255        elif dim == 3: (x, y, z) = boundary
256        if dim >= 1: _x[0] = asBoundaryType(x)
257        if dim >= 2: _y[0] = asBoundaryType(y)
258        if dim >= 3: _z[0] = asBoundaryType(z)
259    return dim
260
261cdef inline object toBoundary(PetscInt dim,
262                              PetscDMBoundaryType x,
263                              PetscDMBoundaryType y,
264                              PetscDMBoundaryType z):
265    if   dim == 0: return ()
266    elif dim == 1: return (x,)
267    elif dim == 2: return (x, y)
268    elif dim == 3: return (x, y, z)
269
270# -----------------------------------------------------------------------------
271
272cdef inline DM ref_DM(PetscDM dm):
273    cdef DM ob = <DM> DM()
274    ob.dm = dm
275    CHKERR(PetscINCREF(ob.obj))
276    return ob
277
278# --------------------------------------------------------------------
279
280cdef PetscErrorCode DM_PyCoarsenHook(
281    PetscDM fine,
282    PetscDM coarse,
283    void    *ctx,
284   ) except PETSC_ERR_PYTHON with gil:
285
286    cdef DM Fine = ref_DM(fine)
287    cdef DM Coarse = ref_DM(coarse)
288    cdef object hooks = Fine.get_attr('__coarsenhooks__')
289    assert hooks is not None and type(hooks) is list
290    for hook in hooks:
291        (hookop, args, kargs) = hook
292        hookop(Fine, Coarse, *args, **kargs)
293    return PETSC_SUCCESS
294
295cdef PetscErrorCode DM_PyRestrictHook(
296    PetscDM  fine,
297    PetscMat mrestrict,
298    PetscVec rscale,
299    PetscMat inject,
300    PetscDM  coarse,
301    void     *ctx,
302   ) except PETSC_ERR_PYTHON with gil:
303
304    cdef DM  Fine = ref_DM(fine)
305    cdef Mat Mrestrict = ref_Mat(mrestrict)
306    cdef Vec Rscale = ref_Vec(rscale)
307    cdef Mat Inject = ref_Mat(inject)
308    cdef DM  Coarse = ref_DM(coarse)
309    cdef object hooks = Fine.get_attr('__restricthooks__')
310    assert hooks is not None and type(hooks) is list
311    for hook in hooks:
312        (hookop, args, kargs) = hook
313        hookop(Fine, Mrestrict, Rscale, Inject, Coarse, *args, **kargs)
314    return PETSC_SUCCESS
315# -----------------------------------------------------------------------------
316