16dd63270SBarry Smith #include <petsc/private/ftnimpl.h>
20ec63f53SRichard Tran Mills #include <petscdmshell.h> /*I "petscdmshell.h" I*/
30ec63f53SRichard Tran Mills
40ec63f53SRichard Tran Mills #if defined(PETSC_HAVE_FORTRAN_CAPS)
50ec63f53SRichard Tran Mills #define dmshellsetcreatematrix_ DMSHELLSETCREATEMATRIX
68a6b6cadSSatish Balay #define dmshellsetcreateglobalvector_ DMSHELLSETCREATEGLOBALVECTOR
78a6b6cadSSatish Balay #define dmshellsetcreatelocalvector_ DMSHELLSETCREATELOCALVECTOR
88a6b6cadSSatish Balay #define dmshellsetglobaltolocal_ DMSHELLSETGLOBALTOLOCAL
98a6b6cadSSatish Balay #define dmshellsetlocaltoglobal_ DMSHELLSETLOCALTOGLOBAL
108a6b6cadSSatish Balay #define dmshellsetlocaltolocal_ DMSHELLSETLOCALTOLOCAL
11*62212064STapashree Pradhan #define dmshellsetcreatefielddecomposition_ DMSHELLSETCREATEFIELDDECOMPOSITION
120ec63f53SRichard Tran Mills #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
130ec63f53SRichard Tran Mills #define dmshellsetcreatematrix_ dmshellsetcreatematrix
140ec63f53SRichard Tran Mills #define dmshellsetcreateglobalvector_ dmshellsetcreateglobalvector
15dc43b69eSJed Brown #define dmshellsetcreatelocalvector_ dmshellsetcreatelocalvector
16f2d7aa87SRichard Tran Mills #define dmshellsetglobaltolocal_ dmshellsetglobaltolocal
17f2d7aa87SRichard Tran Mills #define dmshellsetlocaltoglobal_ dmshellsetlocaltoglobal
188a6b6cadSSatish Balay #define dmshellsetlocaltolocal_ dmshellsetlocaltolocal
19*62212064STapashree Pradhan #define dmshellsetcreatefielddecomposition_ dmshellsetcreatefielddecomposition
200ec63f53SRichard Tran Mills #endif
210ec63f53SRichard Tran Mills
220ec63f53SRichard Tran Mills /*
23dc43b69eSJed Brown * C routines are required for matrix and global vector creation. We define C routines here that call the corresponding
24dc43b69eSJed Brown * Fortran routine (indexed by _cb) that was set by the user.
250ec63f53SRichard Tran Mills */
260ec63f53SRichard Tran Mills
27dc43b69eSJed Brown static struct {
28dc43b69eSJed Brown PetscFortranCallbackId creatematrix;
29dc43b69eSJed Brown PetscFortranCallbackId createglobalvector;
30dc43b69eSJed Brown PetscFortranCallbackId createlocalvector;
31f2d7aa87SRichard Tran Mills PetscFortranCallbackId globaltolocalbegin;
32f2d7aa87SRichard Tran Mills PetscFortranCallbackId globaltolocalend;
33f2d7aa87SRichard Tran Mills PetscFortranCallbackId localtoglobalbegin;
34f2d7aa87SRichard Tran Mills PetscFortranCallbackId localtoglobalend;
35f3db62a7SRichard Tran Mills PetscFortranCallbackId localtolocalbegin;
36f3db62a7SRichard Tran Mills PetscFortranCallbackId localtolocalend;
37*62212064STapashree Pradhan PetscFortranCallbackId createfielddecomposition;
38dc43b69eSJed Brown } _cb;
39dc43b69eSJed Brown
ourcreatematrix(DM dm,Mat * A)40b412c318SBarry Smith static PetscErrorCode ourcreatematrix(DM dm, Mat *A)
410ec63f53SRichard Tran Mills {
42a348f287SBarry Smith PetscObjectUseFortranCallbackSubType(dm, _cb.creatematrix, (DM *, Mat *, PetscErrorCode *), (&dm, A, &ierr));
430ec63f53SRichard Tran Mills }
440ec63f53SRichard Tran Mills
ourcreateglobalvector(DM dm,Vec * v)450ec63f53SRichard Tran Mills static PetscErrorCode ourcreateglobalvector(DM dm, Vec *v)
460ec63f53SRichard Tran Mills {
4717da0f0dSJed Brown PetscObjectUseFortranCallbackSubType(dm, _cb.createglobalvector, (DM *, Vec *, PetscErrorCode *), (&dm, v, &ierr));
480ec63f53SRichard Tran Mills }
490ec63f53SRichard Tran Mills
ourcreatelocalvector(DM dm,Vec * v)50dc43b69eSJed Brown static PetscErrorCode ourcreatelocalvector(DM dm, Vec *v)
510ec63f53SRichard Tran Mills {
5217da0f0dSJed Brown PetscObjectUseFortranCallbackSubType(dm, _cb.createlocalvector, (DM *, Vec *, PetscErrorCode *), (&dm, v, &ierr));
53dc43b69eSJed Brown }
54dc43b69eSJed Brown
ourglobaltolocalbegin(DM dm,Vec g,InsertMode mode,Vec l)55f2d7aa87SRichard Tran Mills static PetscErrorCode ourglobaltolocalbegin(DM dm, Vec g, InsertMode mode, Vec l)
56f2d7aa87SRichard Tran Mills {
5717da0f0dSJed Brown PetscObjectUseFortranCallbackSubType(dm, _cb.globaltolocalbegin, (DM *, Vec *, InsertMode *, Vec *, PetscErrorCode *), (&dm, &g, &mode, &l, &ierr));
58f2d7aa87SRichard Tran Mills }
59f2d7aa87SRichard Tran Mills
ourglobaltolocalend(DM dm,Vec g,InsertMode mode,Vec l)60f2d7aa87SRichard Tran Mills static PetscErrorCode ourglobaltolocalend(DM dm, Vec g, InsertMode mode, Vec l)
61f2d7aa87SRichard Tran Mills {
6217da0f0dSJed Brown PetscObjectUseFortranCallbackSubType(dm, _cb.globaltolocalend, (DM *, Vec *, InsertMode *, Vec *, PetscErrorCode *), (&dm, &g, &mode, &l, &ierr));
63f2d7aa87SRichard Tran Mills }
64f2d7aa87SRichard Tran Mills
ourlocaltoglobalbegin(DM dm,Vec l,InsertMode mode,Vec g)65f2d7aa87SRichard Tran Mills static PetscErrorCode ourlocaltoglobalbegin(DM dm, Vec l, InsertMode mode, Vec g)
66f2d7aa87SRichard Tran Mills {
6717da0f0dSJed Brown PetscObjectUseFortranCallbackSubType(dm, _cb.localtoglobalbegin, (DM *, Vec *, InsertMode *, Vec *, PetscErrorCode *), (&dm, &l, &mode, &g, &ierr));
68f2d7aa87SRichard Tran Mills }
69f2d7aa87SRichard Tran Mills
ourlocaltoglobalend(DM dm,Vec l,InsertMode mode,Vec g)70f2d7aa87SRichard Tran Mills static PetscErrorCode ourlocaltoglobalend(DM dm, Vec l, InsertMode mode, Vec g)
71f2d7aa87SRichard Tran Mills {
7217da0f0dSJed Brown PetscObjectUseFortranCallbackSubType(dm, _cb.localtoglobalend, (DM *, Vec *, InsertMode *, Vec *, PetscErrorCode *), (&dm, &l, &mode, &g, &ierr));
73f2d7aa87SRichard Tran Mills }
74f2d7aa87SRichard Tran Mills
ourlocaltolocalbegin(DM dm,Vec g,InsertMode mode,Vec l)75f3db62a7SRichard Tran Mills static PetscErrorCode ourlocaltolocalbegin(DM dm, Vec g, InsertMode mode, Vec l)
76f3db62a7SRichard Tran Mills {
77f3db62a7SRichard Tran Mills PetscObjectUseFortranCallbackSubType(dm, _cb.localtolocalbegin, (DM *, Vec *, InsertMode *, Vec *, PetscErrorCode *), (&dm, &g, &mode, &l, &ierr));
78f3db62a7SRichard Tran Mills }
79f3db62a7SRichard Tran Mills
ourlocaltolocalend(DM dm,Vec g,InsertMode mode,Vec l)80f3db62a7SRichard Tran Mills static PetscErrorCode ourlocaltolocalend(DM dm, Vec g, InsertMode mode, Vec l)
81f3db62a7SRichard Tran Mills {
82f3db62a7SRichard Tran Mills PetscObjectUseFortranCallbackSubType(dm, _cb.localtolocalend, (DM *, Vec *, InsertMode *, Vec *, PetscErrorCode *), (&dm, &g, &mode, &l, &ierr));
83f3db62a7SRichard Tran Mills }
84f3db62a7SRichard Tran Mills
ourcreatefielddecomposition(DM dm,PetscInt * nfields,char *** names,IS ** is,DM ** subdms)85*62212064STapashree Pradhan static PetscErrorCode ourcreatefielddecomposition(DM dm, PetscInt *nfields, char ***names, IS **is, DM **subdms)
86*62212064STapashree Pradhan {
87*62212064STapashree Pradhan PetscObjectUseFortranCallbackSubType(dm, _cb.createfielddecomposition, (DM *, PetscInt *, char ***, IS **, DM **, PetscErrorCode *), (&dm, nfields, names, is, subdms, &ierr));
88*62212064STapashree Pradhan }
89*62212064STapashree Pradhan
dmshellsetcreatematrix_(DM * dm,void (* func)(DM *,Mat *,PetscErrorCode *,PETSC_FORTRAN_CHARLEN_T len),PetscErrorCode * ierr)9019caf8f3SSatish Balay PETSC_EXTERN void dmshellsetcreatematrix_(DM *dm, void (*func)(DM *, Mat *, PetscErrorCode *, PETSC_FORTRAN_CHARLEN_T len), PetscErrorCode *ierr)
91dc43b69eSJed Brown {
928434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*dm, PETSC_FORTRAN_CALLBACK_SUBTYPE, &_cb.creatematrix, (PetscVoidFn *)func, NULL);
935975b3b6SBarry Smith if (*ierr) return;
940ec63f53SRichard Tran Mills *ierr = DMShellSetCreateMatrix(*dm, ourcreatematrix);
950ec63f53SRichard Tran Mills }
960ec63f53SRichard Tran Mills
dmshellsetcreateglobalvector_(DM * dm,void (* func)(DM *,Vec *,PetscErrorCode *),PetscErrorCode * ierr)9719caf8f3SSatish Balay PETSC_EXTERN void dmshellsetcreateglobalvector_(DM *dm, void (*func)(DM *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
980ec63f53SRichard Tran Mills {
998434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*dm, PETSC_FORTRAN_CALLBACK_SUBTYPE, &_cb.createglobalvector, (PetscVoidFn *)func, NULL);
1005975b3b6SBarry Smith if (*ierr) return;
1010ec63f53SRichard Tran Mills *ierr = DMShellSetCreateGlobalVector(*dm, ourcreateglobalvector);
1020ec63f53SRichard Tran Mills }
103fa59e805SSatish Balay
dmshellsetcreatelocalvector_(DM * dm,void (* func)(DM *,Vec *,PetscErrorCode *),PetscErrorCode * ierr)10419caf8f3SSatish Balay PETSC_EXTERN void dmshellsetcreatelocalvector_(DM *dm, void (*func)(DM *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
105dc43b69eSJed Brown {
1068434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*dm, PETSC_FORTRAN_CALLBACK_SUBTYPE, &_cb.createlocalvector, (PetscVoidFn *)func, NULL);
1075975b3b6SBarry Smith if (*ierr) return;
108dc43b69eSJed Brown *ierr = DMShellSetCreateLocalVector(*dm, ourcreatelocalvector);
109dc43b69eSJed Brown }
110f2d7aa87SRichard Tran Mills
dmshellsetglobaltolocal_(DM * dm,void (* begin)(DM *,Vec *,InsertMode *,Vec *,PetscErrorCode *),void (* end)(DM *,Vec *,InsertMode *,Vec *,PetscErrorCode *),PetscErrorCode * ierr)11119caf8f3SSatish Balay PETSC_EXTERN void dmshellsetglobaltolocal_(DM *dm, void (*begin)(DM *, Vec *, InsertMode *, Vec *, PetscErrorCode *), void (*end)(DM *, Vec *, InsertMode *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
112f2d7aa87SRichard Tran Mills {
1138434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*dm, PETSC_FORTRAN_CALLBACK_SUBTYPE, &_cb.globaltolocalbegin, (PetscVoidFn *)begin, NULL);
1145975b3b6SBarry Smith if (*ierr) return;
1158434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*dm, PETSC_FORTRAN_CALLBACK_SUBTYPE, &_cb.globaltolocalend, (PetscVoidFn *)end, NULL);
1165975b3b6SBarry Smith if (*ierr) return;
117f2d7aa87SRichard Tran Mills *ierr = DMShellSetGlobalToLocal(*dm, ourglobaltolocalbegin, ourglobaltolocalend);
118f2d7aa87SRichard Tran Mills }
119f2d7aa87SRichard Tran Mills
dmshellsetlocaltoglobal_(DM * dm,void (* begin)(DM *,Vec *,InsertMode *,Vec *,PetscErrorCode *),void (* end)(DM *,Vec *,InsertMode *,Vec *,PetscErrorCode *),PetscErrorCode * ierr)12019caf8f3SSatish Balay PETSC_EXTERN void dmshellsetlocaltoglobal_(DM *dm, void (*begin)(DM *, Vec *, InsertMode *, Vec *, PetscErrorCode *), void (*end)(DM *, Vec *, InsertMode *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
121f2d7aa87SRichard Tran Mills {
1228434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*dm, PETSC_FORTRAN_CALLBACK_SUBTYPE, &_cb.localtoglobalbegin, (PetscVoidFn *)begin, NULL);
1235975b3b6SBarry Smith if (*ierr) return;
1248434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*dm, PETSC_FORTRAN_CALLBACK_SUBTYPE, &_cb.localtoglobalend, (PetscVoidFn *)end, NULL);
1255975b3b6SBarry Smith if (*ierr) return;
12617da0f0dSJed Brown *ierr = DMShellSetLocalToGlobal(*dm, ourlocaltoglobalbegin, ourlocaltoglobalend);
127f2d7aa87SRichard Tran Mills }
128f3db62a7SRichard Tran Mills
dmshellsetlocaltolocal_(DM * dm,void (* begin)(DM *,Vec *,InsertMode *,Vec *,PetscErrorCode *),void (* end)(DM *,Vec *,InsertMode *,Vec *,PetscErrorCode *),PetscErrorCode * ierr)12919caf8f3SSatish Balay PETSC_EXTERN void dmshellsetlocaltolocal_(DM *dm, void (*begin)(DM *, Vec *, InsertMode *, Vec *, PetscErrorCode *), void (*end)(DM *, Vec *, InsertMode *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
130f3db62a7SRichard Tran Mills {
1318434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*dm, PETSC_FORTRAN_CALLBACK_SUBTYPE, &_cb.localtolocalbegin, (PetscVoidFn *)begin, NULL);
1325975b3b6SBarry Smith if (*ierr) return;
1338434afd1SBarry Smith *ierr = PetscObjectSetFortranCallback((PetscObject)*dm, PETSC_FORTRAN_CALLBACK_SUBTYPE, &_cb.localtolocalend, (PetscVoidFn *)end, NULL);
1345975b3b6SBarry Smith if (*ierr) return;
135f3db62a7SRichard Tran Mills *ierr = DMShellSetLocalToLocal(*dm, ourlocaltolocalbegin, ourlocaltolocalend);
136f3db62a7SRichard Tran Mills }
137*62212064STapashree Pradhan
dmshellsetcreatefielddecomposition_(DM * dm,void (* func)(DM *,PetscInt *,char ***,IS **,DM **,PetscErrorCode *),PetscErrorCode * ierr)138*62212064STapashree Pradhan PETSC_EXTERN void dmshellsetcreatefielddecomposition_(DM *dm, void (*func)(DM *, PetscInt *, char ***, IS **, DM **, PetscErrorCode *), PetscErrorCode *ierr)
139*62212064STapashree Pradhan {
140*62212064STapashree Pradhan *ierr = PetscObjectSetFortranCallback((PetscObject)*dm, PETSC_FORTRAN_CALLBACK_SUBTYPE, &_cb.createfielddecomposition, (PetscVoidFn *)func, NULL);
141*62212064STapashree Pradhan if (*ierr) return;
142*62212064STapashree Pradhan *ierr = DMShellSetCreateFieldDecomposition(*dm, ourcreatefielddecomposition);
143*62212064STapashree Pradhan }
144