xref: /petsc/src/ksp/ksp/interface/ftn-custom/zitfuncf.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1 #include <petsc/private/ftnimpl.h>
2 #include <petscksp.h>
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5   #define kspmonitorset_                  KSPMONITORSET
6   #define kspconvergeddefaultcreate_      KSPCONVERGEDDEFAULTCREATE
7   #define kspconvergeddefaultdestroycptr_ KSPCONVERGEDDEFAULTDESTROYCPTR
8   #define kspsetconvergencetest_          KSPSETCONVERGENCETEST
9   #define kspconvergeddefault_            KSPCONVERGEDDEFAULT
10   #define kspconvergedskip_               KSPCONVERGEDSKIP
11   #define kspgmresmonitorkrylov_          KSPGMRESMONITORKRYLOV
12   #define kspmonitorresidual_             KSPMONITORRESIDUAL
13   #define kspmonitortrueresidual_         KSPMONITORTRUERESIDUAL
14   #define kspmonitorsolution_             KSPMONITORSOLUTION
15   #define kspmonitorsingularvalue_        KSPMONITORSINGULARVALUE
16   #define kspsetcomputerhs_               KSPSETCOMPUTERHS
17   #define kspsetcomputeinitialguess_      KSPSETCOMPUTEINITIALGUESS
18   #define kspsetcomputeoperators_         KSPSETCOMPUTEOPERATORS
19   #define dmkspsetcomputerhs_             DMKSPSETCOMPUTERHS
20   #define dmkspsetcomputeinitialguess_    DMKSPSETCOMPUTEINITIALGUESS
21   #define dmkspsetcomputeoperators_       DMKSPSETCOMPUTEOPERATORS
22 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
23   #define kspmonitorset_                  kspmonitorset
24   #define kspconvergeddefaultcreate_      kspconvergeddefaultcreate
25   #define kspconvergeddefaultdestroycptr_ kspconvergeddefaultdestroycptr
26   #define kspsetconvergencetest_          kspsetconvergencetest
27   #define kspconvergeddefault_            kspconvergeddefault
28   #define kspconvergedskip_               kspconvergedskip
29   #define kspgmresmonitorkrylov_          kspgmresmonitorkrylov
30   #define kspmonitorresidual_             kspmonitorresidual
31   #define kspmonitortrueresidual_         kspmonitortrueresidual
32   #define kspmonitorsolution_             kspmonitorsolution
33   #define kspmonitorsingularvalue_        kspmonitorsingularvalue
34   #define kspsetcomputerhs_               kspsetcomputerhs
35   #define kspsetcomputeinitialguess_      kspsetcomputeinitialguess
36   #define kspsetcomputeoperators_         kspsetcomputeoperators
37   #define dmkspsetcomputerhs_             dmkspsetcomputerhs
38   #define dmkspsetcomputeinitialguess_    dmkspsetcomputeinitialguess
39   #define dmkspsetcomputeoperators_       dmkspsetcomputeoperators
40 #endif
41 
42 /* These are defined in zdmkspf.c */
43 PETSC_EXTERN void dmkspsetcomputerhs_(DM *dm, void (*func)(KSP *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr);
44 PETSC_EXTERN void dmkspsetcomputeinitialguess_(DM *dm, void (*func)(KSP *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr);
45 PETSC_EXTERN void dmkspsetcomputeoperators_(DM *dm, void (*func)(KSP *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr);
46 
47 /*
48         These cannot be called from Fortran but allow Fortran users to transparently set these monitors from .F code
49 */
50 
51 PETSC_EXTERN void kspconvergeddefault_(KSP *, PetscInt *, PetscReal *, KSPConvergedReason *, void *, PetscErrorCode *);
52 PETSC_EXTERN void kspconvergedskip_(KSP *, PetscInt *, PetscReal *, KSPConvergedReason *, void *, PetscErrorCode *);
53 PETSC_EXTERN void kspgmresmonitorkrylov_(KSP *, PetscInt *, PetscReal *, PetscViewerAndFormat *, PetscErrorCode *);
54 PETSC_EXTERN void kspmonitorresidual_(KSP *, PetscInt *, PetscReal *, PetscViewerAndFormat *, PetscErrorCode *);
55 PETSC_EXTERN void kspmonitorsingularvalue_(KSP *, PetscInt *, PetscReal *, PetscViewerAndFormat *, PetscErrorCode *);
56 PETSC_EXTERN void kspmonitortrueresidual_(KSP *, PetscInt *, PetscReal *, PetscViewerAndFormat *, PetscErrorCode *);
57 PETSC_EXTERN void kspmonitorsolution_(KSP *, PetscInt *, PetscReal *, PetscViewerAndFormat *, PetscErrorCode *);
58 
59 static struct {
60   PetscFortranCallbackId monitor;
61   PetscFortranCallbackId monitordestroy;
62   PetscFortranCallbackId test;
63   PetscFortranCallbackId testdestroy;
64 } _cb;
65 
ourmonitor(KSP ksp,PetscInt i,PetscReal d,PetscCtx ctx)66 static PetscErrorCode ourmonitor(KSP ksp, PetscInt i, PetscReal d, PetscCtx ctx)
67 {
68   PetscObjectUseFortranCallback(ksp, _cb.monitor, (KSP *, PetscInt *, PetscReal *, void *, PetscErrorCode *), (&ksp, &i, &d, _ctx, &ierr));
69 }
70 
ourdestroy(PetscCtxRt ctx)71 static PetscErrorCode ourdestroy(PetscCtxRt ctx)
72 {
73   KSP ksp = *(KSP *)ctx;
74   PetscObjectUseFortranCallback(ksp, _cb.monitordestroy, (void *, PetscErrorCode *), (_ctx, &ierr));
75 }
76 
77 /* These are not extern C because they are passed into non-extern C user level functions */
ourtest(KSP ksp,PetscInt i,PetscReal d,KSPConvergedReason * reason,PetscCtx ctx)78 static PetscErrorCode ourtest(KSP ksp, PetscInt i, PetscReal d, KSPConvergedReason *reason, PetscCtx ctx)
79 {
80   PetscObjectUseFortranCallback(ksp, _cb.test, (KSP *, PetscInt *, PetscReal *, KSPConvergedReason *, void *, PetscErrorCode *), (&ksp, &i, &d, reason, _ctx, &ierr));
81 }
82 
ourtestdestroy(PetscCtxRt ctx)83 static PetscErrorCode ourtestdestroy(PetscCtxRt ctx)
84 {
85   KSP ksp = *(KSP *)ctx;
86   PetscObjectUseFortranCallback(ksp, _cb.testdestroy, (void **, PetscErrorCode *), (&_ctx, &ierr));
87 }
88 
89 /*
90    For the built in monitors we ignore the monitordestroy that is passed in and use PetscViewerAndFormatDestroy()
91 */
kspmonitorset_(KSP * ksp,void (* monitor)(KSP *,PetscInt *,PetscReal *,void *,PetscErrorCode *),void * mctx,void (* monitordestroy)(void *,PetscErrorCode *),PetscErrorCode * ierr)92 PETSC_EXTERN void kspmonitorset_(KSP *ksp, void (*monitor)(KSP *, PetscInt *, PetscReal *, void *, PetscErrorCode *), void *mctx, void (*monitordestroy)(void *, PetscErrorCode *), PetscErrorCode *ierr)
93 {
94   CHKFORTRANNULLFUNCTION(monitordestroy);
95 
96   if ((PetscFortranCallbackFn *)monitor == (PetscFortranCallbackFn *)kspmonitorresidual_) {
97     *ierr = KSPMonitorSet(*ksp, (KSPMonitorFn *)KSPMonitorResidual, *(PetscViewerAndFormat **)mctx, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy);
98   } else if ((PetscFortranCallbackFn *)monitor == (PetscFortranCallbackFn *)kspmonitorsolution_) {
99     *ierr = KSPMonitorSet(*ksp, (KSPMonitorFn *)KSPMonitorSolution, *(PetscViewerAndFormat **)mctx, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy);
100   } else if ((PetscFortranCallbackFn *)monitor == (PetscFortranCallbackFn *)kspmonitortrueresidual_) {
101     *ierr = KSPMonitorSet(*ksp, (KSPMonitorFn *)KSPMonitorTrueResidual, *(PetscViewerAndFormat **)mctx, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy);
102   } else if ((PetscFortranCallbackFn *)monitor == (PetscFortranCallbackFn *)kspmonitorsingularvalue_) {
103     *ierr = KSPMonitorSet(*ksp, (KSPMonitorFn *)KSPMonitorSingularValue, *(PetscViewerAndFormat **)mctx, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy);
104   } else if ((PetscFortranCallbackFn *)monitor == (PetscFortranCallbackFn *)kspgmresmonitorkrylov_) {
105     *ierr = KSPMonitorSet(*ksp, (KSPMonitorFn *)KSPGMRESMonitorKrylov, *(PetscViewerAndFormat **)mctx, (PetscCtxDestroyFn *)PetscViewerAndFormatDestroy);
106   } else {
107     *ierr = PetscObjectSetFortranCallback((PetscObject)*ksp, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.monitor, (PetscFortranCallbackFn *)monitor, mctx);
108     if (*ierr) return;
109     *ierr = PetscObjectSetFortranCallback((PetscObject)*ksp, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.monitordestroy, (PetscFortranCallbackFn *)monitordestroy, mctx);
110     if (*ierr) return;
111     *ierr = KSPMonitorSet(*ksp, ourmonitor, *ksp, ourdestroy);
112   }
113 }
114 
115 PETSC_EXTERN void kspconvergeddefaultdestroycptr_(PetscCtxRt, PetscErrorCode *);
116 
kspsetconvergencetest_(KSP * ksp,void (* converge)(KSP *,PetscInt *,PetscReal *,KSPConvergedReason *,void *,PetscErrorCode *),void ** cctx,void (* destroy)(PetscCtxRt,PetscErrorCode *),PetscErrorCode * ierr)117 PETSC_EXTERN void kspsetconvergencetest_(KSP *ksp, void (*converge)(KSP *, PetscInt *, PetscReal *, KSPConvergedReason *, void *, PetscErrorCode *), void **cctx, void (*destroy)(PetscCtxRt, PetscErrorCode *), PetscErrorCode *ierr)
118 {
119   CHKFORTRANNULLFUNCTION(destroy);
120 
121   if (converge == kspconvergeddefault_) {
122     *ierr = KSPSetConvergenceTest(*ksp, KSPConvergedDefault, &cctx, KSPConvergedDefaultDestroy);
123   } else if (converge == kspconvergedskip_) {
124     *ierr = KSPSetConvergenceTest(*ksp, KSPConvergedSkip, NULL, NULL);
125   } else {
126     if (destroy == kspconvergeddefaultdestroycptr_) cctx = *(void ***)cctx;
127     *ierr = PetscObjectSetFortranCallback((PetscObject)*ksp, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.test, (PetscFortranCallbackFn *)converge, cctx);
128     if (*ierr) return;
129     *ierr = PetscObjectSetFortranCallback((PetscObject)*ksp, PETSC_FORTRAN_CALLBACK_CLASS, &_cb.testdestroy, (PetscFortranCallbackFn *)destroy, cctx);
130     if (*ierr) return;
131     *ierr = KSPSetConvergenceTest(*ksp, ourtest, *ksp, ourtestdestroy);
132   }
133 }
134 
kspsetcomputerhs_(KSP * ksp,void (* func)(KSP *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)135 PETSC_EXTERN void kspsetcomputerhs_(KSP *ksp, void (*func)(KSP *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
136 {
137   DM dm;
138   *ierr = KSPGetDM(*ksp, &dm);
139   if (!*ierr) dmkspsetcomputerhs_(&dm, func, ctx, ierr);
140 }
141 
kspsetcomputeinitialguess_(KSP * ksp,void (* func)(KSP *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)142 PETSC_EXTERN void kspsetcomputeinitialguess_(KSP *ksp, void (*func)(KSP *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
143 {
144   DM dm;
145   *ierr = KSPGetDM(*ksp, &dm);
146   if (!*ierr) dmkspsetcomputeinitialguess_(&dm, func, ctx, ierr);
147 }
148 
kspsetcomputeoperators_(KSP * ksp,void (* func)(KSP *,Vec *,void *,PetscErrorCode *),PetscCtx ctx,PetscErrorCode * ierr)149 PETSC_EXTERN void kspsetcomputeoperators_(KSP *ksp, void (*func)(KSP *, Vec *, void *, PetscErrorCode *), PetscCtx ctx, PetscErrorCode *ierr)
150 {
151   DM dm;
152   *ierr = KSPGetDM(*ksp, &dm);
153   if (!*ierr) dmkspsetcomputeoperators_(&dm, func, ctx, ierr);
154 }
155 
kspconvergeddefaultcreate_(PetscFortranAddr * ctx,PetscErrorCode * ierr)156 PETSC_EXTERN void kspconvergeddefaultcreate_(PetscFortranAddr *ctx, PetscErrorCode *ierr)
157 {
158   *ierr = KSPConvergedDefaultCreate((void **)ctx);
159 }
160