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