xref: /petsc/src/binding/petsc4py/src/petsc4py/PETSc/petscdmshell.pxi (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
1ctypedef PetscErrorCode (*PetscDMShellXToYFunction)(PetscDM,
2                                                    PetscVec,
3                                                    PetscInsertMode,
4                                                    PetscVec) except PETSC_ERR_PYTHON
5cdef extern from * nogil:
6    ctypedef PetscErrorCode (*PetscDMShellCreateVectorFunction)(PetscDM,
7                                                                PetscVec*) except PETSC_ERR_PYTHON
8    ctypedef PetscErrorCode (*PetscDMShellCreateMatrixFunction)(PetscDM,
9                                                                PetscMat*) except PETSC_ERR_PYTHON
10    ctypedef PetscErrorCode (*PetscDMShellTransferFunction)(PetscDM,
11                                                            MPI_Comm,
12                                                            PetscDM*) except PETSC_ERR_PYTHON
13    ctypedef PetscErrorCode (*PetscDMShellCreateInterpolationFunction)(PetscDM,
14                                                                       PetscDM,
15                                                                       PetscMat*,
16                                                                       PetscVec*) except PETSC_ERR_PYTHON
17    ctypedef PetscErrorCode (*PetscDMShellCreateInjectionFunction)(PetscDM,
18                                                                   PetscDM,
19                                                                   PetscMat*) except PETSC_ERR_PYTHON
20    ctypedef PetscErrorCode (*PetscDMShellCreateRestrictionFunction)(PetscDM,
21                                                                     PetscDM,
22                                                                     PetscMat*) except PETSC_ERR_PYTHON
23    ctypedef PetscErrorCode (*PetscDMShellCreateFieldDecompositionFunction)(PetscDM,
24                                                                            PetscInt*,
25                                                                            char***,
26                                                                            PetscIS**,
27                                                                            PetscDM**) except PETSC_ERR_PYTHON
28    ctypedef PetscErrorCode (*PetscDMShellCreateDomainDecompositionFunction)(PetscDM,
29                                                                             PetscInt*,
30                                                                             char***,
31                                                                             PetscIS**,
32                                                                             PetscIS**,
33                                                                             PetscDM**) except PETSC_ERR_PYTHON
34    ctypedef PetscErrorCode (*PetscDMShellCreateDomainDecompositionScattersFunction)(PetscDM,
35                                                                                     PetscInt,
36                                                                                     PetscDM*,
37                                                                                     PetscScatter**,
38                                                                                     PetscScatter**,
39                                                                                     PetscScatter**) except PETSC_ERR_PYTHON
40    ctypedef PetscErrorCode (*PetscDMShellCreateSubDM)(PetscDM,
41                                                       PetscInt,
42                                                       PetscInt[],
43                                                       PetscIS*,
44                                                       PetscDM*) except PETSC_ERR_PYTHON
45    PetscErrorCode DMShellCreate(MPI_Comm, PetscDM*)
46    PetscErrorCode DMShellSetMatrix(PetscDM, PetscMat)
47    PetscErrorCode DMShellSetGlobalVector(PetscDM, PetscVec)
48    PetscErrorCode DMShellSetLocalVector(PetscDM, PetscVec)
49    PetscErrorCode DMShellSetCreateGlobalVector(PetscDM, PetscDMShellCreateVectorFunction)
50    PetscErrorCode DMShellSetCreateLocalVector(PetscDM, PetscDMShellCreateVectorFunction)
51    PetscErrorCode DMShellSetGlobalToLocal(PetscDM, PetscDMShellXToYFunction, PetscDMShellXToYFunction)
52    PetscErrorCode DMShellSetGlobalToLocalVecScatter(PetscDM, PetscScatter)
53    PetscErrorCode DMShellSetLocalToGlobal(PetscDM, PetscDMShellXToYFunction, PetscDMShellXToYFunction)
54    PetscErrorCode DMShellSetLocalToGlobalVecScatter(PetscDM, PetscScatter)
55    PetscErrorCode DMShellSetLocalToLocal(PetscDM, PetscDMShellXToYFunction, PetscDMShellXToYFunction)
56    PetscErrorCode DMShellSetLocalToLocalVecScatter(PetscDM, PetscScatter)
57    PetscErrorCode DMShellSetCreateMatrix(PetscDM, PetscDMShellCreateMatrixFunction)
58    PetscErrorCode DMShellSetCoarsen(PetscDM, PetscDMShellTransferFunction)
59    PetscErrorCode DMShellSetRefine(PetscDM, PetscDMShellTransferFunction)
60    PetscErrorCode DMShellSetCreateInterpolation(PetscDM, PetscDMShellCreateInterpolationFunction)
61    PetscErrorCode DMShellSetCreateInjection(PetscDM, PetscDMShellCreateInjectionFunction)
62    PetscErrorCode DMShellSetCreateRestriction(PetscDM, PetscDMShellCreateRestrictionFunction)
63    PetscErrorCode DMShellSetCreateFieldDecomposition(PetscDM, PetscDMShellCreateFieldDecompositionFunction)
64    PetscErrorCode DMShellSetCreateDomainDecomposition(PetscDM, PetscDMShellCreateDomainDecompositionFunction)
65    PetscErrorCode DMShellSetCreateDomainDecompositionScatters(PetscDM, PetscDMShellCreateDomainDecompositionScattersFunction)
66    PetscErrorCode DMShellSetCreateSubDM(PetscDM, PetscDMShellCreateSubDM)
67
68    PetscErrorCode VecGetDM(PetscVec, PetscDM*)
69    PetscErrorCode VecSetDM(PetscVec, PetscDM)
70    PetscErrorCode MatGetDM(PetscMat, PetscDM*)
71    PetscErrorCode MatSetDM(PetscMat, PetscDM)
72
73
74cdef PetscErrorCode DMSHELL_CreateGlobalVector(
75    PetscDM dm,
76    PetscVec *v) except PETSC_ERR_PYTHON with gil:
77    cdef DM Dm = subtype_DM(dm)()
78    cdef Vec vec
79    Dm.dm = dm
80    CHKERR(PetscINCREF(Dm.obj))
81    context = Dm.get_attr('__create_global_vector__')
82    assert context is not None and type(context) is tuple
83    (create_gvec, args, kargs) = context
84    vec = create_gvec(Dm, *args, **kargs)
85    CHKERR(PetscINCREF(vec.obj))
86    v[0] = vec.vec
87    cdef PetscDM odm = NULL
88    CHKERR(VecGetDM(v[0], &odm))
89    if odm == NULL:
90        CHKERR(VecSetDM(v[0], dm))
91    return PETSC_SUCCESS
92
93cdef PetscErrorCode DMSHELL_CreateLocalVector(
94    PetscDM dm,
95    PetscVec *v) except PETSC_ERR_PYTHON with gil:
96    cdef DM Dm = subtype_DM(dm)()
97    cdef Vec vec
98    Dm.dm = dm
99    CHKERR(PetscINCREF(Dm.obj))
100    context = Dm.get_attr('__create_local_vector__')
101    assert context is not None and type(context) is tuple
102    (create_lvec, args, kargs) = context
103    vec = create_lvec(Dm, *args, **kargs)
104    CHKERR(PetscINCREF(vec.obj))
105    v[0] = vec.vec
106    cdef PetscDM odm = NULL
107    CHKERR(VecGetDM(v[0], &odm))
108    if odm == NULL:
109        CHKERR(VecSetDM(v[0], dm))
110    return PETSC_SUCCESS
111
112cdef PetscErrorCode DMSHELL_GlobalToLocalBegin(
113    PetscDM dm,
114    PetscVec g,
115    PetscInsertMode mode,
116    PetscVec l) except PETSC_ERR_PYTHON with gil:
117    cdef DM Dm = subtype_DM(dm)()
118    cdef Vec gvec = ref_Vec(g)
119    cdef Vec lvec = ref_Vec(l)
120    Dm.dm = dm
121    CHKERR(PetscINCREF(Dm.obj))
122    context = Dm.get_attr('__g2l_begin__')
123    assert context is not None and type(context) is tuple
124    (begin, args, kargs) = context
125    begin(Dm, gvec, mode, lvec, *args, **kargs)
126    return PETSC_SUCCESS
127
128cdef PetscErrorCode DMSHELL_GlobalToLocalEnd(
129    PetscDM dm,
130    PetscVec g,
131    PetscInsertMode mode,
132    PetscVec l) except PETSC_ERR_PYTHON with gil:
133    cdef DM Dm = subtype_DM(dm)()
134    cdef Vec gvec = ref_Vec(g)
135    cdef Vec lvec = ref_Vec(l)
136    Dm.dm = dm
137    CHKERR(PetscINCREF(Dm.obj))
138    context = Dm.get_attr('__g2l_end__')
139    assert context is not None and type(context) is tuple
140    (end, args, kargs) = context
141    end(Dm, gvec, mode, lvec, *args, **kargs)
142    return PETSC_SUCCESS
143
144cdef PetscErrorCode DMSHELL_LocalToGlobalBegin(
145    PetscDM dm,
146    PetscVec g,
147    PetscInsertMode mode,
148    PetscVec l) except PETSC_ERR_PYTHON with gil:
149    cdef DM Dm = subtype_DM(dm)()
150    cdef Vec gvec = ref_Vec(g)
151    cdef Vec lvec = ref_Vec(l)
152    Dm.dm = dm
153    CHKERR(PetscINCREF(Dm.obj))
154    context = Dm.get_attr('__l2g_begin__')
155    assert context is not None and type(context) is tuple
156    (begin, args, kargs) = context
157    begin(Dm, gvec, mode, lvec, *args, **kargs)
158    return PETSC_SUCCESS
159
160cdef PetscErrorCode DMSHELL_LocalToGlobalEnd(
161    PetscDM dm,
162    PetscVec g,
163    PetscInsertMode mode,
164    PetscVec l) except PETSC_ERR_PYTHON with gil:
165    cdef DM Dm = subtype_DM(dm)()
166    cdef Vec gvec = ref_Vec(g)
167    cdef Vec lvec = ref_Vec(l)
168    Dm.dm = dm
169    CHKERR(PetscINCREF(Dm.obj))
170    context = Dm.get_attr('__l2g_end__')
171    assert context is not None and type(context) is tuple
172    (end, args, kargs) = context
173    end(Dm, gvec, mode, lvec, *args, **kargs)
174    return PETSC_SUCCESS
175
176cdef PetscErrorCode DMSHELL_LocalToLocalBegin(
177    PetscDM dm,
178    PetscVec g,
179    PetscInsertMode mode,
180    PetscVec l) except PETSC_ERR_PYTHON with gil:
181    cdef DM Dm = subtype_DM(dm)()
182    cdef Vec gvec = ref_Vec(g)
183    cdef Vec lvec = ref_Vec(l)
184    Dm.dm = dm
185    CHKERR(PetscINCREF(Dm.obj))
186    context = Dm.get_attr('__l2l_begin__')
187    assert context is not None and type(context) is tuple
188    (begin, args, kargs) = context
189    begin(Dm, gvec, mode, lvec, *args, **kargs)
190    return PETSC_SUCCESS
191
192cdef PetscErrorCode DMSHELL_LocalToLocalEnd(
193    PetscDM dm,
194    PetscVec g,
195    PetscInsertMode mode,
196    PetscVec l) except PETSC_ERR_PYTHON with gil:
197    cdef DM Dm = subtype_DM(dm)()
198    cdef Vec gvec = ref_Vec(g)
199    cdef Vec lvec = ref_Vec(l)
200    Dm.dm = dm
201    CHKERR(PetscINCREF(Dm.obj))
202    context = Dm.get_attr('__l2l_end__')
203    assert context is not None and type(context) is tuple
204    (end, args, kargs) = context
205    end(Dm, gvec, mode, lvec, *args, **kargs)
206    return PETSC_SUCCESS
207
208cdef PetscErrorCode DMSHELL_CreateMatrix(
209    PetscDM dm,
210    PetscMat *cmat) except PETSC_ERR_PYTHON with gil:
211    cdef DM Dm = subtype_DM(dm)()
212    cdef Mat mat
213    Dm.dm = dm
214    CHKERR(PetscINCREF(Dm.obj))
215    context = Dm.get_attr('__create_matrix__')
216    assert context is not None and type(context) is tuple
217    (matrix, args, kargs) = context
218    mat = matrix(Dm, *args, **kargs)
219    CHKERR(PetscINCREF(mat.obj))
220    cmat[0] = mat.mat
221    cdef PetscDM odm = NULL
222    CHKERR(MatGetDM(cmat[0], &odm))
223    if odm == NULL:
224        CHKERR(MatSetDM(cmat[0], dm))
225    return PETSC_SUCCESS
226
227cdef PetscErrorCode DMSHELL_Coarsen(
228    PetscDM dm,
229    MPI_Comm comm,
230    PetscDM *dmc) except PETSC_ERR_PYTHON with gil:
231    cdef DM Dm = subtype_DM(dm)()
232    cdef DM Dmc
233    cdef Comm Comm = new_Comm(comm)
234    Dm.dm = dm
235    CHKERR(PetscINCREF(Dm.obj))
236    context = Dm.get_attr('__coarsen__')
237    assert context is not None and type(context) is tuple
238    (coarsen, args, kargs) = context
239    Dmc = coarsen(Dm, Comm, *args, **kargs)
240    CHKERR(PetscINCREF(Dmc.obj))
241    dmc[0] = Dmc.dm
242    return PETSC_SUCCESS
243
244cdef PetscErrorCode DMSHELL_Refine(
245    PetscDM dm,
246    MPI_Comm comm,
247    PetscDM *dmf) except PETSC_ERR_PYTHON with gil:
248    cdef DM Dm = subtype_DM(dm)()
249    cdef DM Dmf
250    cdef Comm Comm = new_Comm(comm)
251    Dm.dm = dm
252    CHKERR(PetscINCREF(Dm.obj))
253    context = Dm.get_attr('__refine__')
254    assert context is not None and type(context) is tuple
255    (refine, args, kargs) = context
256    Dmf = refine(Dm, Comm, *args, **kargs)
257    CHKERR(PetscINCREF(Dmf.obj))
258    dmf[0] = Dmf.dm
259    return PETSC_SUCCESS
260
261cdef PetscErrorCode DMSHELL_CreateInterpolation(
262    PetscDM dmc,
263    PetscDM dmf,
264    PetscMat *cmat,
265    PetscVec *cvec) except PETSC_ERR_PYTHON with gil:
266    cdef DM Dmc = subtype_DM(dmc)()
267    cdef DM Dmf = subtype_DM(dmf)()
268    cdef Mat mat
269    cdef Vec vec
270    Dmc.dm = dmc
271    CHKERR(PetscINCREF(Dmc.obj))
272    Dmf.dm = dmf
273    CHKERR(PetscINCREF(Dmf.obj))
274    context = Dmc.get_attr('__create_interpolation__')
275    assert context is not None and type(context) is tuple
276    (interpolation, args, kargs) = context
277    mat, vec = interpolation(Dmc, Dmf, *args, **kargs)
278    CHKERR(PetscINCREF(mat.obj))
279    cmat[0] = mat.mat
280    if cvec == NULL:
281        return PETSC_SUCCESS
282    if vec is None:
283        cvec[0] = NULL
284    else:
285        CHKERR(PetscINCREF(vec.obj))
286        cvec[0] = vec.vec
287    return PETSC_SUCCESS
288
289cdef PetscErrorCode DMSHELL_CreateInjection(
290    PetscDM dmc,
291    PetscDM dmf,
292    PetscMat *cmat) except PETSC_ERR_PYTHON with gil:
293    cdef DM Dmc = subtype_DM(dmc)()
294    cdef DM Dmf = subtype_DM(dmf)()
295    cdef Mat mat
296    Dmc.dm = dmc
297    CHKERR(PetscINCREF(Dmc.obj))
298    Dmf.dm = dmf
299    CHKERR(PetscINCREF(Dmf.obj))
300    context = Dmc.get_attr('__create_injection__')
301    assert context is not None and type(context) is tuple
302    (injection, args, kargs) = context
303    mat = injection(Dmc, Dmf, *args, **kargs)
304    CHKERR(PetscINCREF(mat.obj))
305    cmat[0] = mat.mat
306    return PETSC_SUCCESS
307
308cdef PetscErrorCode DMSHELL_CreateRestriction(
309    PetscDM dmf,
310    PetscDM dmc,
311    PetscMat *cmat) except PETSC_ERR_PYTHON with gil:
312    cdef DM Dmf = subtype_DM(dmf)()
313    cdef DM Dmc = subtype_DM(dmc)()
314    cdef Mat mat
315    Dmf.dm = dmf
316    CHKERR(PetscINCREF(Dmf.obj))
317    Dmc.dm = dmc
318    CHKERR(PetscINCREF(Dmc.obj))
319    context = Dmf.get_attr('__create_restriction__')
320    assert context is not None and type(context) is tuple
321    (restriction, args, kargs) = context
322    mat = restriction(Dmf, Dmc, *args, **kargs)
323    CHKERR(PetscINCREF(mat.obj))
324    cmat[0] = mat.mat
325    return PETSC_SUCCESS
326
327cdef PetscErrorCode DMSHELL_CreateFieldDecomposition(
328    PetscDM dm,
329    PetscInt *clen,
330    char ***namelist,
331    PetscIS **islist,
332    PetscDM **dmlist) except PETSC_ERR_PYTHON with gil:
333    cdef DM Dm = subtype_DM(dm)()
334    cdef int i
335    cdef const char *cname = NULL
336    Dm.dm = dm
337    CHKERR(PetscINCREF(Dm.obj))
338    context = Dm.get_attr('__create_field_decomp__')
339    assert context is not None and type(context) is tuple
340    (decomp, args, kargs) = context
341    names, ises, dms = decomp(Dm, *args, **kargs)
342
343    if clen != NULL:
344        if names is not None:
345            clen[0] = <PetscInt>len(names)
346        elif ises is not None:
347            clen[0] = <PetscInt>len(ises)
348        elif dms is not None:
349            clen[0] = <PetscInt>len(dms)
350        else:
351            clen[0] = 0
352
353    if namelist != NULL and names is not None:
354        CHKERR(PetscMalloc(len(names)*sizeof(char**), namelist))
355        for i in range(len(names)):
356            names[i] = str2bytes(names[i], &cname)
357            CHKERR(PetscStrallocpy(cname, &namelist[0][i]))
358
359    if islist != NULL and ises is not None:
360        CHKERR(PetscMalloc(len(ises)*sizeof(PetscIS), islist))
361        for i in range(len(ises)):
362            islist[0][i] = (<IS?>ises[i]).iset
363            CHKERR(PetscINCREF((<IS?>ises[i]).obj))
364
365    if dmlist != NULL and dms is not None:
366        CHKERR(PetscMalloc(len(dms)*sizeof(PetscDM), dmlist))
367        for i in range(len(dms)):
368            dmlist[0][i] = (<DM?>dms[i]).dm
369            CHKERR(PetscINCREF((<DM?>dms[i]).obj))
370    return PETSC_SUCCESS
371
372cdef PetscErrorCode DMSHELL_CreateDomainDecomposition(
373    PetscDM dm,
374    PetscInt *clen,
375    char ***namelist,
376    PetscIS **innerislist,
377    PetscIS **outerislist,
378    PetscDM **dmlist) except PETSC_ERR_PYTHON with gil:
379    cdef DM Dm = subtype_DM(dm)()
380    cdef int i
381    cdef const char *cname = NULL
382    Dm.dm = dm
383    CHKERR(PetscINCREF(Dm.obj))
384    context = Dm.get_attr('__create_domain_decomp__')
385    assert context is not None and type(context) is tuple
386    (decomp, args, kargs) = context
387    names, innerises, outerises, dms = decomp(Dm, *args, **kargs)
388
389    if clen != NULL:
390        if names is not None:
391            clen[0] = <PetscInt>len(names)
392        elif innerises is not None:
393            clen[0] = <PetscInt>len(innerises)
394        elif outerises is not None:
395            clen[0] = <PetscInt>len(outerises)
396        elif dms is not None:
397            clen[0] = <PetscInt>len(dms)
398        else:
399            clen[0] = 0
400
401    if namelist != NULL and names is not None:
402        CHKERR(PetscMalloc(len(names)*sizeof(char**), namelist))
403        for i in range(len(names)):
404            names[i] = str2bytes(names[i], &cname)
405            CHKERR(PetscStrallocpy(cname, &namelist[0][i]))
406
407    if innerislist != NULL and innerises is not None:
408        CHKERR(PetscMalloc(len(innerises)*sizeof(PetscIS), innerislist))
409        for i in range(len(innerises)):
410            innerislist[0][i] = (<IS?>innerises[i]).iset
411            CHKERR(PetscINCREF((<IS?>innerises[i]).obj))
412
413    if outerislist != NULL and outerises is not None:
414        CHKERR(PetscMalloc(len(outerises)*sizeof(PetscIS), outerislist))
415        for i in range(len(outerises)):
416            outerislist[0][i] = (<IS?>outerises[i]).iset
417            CHKERR(PetscINCREF((<IS?>outerises[i]).obj))
418
419    if dmlist != NULL and dms is not None:
420        CHKERR(PetscMalloc(len(dms)*sizeof(PetscDM), dmlist))
421        for i in range(len(dms)):
422            dmlist[0][i] = (<DM?>dms[i]).dm
423            CHKERR(PetscINCREF((<DM?>dms[i]).obj))
424    return PETSC_SUCCESS
425
426cdef PetscErrorCode DMSHELL_CreateDomainDecompositionScatters(
427    PetscDM dm,
428    PetscInt clen,
429    PetscDM *subdms,
430    PetscScatter** iscat,
431    PetscScatter** oscat,
432    PetscScatter** gscat) except PETSC_ERR_PYTHON with gil:
433
434    cdef DM Dm = subtype_DM(dm)()
435    cdef int i
436    cdef DM subdm = None
437
438    Dm.dm = dm
439    CHKERR(PetscINCREF(Dm.obj))
440
441    psubdms = []
442    for i from 0 <= i < clen:
443        subdm = subtype_DM(subdms[i])()
444        subdm.dm = subdms[i]
445        CHKERR(PetscINCREF(subdm.obj))
446        psubdms.append(subdm)
447
448    context = Dm.get_attr('__create_domain_decomp_scatters__')
449    assert context is not None and type(context) is tuple
450    (scatters, args, kargs) = context
451    (iscatter, oscatter, gscatter) = scatters(Dm, psubdms, *args, **kargs)
452
453    assert len(iscatter) == clen
454    assert len(oscatter) == clen
455    assert len(gscatter) == clen
456
457    CHKERR (PetscMalloc(clen*sizeof(PetscScatter), iscat))
458    CHKERR (PetscMalloc(clen*sizeof(PetscScatter), oscat))
459    CHKERR (PetscMalloc(clen*sizeof(PetscScatter), gscat))
460
461    for i in range(clen):
462        iscat[0][i] = (<Scatter?>iscatter[i]).sct
463        CHKERR(PetscINCREF((<Scatter?>iscatter[i]).obj))
464
465        oscat[0][i] = (<Scatter?>oscatter[i]).sct
466        CHKERR(PetscINCREF((<Scatter?>oscatter[i]).obj))
467
468        gscat[0][i] = (<Scatter?>gscatter[i]).sct
469        CHKERR(PetscINCREF((<Scatter?>gscatter[i]).obj))
470
471    return PETSC_SUCCESS
472
473cdef PetscErrorCode DMSHELL_CreateSubDM(
474    PetscDM cdm,
475    PetscInt numFields,
476    const PetscInt cfields[],
477    PetscIS *ciset,
478    PetscDM *csubdm) except PETSC_ERR_PYTHON with gil:
479    cdef DM dm = subtype_DM(cdm)()
480    cdef IS iset
481    cdef DM subdm
482    dm.dm = cdm
483    CHKERR(PetscINCREF(dm.obj))
484    context = dm.get_attr('__create_subdm__')
485    assert context is not None and type(context) is tuple
486    (create_subdm, args, kargs) = context
487
488    fields = array_i(numFields, cfields)
489
490    iset, subdm = create_subdm(dm, fields, *args, **kargs)
491
492    CHKERR(PetscINCREF(iset.obj))
493    CHKERR(PetscINCREF(subdm.obj))
494    ciset[0] = iset.iset
495    csubdm[0] = subdm.dm
496    return PETSC_SUCCESS
497