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