xref: /petsc/src/snes/interface/ftn-custom/zsnesf.c (revision 025f1a0480dcff154a24a80c8d03877a2ba7e925)
1 #include "private/fortranimpl.h"
2 #include "petscsnes.h"
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5 #define snessolve_                       SNESSOLVE
6 #define snesdefaultcomputejacobian_      SNESDEFAULTCOMPUTEJACOBIAN
7 #define snesdefaultcomputejacobiancolor_ SNESDEFAULTCOMPUTEJACOBIANCOLOR
8 #define snesdacomputejacobian_           SNESDACOMPUTEJACOBIAN
9 #define snesdacomputejacobianwithadifor_ SNESDACOMPUTEJACOBIANWITHADIFOR
10 #define snessetjacobian_                 SNESSETJACOBIAN
11 #define snesgetoptionsprefix_            SNESGETOPTIONSPREFIX
12 #define snesgettype_                     SNESGETTYPE
13 #define snesdaformfunction_              SNESDAFORMFUNCTION
14 #define snessetfunction_                 SNESSETFUNCTION
15 #define snesgetfunction_                 SNESGETFUNCTION
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 snesgetapplicationcontext_       SNESGETAPPLICATIONCONTEXT
31 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
32 #define snessolve_                       snessolve
33 #define snesdefaultcomputejacobian_      snesdefaultcomputejacobian
34 #define snesdefaultcomputejacobiancolor_ snesdefaultcomputejacobiancolor
35 #define snesdacomputejacobian_           snesdacomputejacobian
36 #define snesdacomputejacobianwithadifor_ snesdacomputejacobianwithadifor
37 #define snessetjacobian_                 snessetjacobian
38 #define snesgetoptionsprefix_            snesgetoptionsprefix
39 #define snesgettype_                     snesgettype
40 #define snesdaformfunction_              snesdaformfunction
41 #define snessetfunction_                 snessetfunction
42 #define snesgetfunction_                 snesgetfunction
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 snesgetapplicationcontext_       snesgetapplicationcontext
58 #endif
59 
60 static PetscErrorCode oursnesfunction(SNES snes,Vec x,Vec f,void *ctx)
61 {
62   PetscErrorCode ierr = 0;
63   (*(void (PETSC_STDCALL *)(SNES*,Vec*,Vec*,void*,PetscErrorCode*))(((PetscObject)snes)->fortran_func_pointers[0]))(&snes,&x,&f,ctx,&ierr);CHKERRQ(ierr);
64   return 0;
65 }
66 
67 static PetscErrorCode oursnestest(SNES snes,PetscInt it,PetscReal a,PetscReal d,PetscReal c,SNESConvergedReason*reason,void*ctx)
68 {
69   PetscErrorCode ierr = 0;
70   (*(void (PETSC_STDCALL *)(SNES*,PetscInt*,PetscReal*,PetscReal*,PetscReal*,SNESConvergedReason*,void*,PetscErrorCode*))(((PetscObject)snes)->fortran_func_pointers[1]))(&snes,&it,&a,&d,&c,reason,ctx,&ierr);CHKERRQ(ierr);
71   return 0;
72 }
73 
74 static PetscErrorCode oursnesjacobian(SNES snes,Vec x,Mat* m,Mat* p,MatStructure* type,void*ctx)
75 {
76   PetscErrorCode ierr = 0;
77   (*(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);
78   return 0;
79 }
80 static PetscErrorCode oursnesmonitor(SNES snes,PetscInt i,PetscReal d,void*ctx)
81 {
82   PetscErrorCode ierr = 0;
83 
84   void           (*mctx)(void) = ((PetscObject)snes)->fortran_func_pointers[4];
85   (*(void (PETSC_STDCALL *)(SNES*,PetscInt*,PetscReal*,void*,PetscErrorCode*))(((PetscObject)snes)->fortran_func_pointers[3]))(&snes,&i,&d,(void*)mctx,&ierr);CHKERRQ(ierr);
86   return 0;
87 }
88 static PetscErrorCode ourmondestroy(void* ctx)
89 {
90   PetscErrorCode ierr = 0;
91   SNES           snes = (SNES)ctx;
92   void           (*mctx)(void) = ((PetscObject)snes)->fortran_func_pointers[4];
93   (*(void (PETSC_STDCALL *)(PetscVoidFunction,PetscErrorCode*))(((PetscObject)snes)->fortran_func_pointers[5]))(mctx,&ierr);CHKERRQ(ierr);
94   return 0;
95 }
96 
97 EXTERN_C_BEGIN
98 /* ---------------------------------------------------------*/
99 /*
100      snesdefaultcomputejacobian() and snesdefaultcomputejacobiancolor()
101   These can be used directly from Fortran but are mostly so that
102   Fortran SNESSetJacobian() will properly handle the defaults being passed in.
103 
104   functions, hence no STDCALL
105 */
106 void snesdefaultcomputejacobian_(SNES *snes,Vec *x,Mat *m,Mat *p,MatStructure* type,void *ctx,PetscErrorCode *ierr)
107 {
108   *ierr = SNESDefaultComputeJacobian(*snes,*x,m,p,type,ctx);
109 }
110 void  snesdefaultcomputejacobiancolor_(SNES *snes,Vec *x,Mat *m,Mat *p,MatStructure* type,void *ctx,PetscErrorCode *ierr)
111 {
112   *ierr = SNESDefaultComputeJacobianColor(*snes,*x,m,p,type,*(MatFDColoring*)ctx);
113 }
114 
115 void  snesdacomputejacobianwithadifor_(SNES *snes,Vec *X,Mat *m,Mat *p,MatStructure* type,void *ctx,PetscErrorCode *ierr)
116 {
117   (*PetscErrorPrintf)("Cannot call this function from Fortran");
118   *ierr = 1;
119 }
120 
121 void  snesdacomputejacobian_(SNES *snes,Vec *X,Mat *m,Mat *p,MatStructure* type,void *ctx,PetscErrorCode *ierr)
122 {
123   (*PetscErrorPrintf)("Cannot call this function from Fortran");
124   *ierr = 1;
125 }
126 
127 void PETSC_STDCALL snessetjacobian_(SNES *snes,Mat *A,Mat *B,void (PETSC_STDCALL *func)(SNES*,Vec*,Mat*,Mat*,
128             MatStructure*,void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr)
129 {
130   CHKFORTRANNULLOBJECT(ctx);
131   CHKFORTRANNULLFUNCTION(func);
132   PetscObjectAllocateFortranPointers(*snes,10);
133   if ((PetscVoidFunction)func == (PetscVoidFunction)snesdefaultcomputejacobian_) {
134     *ierr = SNESSetJacobian(*snes,*A,*B,SNESDefaultComputeJacobian,ctx);
135   } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesdefaultcomputejacobiancolor_) {
136     *ierr = SNESSetJacobian(*snes,*A,*B,SNESDefaultComputeJacobianColor,*(MatFDColoring*)ctx);
137   } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesdacomputejacobianwithadifor_) {
138     *ierr = SNESSetJacobian(*snes,*A,*B,SNESDAComputeJacobianWithAdifor,ctx);
139   } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesdacomputejacobian_) {
140     *ierr = SNESSetJacobian(*snes,*A,*B,SNESDAComputeJacobian,ctx);
141   } else if (!func) {
142     *ierr = SNESSetJacobian(*snes,*A,*B,0,ctx);
143   } else {
144     ((PetscObject)*snes)->fortran_func_pointers[2] = (PetscVoidFunction)func;
145     *ierr = SNESSetJacobian(*snes,*A,*B,oursnesjacobian,ctx);
146   }
147 }
148 /* -------------------------------------------------------------*/
149 
150 void PETSC_STDCALL   snessolve_(SNES *snes,Vec *b,Vec *x, int *__ierr )
151 {
152   Vec B = *b;
153   if (*b == PETSC_NULL_OBJECT_Fortran) B = PETSC_NULL;
154   *__ierr = SNESSolve(*snes,B,*x);
155 }
156 
157 void PETSC_STDCALL snesgetoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
158 {
159   const char *tname;
160 
161   *ierr = SNESGetOptionsPrefix(*snes,&tname);
162   *ierr = PetscStrncpy(prefix,tname,len);if (*ierr) return;
163 }
164 
165 void PETSC_STDCALL snesgettype_(SNES *snes,CHAR name PETSC_MIXED_LEN(len),
166                                 PetscErrorCode *ierr PETSC_END_LEN(len))
167 {
168   const char *tname;
169 
170   *ierr = SNESGetType(*snes,&tname);
171   *ierr = PetscStrncpy(name,tname,len);if (*ierr) return;
172   FIXRETURNCHAR(PETSC_TRUE,name,len);
173 }
174 
175 void PETSC_STDCALL snesgetapplicationcontext_(SNES *snes,void **ctx,PetscErrorCode *ierr)
176 {
177   *ierr = SNESGetApplicationContext(*snes,ctx);
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 void  snesdaformfunction_(SNES *snes,Vec *X, Vec *F,void *ptr,PetscErrorCode *ierr)
188 {
189   *ierr = SNESDAFormFunction(*snes,*X,*F,ptr);
190 }
191 
192 void PETSC_STDCALL snessetfunction_(SNES *snes,Vec *r,void (PETSC_STDCALL *func)(SNES*,Vec*,Vec*,void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr)
193 {
194   CHKFORTRANNULLOBJECT(ctx);
195   PetscObjectAllocateFortranPointers(*snes,10);
196   if ((PetscVoidFunction)func == (PetscVoidFunction)snesdaformfunction_) {
197     *ierr = SNESSetFunction(*snes,*r,SNESDAFormFunction,ctx);
198   } else {
199     ((PetscObject)*snes)->fortran_func_pointers[0] = (PetscVoidFunction)func;
200     *ierr = SNESSetFunction(*snes,*r,oursnesfunction,ctx);
201   }
202 }
203 /* ---------------------------------------------------------*/
204 
205 /* the func argument is ignored */
206 void PETSC_STDCALL snesgetfunction_(SNES *snes,Vec *r,void *func,void **ctx,PetscErrorCode *ierr)
207 {
208   CHKFORTRANNULLINTEGER(ctx);
209   CHKFORTRANNULLOBJECT(r);
210   *ierr = SNESGetFunction(*snes,r,PETSC_NULL,ctx);
211 }
212 /*----------------------------------------------------------------------*/
213 
214 void snesdefaultconverged_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r,
215                                        void *ct,PetscErrorCode *ierr)
216 {
217   *ierr = SNESDefaultConverged(*snes,*it,*a,*b,*c,r,ct);
218 }
219 
220 void snesskipconverged_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r,
221                                        void *ct,PetscErrorCode *ierr)
222 {
223   *ierr = SNESSkipConverged(*snes,*it,*a,*b,*c,r,ct);
224 }
225 
226 void PETSC_STDCALL snessetconvergencetest_(SNES *snes,
227        void (PETSC_STDCALL *func)(SNES*,PetscInt*,PetscReal*,PetscReal*,PetscReal*,SNESConvergedReason*,void*,PetscErrorCode*),
228        void *cctx,PetscErrorCode *ierr)
229 {
230   CHKFORTRANNULLOBJECT(cctx);
231   PetscObjectAllocateFortranPointers(*snes,10);
232 
233   if ((PetscVoidFunction)func == (PetscVoidFunction)snesdefaultconverged_){
234     *ierr = SNESSetConvergenceTest(*snes,SNESDefaultConverged,0);
235   } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesskipconverged_){
236     *ierr = SNESSetConvergenceTest(*snes,SNESSkipConverged,0);
237   } else {
238     ((PetscObject)*snes)->fortran_func_pointers[1] = (PetscVoidFunction)func;
239     *ierr = SNESSetConvergenceTest(*snes,oursnestest,cctx);
240   }
241 }
242 /*----------------------------------------------------------------------*/
243 
244 void PETSC_STDCALL snesview_(SNES *snes,PetscViewer *viewer, PetscErrorCode *ierr)
245 {
246   PetscViewer v;
247   PetscPatchDefaultViewers_Fortran(viewer,v);
248   *ierr = SNESView(*snes,v);
249 }
250 
251 /*  func is currently ignored from Fortran */
252 void PETSC_STDCALL snesgetjacobian_(SNES *snes,Mat *A,Mat *B,int *func,void **ctx,PetscErrorCode *ierr)
253 {
254   CHKFORTRANNULLINTEGER(ctx);
255   CHKFORTRANNULLOBJECT(A);
256   CHKFORTRANNULLOBJECT(B);
257   *ierr = SNESGetJacobian(*snes,A,B,0,ctx);
258 }
259 
260 void PETSC_STDCALL snesgetconvergencehistory_(SNES *snes,PetscInt *na,PetscErrorCode *ierr)
261 {
262   *ierr = SNESGetConvergenceHistory(*snes,PETSC_NULL,PETSC_NULL,na);
263 }
264 
265 void PETSC_STDCALL snessettype_(SNES *snes,CHAR type PETSC_MIXED_LEN(len),
266                                 PetscErrorCode *ierr PETSC_END_LEN(len))
267 {
268   char *t;
269 
270   FIXCHAR(type,len,t);
271   *ierr = SNESSetType(*snes,t);
272   FREECHAR(type,t);
273 }
274 
275 void PETSC_STDCALL snesappendoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len),
276                                             PetscErrorCode *ierr PETSC_END_LEN(len))
277 {
278   char *t;
279 
280   FIXCHAR(prefix,len,t);
281   *ierr = SNESAppendOptionsPrefix(*snes,t);
282   FREECHAR(prefix,t);
283 }
284 
285 void PETSC_STDCALL snessetoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len),
286                                         PetscErrorCode *ierr PETSC_END_LEN(len))
287 {
288   char *t;
289 
290   FIXCHAR(prefix,len,t);
291   *ierr = SNESSetOptionsPrefix(*snes,t);
292   FREECHAR(prefix,t);
293 }
294 
295 /*----------------------------------------------------------------------*/
296 /* functions, hence no STDCALL */
297 
298 void snesmonitorlg_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr)
299 {
300   *ierr = SNESMonitorLG(*snes,*its,*fgnorm,dummy);
301 }
302 
303 void snesmonitordefault_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr)
304 {
305   *ierr = SNESMonitorDefault(*snes,*its,*fgnorm,dummy);
306 }
307 
308 void snesmonitorsolution_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr)
309 {
310   *ierr = SNESMonitorSolution(*snes,*its,*fgnorm,dummy);
311 }
312 
313 void snesmonitorsolutionupdate_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr)
314 {
315   *ierr = SNESMonitorSolutionUpdate(*snes,*its,*fgnorm,dummy);
316 }
317 
318 
319 void PETSC_STDCALL snesmonitorset_(SNES *snes,void (PETSC_STDCALL *func)(SNES*,PetscInt*,PetscReal*,void*,PetscErrorCode*),
320                     void *mctx,void (PETSC_STDCALL *mondestroy)(void*,PetscErrorCode*),PetscErrorCode *ierr)
321 {
322   CHKFORTRANNULLOBJECT(mctx);
323   PetscObjectAllocateFortranPointers(*snes,10);
324   if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitordefault_) {
325     *ierr = SNESMonitorSet(*snes,SNESMonitorDefault,0,0);
326   } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorsolution_) {
327     *ierr = SNESMonitorSet(*snes,SNESMonitorSolution,0,0);
328   } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorsolutionupdate_) {
329     *ierr = SNESMonitorSet(*snes,SNESMonitorSolutionUpdate,0,0);
330   } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorlg_) {
331     *ierr = SNESMonitorSet(*snes,SNESMonitorLG,0,0);
332   } else {
333     ((PetscObject)*snes)->fortran_func_pointers[3] = (PetscVoidFunction)func;
334     ((PetscObject)*snes)->fortran_func_pointers[4] = (PetscVoidFunction)mctx;
335 
336     if (FORTRANNULLFUNCTION(mondestroy)){
337       *ierr = SNESMonitorSet(*snes,oursnesmonitor,*snes,0);
338     } else {
339       ((PetscObject)*snes)->fortran_func_pointers[5] = (PetscVoidFunction)mondestroy;
340       *ierr = SNESMonitorSet(*snes,oursnesmonitor,*snes,ourmondestroy);
341     }
342   }
343 }
344 
345 
346 
347 EXTERN_C_END
348