xref: /petsc/src/snes/interface/ftn-custom/zsnesf.c (revision fb8e56e08d4d0bfe9fc63603ca1f7fddd68abbdb)
1 #include <petsc-private/fortranimpl.h>
2 #include <petscsnes.h>
3 #include <petscviewer.h>
4 
5 #if defined(PETSC_HAVE_FORTRAN_CAPS)
6 #define matmffdcomputejacobian_          MATMFFDCOMPUTEJACOBIAN
7 #define snessolve_                       SNESSOLVE
8 #define snescomputejacobiandefault_      SNESCOMPUTEJACOBIANDEFAULT
9 #define snescomputejacobiandefaultcolor_ SNESCOMPUTEJACOBIANDEFAULTCOLOR
10 #define snessetjacobian_                 SNESSETJACOBIAN
11 #define snesgetoptionsprefix_            SNESGETOPTIONSPREFIX
12 #define snesgettype_                     SNESGETTYPE
13 #define snessetfunction_                 SNESSETFUNCTION
14 #define snessetgs_                       SNESSETGS
15 #define snesgetfunction_                 SNESGETFUNCTION
16 #define snesgetgs_                       SNESGETGS
17 #define snessetconvergencetest_          SNESSETCONVERGENCETEST
18 #define snesconvergeddefault_            SNESCONVERGEDDEFAULT
19 #define snesskipconverged_               SNESSKIPCONVERGED
20 #define snesview_                        SNESVIEW
21 #define snesgetconvergencehistory_       SNESGETCONVERGENCEHISTORY
22 #define snesgetjacobian_                 SNESGETJACOBIAN
23 #define snessettype_                     SNESSETTYPE
24 #define snesappendoptionsprefix_         SNESAPPENDOPTIONSPREFIX
25 #define snessetoptionsprefix_            SNESSETOPTIONSPREFIX
26 #define snesmonitordefault_              SNESMONITORDEFAULT
27 #define snesmonitorsolution_             SNESMONITORSOLUTION
28 #define snesmonitorlgresidualnorm_       SNESMONITORLGRESIDUALNORM
29 #define snesmonitorsolutionupdate_       SNESMONITORSOLUTIONUPDATE
30 #define snesmonitorset_                  SNESMONITORSET
31 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
32 #define matmffdcomputejacobian_          matmffdcomputejacobian
33 #define snessolve_                       snessolve
34 #define snescomputejacobiandefault_      snescomputejacobiandefault
35 #define snescomputejacobiandefaultcolor_ snescomputejacobiandefaultcolor
36 #define snessetjacobian_                 snessetjacobian
37 #define snesgetoptionsprefix_            snesgetoptionsprefix
38 #define snesgettype_                     snesgettype
39 #define snessetfunction_                 snessetfunction
40 #define snessetgs_                       snessetgs
41 #define snesgetfunction_                 snesgetfunction
42 #define snesgetgs_                       snesgetgs
43 #define snessetconvergencetest_          snessetconvergencetest
44 #define snesconvergeddefault_            snesconvergeddefault
45 #define snesskipconverged_               snesskipconverged
46 #define snesview_                        snesview
47 #define snesgetjacobian_                 snesgetjacobian
48 #define snesgetconvergencehistory_       snesgetconvergencehistory
49 #define snessettype_                     snessettype
50 #define snesappendoptionsprefix_         snesappendoptionsprefix
51 #define snessetoptionsprefix_            snessetoptionsprefix
52 #define snesmonitorlgresidualnorm_       snesmonitorlgresidualnorm
53 #define snesmonitordefault_              snesmonitordefault
54 #define snesmonitorsolution_             snesmonitorsolution
55 #define snesmonitorsolutionupdate_       snesmonitorsolutionupdate
56 #define snesmonitorset_                  snesmonitorset
57 #endif
58 
59 static struct {
60   PetscFortranCallbackId function;
61   PetscFortranCallbackId test;
62   PetscFortranCallbackId destroy;
63   PetscFortranCallbackId jacobian;
64   PetscFortranCallbackId monitor;
65   PetscFortranCallbackId mondestroy;
66   PetscFortranCallbackId gs;
67 } _cb;
68 
69 #undef __FUNCT__
70 #define __FUNCT__ "oursnesfunction"
71 static PetscErrorCode oursnesfunction(SNES snes,Vec x,Vec f,void *ctx)
72 {
73   PetscObjectUseFortranCallback(snes,_cb.function,(SNES*,Vec*,Vec*,void*,PetscErrorCode*),(&snes,&x,&f,_ctx,&ierr));
74   return 0;
75 }
76 
77 #undef __FUNCT__
78 #define __FUNCT__ "oursnestest"
79 static PetscErrorCode oursnestest(SNES snes,PetscInt it,PetscReal a,PetscReal d,PetscReal c,SNESConvergedReason *reason,void *ctx)
80 {
81   PetscObjectUseFortranCallback(snes,_cb.test,(SNES*,PetscInt*,PetscReal*,PetscReal*,PetscReal*,SNESConvergedReason*,void*,PetscErrorCode*),(&snes,&it,&a,&d,&c,reason,_ctx,&ierr));
82   return 0;
83 }
84 
85 #undef __FUNCT__
86 #define __FUNCT__ "ourdestroy"
87 static PetscErrorCode ourdestroy(void *ctx)
88 {
89   PetscObjectUseFortranCallback(ctx,_cb.destroy,(void*,PetscErrorCode*),(_ctx,&ierr));
90   return 0;
91 }
92 
93 #undef __FUNCT__
94 #define __FUNCT__ "oursnesjacobian"
95 static PetscErrorCode oursnesjacobian(SNES snes,Vec x,Mat *m,Mat *p,MatStructure *type,void *ctx)
96 {
97   PetscObjectUseFortranCallback(snes,_cb.jacobian,(SNES*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*),(&snes,&x,m,p,type,_ctx,&ierr));
98   return 0;
99 }
100 
101 #undef __FUNCT__
102 #define __FUNCT__ "oursnesgs"
103 static PetscErrorCode oursnesgs(SNES snes,Vec x,Vec b,void *ctx)
104 {
105   PetscObjectUseFortranCallback(snes,_cb.gs,(SNES*,Vec*,Vec*,void*,PetscErrorCode*),(&snes,&x,&b,_ctx,&ierr));
106   return 0;
107 }
108 #undef __FUNCT__
109 #define __FUNCT__ "oursnesmonitor"
110 static PetscErrorCode oursnesmonitor(SNES snes,PetscInt i,PetscReal d,void *ctx)
111 {
112   PetscObjectUseFortranCallback(snes,_cb.monitor,(SNES*,PetscInt*,PetscReal*,void*,PetscErrorCode*),(&snes,&i,&d,_ctx,&ierr));
113   return 0;
114 }
115 #undef __FUNCT__
116 #define __FUNCT__ "ourmondestroy"
117 static PetscErrorCode ourmondestroy(void **ctx)
118 {
119   SNES snes = (SNES)*ctx;
120   PetscObjectUseFortranCallback(snes,_cb.mondestroy,(void*,PetscErrorCode*),(_ctx,&ierr));
121   return 0;
122 }
123 
124 /* ---------------------------------------------------------*/
125 /*
126      snescomputejacobiandefault() and snescomputejacobiandefaultcolor()
127   These can be used directly from Fortran but are mostly so that
128   Fortran SNESSetJacobian() will properly handle the defaults being passed in.
129 
130   functions, hence no STDCALL
131 */
132 PETSC_EXTERN void matmffdcomputejacobian_(SNES *snes,Vec *x,Mat *m,Mat *p,MatStructure *type,void *ctx,PetscErrorCode *ierr)
133 {
134   *ierr = MatMFFDComputeJacobian(*snes,*x,m,p,type,ctx);
135 }
136 PETSC_EXTERN void snescomputejacobiandefault_(SNES *snes,Vec *x,Mat *m,Mat *p,MatStructure *type,void *ctx,PetscErrorCode *ierr)
137 {
138   *ierr = SNESComputeJacobianDefault(*snes,*x,m,p,type,ctx);
139 }
140 PETSC_EXTERN void  snescomputejacobiandefaultcolor_(SNES *snes,Vec *x,Mat *m,Mat *p,MatStructure *type,void *ctx,PetscErrorCode *ierr)
141 {
142   *ierr = SNESComputeJacobianDefaultColor(*snes,*x,m,p,type,*(MatFDColoring*)ctx);
143 }
144 
145 PETSC_EXTERN void PETSC_STDCALL snessetjacobian_(SNES *snes,Mat *A,Mat *B,
146                                     void (PETSC_STDCALL *func)(SNES*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*),
147                                     void *ctx,PetscErrorCode *ierr)
148 {
149   CHKFORTRANNULLOBJECT(ctx);
150   CHKFORTRANNULLFUNCTION(func);
151   if ((PetscVoidFunction)func == (PetscVoidFunction)snescomputejacobiandefault_) {
152     *ierr = SNESSetJacobian(*snes,*A,*B,SNESComputeJacobianDefault,ctx);
153   } else if ((PetscVoidFunction)func == (PetscVoidFunction)snescomputejacobiandefaultcolor_) {
154     *ierr = SNESSetJacobian(*snes,*A,*B,SNESComputeJacobianDefaultColor,*(MatFDColoring*)ctx);
155   } else if ((PetscVoidFunction)func == (PetscVoidFunction)matmffdcomputejacobian_) {
156     *ierr = SNESSetJacobian(*snes,*A,*B,MatMFFDComputeJacobian,ctx);
157   } else if (!func) {
158     *ierr = SNESSetJacobian(*snes,*A,*B,0,ctx);
159   } else {
160     *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.jacobian,(PetscVoidFunction)func,ctx);
161     if (!*ierr) *ierr = SNESSetJacobian(*snes,*A,*B,oursnesjacobian,NULL);
162   }
163 }
164 /* -------------------------------------------------------------*/
165 
166 PETSC_EXTERN void PETSC_STDCALL snessolve_(SNES *snes,Vec *b,Vec *x, int *__ierr)
167 {
168   Vec B = *b,X = *x;
169   if (FORTRANNULLOBJECT(b)) B = NULL;
170   if (FORTRANNULLOBJECT(x)) X = NULL;
171   *__ierr = SNESSolve(*snes,B,X);
172 }
173 
174 PETSC_EXTERN void PETSC_STDCALL snesgetoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
175 {
176   const char *tname;
177 
178   *ierr = SNESGetOptionsPrefix(*snes,&tname);
179   *ierr = PetscStrncpy(prefix,tname,len);if (*ierr) return;
180 }
181 
182 PETSC_EXTERN void PETSC_STDCALL snesgettype_(SNES *snes,CHAR name PETSC_MIXED_LEN(len), PetscErrorCode *ierr PETSC_END_LEN(len))
183 {
184   const char *tname;
185 
186   *ierr = SNESGetType(*snes,&tname);
187   *ierr = PetscStrncpy(name,tname,len);if (*ierr) return;
188   FIXRETURNCHAR(PETSC_TRUE,name,len);
189 }
190 
191 /* ---------------------------------------------------------*/
192 
193 /*
194    These are not usually called from Fortran but allow Fortran users
195    to transparently set these monitors from .F code
196 
197    functions, hence no STDCALL
198 */
199 
200 PETSC_EXTERN void PETSC_STDCALL snessetfunction_(SNES *snes,Vec *r,void (PETSC_STDCALL *func)(SNES*,Vec*,Vec*,void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr)
201 {
202   CHKFORTRANNULLOBJECT(ctx);
203   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.function,(PetscVoidFunction)func,ctx);
204   if (!*ierr) *ierr = SNESSetFunction(*snes,*r,oursnesfunction,NULL);
205 }
206 
207 
208 PETSC_EXTERN void PETSC_STDCALL snessetgs_(SNES *snes,void (PETSC_STDCALL *func)(SNES*,Vec*,Vec*,void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr)
209 {
210   CHKFORTRANNULLOBJECT(ctx);
211   *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.gs,(PetscVoidFunction)func,ctx);
212   if (!*ierr) *ierr = SNESSetGS(*snes,oursnesgs,NULL);
213 }
214 /* ---------------------------------------------------------*/
215 
216 /* the func argument is ignored */
217 PETSC_EXTERN void PETSC_STDCALL snesgetfunction_(SNES *snes,Vec *r,void *func,void **ctx,PetscErrorCode *ierr)
218 {
219   CHKFORTRANNULLINTEGER(ctx);
220   CHKFORTRANNULLOBJECT(r);
221   *ierr = SNESGetFunction(*snes,r,NULL,NULL); if (*ierr) return;
222   *ierr = PetscObjectGetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,_cb.function,NULL,ctx);
223 }
224 
225 PETSC_EXTERN void PETSC_STDCALL snesgetgs_(SNES *snes,void *func,void **ctx,PetscErrorCode *ierr)
226 {
227   CHKFORTRANNULLINTEGER(ctx);
228   *ierr = PetscObjectGetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,_cb.gs,NULL,ctx);
229 }
230 
231 /*----------------------------------------------------------------------*/
232 
233 PETSC_EXTERN void snesconvergeddefault_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r, void *ct,PetscErrorCode *ierr)
234 {
235   *ierr = SNESConvergedDefault(*snes,*it,*a,*b,*c,r,ct);
236 }
237 
238 PETSC_EXTERN void snesskipconverged_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r,void *ct,PetscErrorCode *ierr)
239 {
240   *ierr = SNESSkipConverged(*snes,*it,*a,*b,*c,r,ct);
241 }
242 
243 PETSC_EXTERN void PETSC_STDCALL snessetconvergencetest_(SNES *snes,void (PETSC_STDCALL *func)(SNES*,PetscInt*,PetscReal*,PetscReal*,PetscReal*,SNESConvergedReason*,void*,PetscErrorCode*), void *cctx,void (PETSC_STDCALL *destroy)(void*),PetscErrorCode *ierr)
244 {
245   CHKFORTRANNULLOBJECT(cctx);
246   CHKFORTRANNULLFUNCTION(destroy);
247 
248   if ((PetscVoidFunction)func == (PetscVoidFunction)snesconvergeddefault_) {
249     *ierr = SNESSetConvergenceTest(*snes,SNESConvergedDefault,0,0);
250   } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesskipconverged_) {
251     *ierr = SNESSetConvergenceTest(*snes,SNESSkipConverged,0,0);
252   } else {
253     *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.test,(PetscVoidFunction)func,cctx);
254     if (*ierr) return;
255     if (!destroy) {
256       *ierr = SNESSetConvergenceTest(*snes,oursnestest,*snes,NULL);
257     } else {
258       *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.destroy,(PetscVoidFunction)destroy,cctx);
259       if (!*ierr) *ierr = SNESSetConvergenceTest(*snes,oursnestest,*snes,ourdestroy);
260     }
261   }
262 }
263 /*----------------------------------------------------------------------*/
264 
265 PETSC_EXTERN void PETSC_STDCALL snesview_(SNES *snes,PetscViewer *viewer, PetscErrorCode *ierr)
266 {
267   PetscViewer v;
268   PetscPatchDefaultViewers_Fortran(viewer,v);
269   *ierr = SNESView(*snes,v);
270 }
271 
272 /*  func is currently ignored from Fortran */
273 PETSC_EXTERN void PETSC_STDCALL snesgetjacobian_(SNES *snes,Mat *A,Mat *B,int *func,void **ctx,PetscErrorCode *ierr)
274 {
275   CHKFORTRANNULLINTEGER(ctx);
276   CHKFORTRANNULLOBJECT(A);
277   CHKFORTRANNULLOBJECT(B);
278   *ierr = SNESGetJacobian(*snes,A,B,0,NULL); if (*ierr) return;
279   *ierr = PetscObjectGetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,_cb.jacobian,NULL,ctx);
280 
281 }
282 
283 PETSC_EXTERN void PETSC_STDCALL snesgetconvergencehistory_(SNES *snes,PetscInt *na,PetscErrorCode *ierr)
284 {
285   *ierr = SNESGetConvergenceHistory(*snes,NULL,NULL,na);
286 }
287 
288 PETSC_EXTERN void PETSC_STDCALL snessettype_(SNES *snes,CHAR type PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
289 {
290   char *t;
291 
292   FIXCHAR(type,len,t);
293   *ierr = SNESSetType(*snes,t);
294   FREECHAR(type,t);
295 }
296 
297 PETSC_EXTERN void PETSC_STDCALL snesappendoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
298 {
299   char *t;
300 
301   FIXCHAR(prefix,len,t);
302   *ierr = SNESAppendOptionsPrefix(*snes,t);
303   FREECHAR(prefix,t);
304 }
305 
306 PETSC_EXTERN void PETSC_STDCALL snessetoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
307 {
308   char *t;
309 
310   FIXCHAR(prefix,len,t);
311   *ierr = SNESSetOptionsPrefix(*snes,t);
312   FREECHAR(prefix,t);
313 }
314 
315 /*----------------------------------------------------------------------*/
316 /* functions, hence no STDCALL */
317 
318 PETSC_EXTERN void snesmonitorlgresidualnorm_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr)
319 {
320   *ierr = SNESMonitorLGResidualNorm(*snes,*its,*fgnorm,dummy);
321 }
322 
323 PETSC_EXTERN void snesmonitordefault_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr)
324 {
325   *ierr = SNESMonitorDefault(*snes,*its,*fgnorm,dummy);
326 }
327 
328 PETSC_EXTERN void snesmonitorsolution_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr)
329 {
330   *ierr = SNESMonitorSolution(*snes,*its,*fgnorm,dummy);
331 }
332 
333 PETSC_EXTERN void snesmonitorsolutionupdate_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr)
334 {
335   *ierr = SNESMonitorSolutionUpdate(*snes,*its,*fgnorm,dummy);
336 }
337 
338 
339 PETSC_EXTERN void PETSC_STDCALL snesmonitorset_(SNES *snes,void (PETSC_STDCALL *func)(SNES*,PetscInt*,PetscReal*,void*,PetscErrorCode*),void *mctx,void (PETSC_STDCALL *mondestroy)(void*,PetscErrorCode*),PetscErrorCode *ierr)
340 {
341   CHKFORTRANNULLOBJECT(mctx);
342   if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitordefault_) {
343     *ierr = SNESMonitorSet(*snes,SNESMonitorDefault,0,0);
344   } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorsolution_) {
345     *ierr = SNESMonitorSet(*snes,SNESMonitorSolution,0,0);
346   } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorsolutionupdate_) {
347     *ierr = SNESMonitorSet(*snes,SNESMonitorSolutionUpdate,0,0);
348   } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorlgresidualnorm_) {
349     *ierr = SNESMonitorSet(*snes,SNESMonitorLGResidualNorm,0,0);
350   } else {
351     *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.monitor,(PetscVoidFunction)func,mctx);
352     if (*ierr) return;
353     if (FORTRANNULLFUNCTION(mondestroy)) {
354       *ierr = SNESMonitorSet(*snes,oursnesmonitor,*snes,NULL);
355     } else {
356       CHKFORTRANNULLFUNCTION(mondestroy);
357       *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.mondestroy,(PetscVoidFunction)mondestroy,mctx);
358       if (!*ierr) *ierr = SNESMonitorSet(*snes,oursnesmonitor,*snes,ourmondestroy);
359     }
360   }
361 }
362 
363