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