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