1 #include <petsc/private/fortranimpl.h> 2 #include <petscsnes.h> 3 #include <petscviewer.h> 4 #include <../src/sys/f90-src/f90impl.h> 5 6 #if defined(PETSC_HAVE_FORTRAN_CAPS) 7 #define matmffdcomputejacobian_ MATMFFDCOMPUTEJACOBIAN 8 #define snessolve_ SNESSOLVE 9 #define snescomputejacobiandefault_ SNESCOMPUTEJACOBIANDEFAULT 10 #define snescomputejacobiandefaultcolor_ SNESCOMPUTEJACOBIANDEFAULTCOLOR 11 #define snessetjacobian_ SNESSETJACOBIAN 12 #define snesgetoptionsprefix_ SNESGETOPTIONSPREFIX 13 #define snesgettype_ SNESGETTYPE 14 #define snessetfunction_ SNESSETFUNCTION 15 #define snessetngs_ SNESSETNGS 16 #define snessetupdate_ SNESSETUPDATE 17 #define snesgetfunction_ SNESGETFUNCTION 18 #define snesgetngs_ SNESGETNGS 19 #define snessetconvergencetest_ SNESSETCONVERGENCETEST 20 #define snesconvergeddefault_ SNESCONVERGEDDEFAULT 21 #define snesconvergedskip_ SNESCONVERGEDSKIP 22 #define snesview_ SNESVIEW 23 #define snesgetconvergencehistory_ SNESGETCONVERGENCEHISTORY 24 #define snesgetjacobian_ SNESGETJACOBIAN 25 #define snessettype_ SNESSETTYPE 26 #define snesappendoptionsprefix_ SNESAPPENDOPTIONSPREFIX 27 #define snessetoptionsprefix_ SNESSETOPTIONSPREFIX 28 #define snesmonitordefault_ SNESMONITORDEFAULT 29 #define snesmonitorsolution_ SNESMONITORSOLUTION 30 #define snesmonitorlgresidualnorm_ SNESMONITORLGRESIDUALNORM 31 #define snesmonitorsolutionupdate_ SNESMONITORSOLUTIONUPDATE 32 #define snesmonitorset_ SNESMONITORSET 33 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 34 #define matmffdcomputejacobian_ matmffdcomputejacobian 35 #define snessolve_ snessolve 36 #define snescomputejacobiandefault_ snescomputejacobiandefault 37 #define snescomputejacobiandefaultcolor_ snescomputejacobiandefaultcolor 38 #define snessetjacobian_ snessetjacobian 39 #define snesgetoptionsprefix_ snesgetoptionsprefix 40 #define snesgettype_ snesgettype 41 #define snessetfunction_ snessetfunction 42 #define snessetngs_ snessetngs 43 #define snessetupdate_ snessetupdate 44 #define snesgetfunction_ snesgetfunction 45 #define snesgetngs_ snesgetngs 46 #define snessetconvergencetest_ snessetconvergencetest 47 #define snesconvergeddefault_ snesconvergeddefault 48 #define snesconvergedskip_ snesconvergedskip 49 #define snesview_ snesview 50 #define snesgetjacobian_ snesgetjacobian 51 #define snesgetconvergencehistory_ snesgetconvergencehistory 52 #define snessettype_ snessettype 53 #define snesappendoptionsprefix_ snesappendoptionsprefix 54 #define snessetoptionsprefix_ snessetoptionsprefix 55 #define snesmonitorlgresidualnorm_ snesmonitorlgresidualnorm 56 #define snesmonitordefault_ snesmonitordefault 57 #define snesmonitorsolution_ snesmonitorsolution 58 #define snesmonitorsolutionupdate_ snesmonitorsolutionupdate 59 #define snesmonitorset_ snesmonitorset 60 #endif 61 62 static struct { 63 PetscFortranCallbackId function; 64 PetscFortranCallbackId test; 65 PetscFortranCallbackId destroy; 66 PetscFortranCallbackId jacobian; 67 PetscFortranCallbackId monitor; 68 PetscFortranCallbackId mondestroy; 69 PetscFortranCallbackId ngs; 70 PetscFortranCallbackId update; 71 #if defined(PETSC_HAVE_F90_2PTR_ARG) 72 PetscFortranCallbackId function_pgiptr; 73 #endif 74 } _cb; 75 76 #undef __FUNCT__ 77 #define __FUNCT__ "oursnesfunction" 78 static PetscErrorCode oursnesfunction(SNES snes,Vec x,Vec f,void *ctx) 79 { 80 #if defined(PETSC_HAVE_F90_2PTR_ARG) 81 void* ptr; 82 PetscObjectGetFortranCallback((PetscObject)snes,PETSC_FORTRAN_CALLBACK_CLASS,_cb.function_pgiptr,NULL,&ptr); 83 #endif 84 PetscObjectUseFortranCallback(snes,_cb.function,(SNES*,Vec*,Vec*,void*,PetscErrorCode* PETSC_F90_2PTR_PROTO_NOVAR),(&snes,&x,&f,_ctx,&ierr PETSC_F90_2PTR_PARAM(ptr))); 85 } 86 87 #undef __FUNCT__ 88 #define __FUNCT__ "oursnestest" 89 static PetscErrorCode oursnestest(SNES snes,PetscInt it,PetscReal a,PetscReal d,PetscReal c,SNESConvergedReason *reason,void *ctx) 90 { 91 PetscObjectUseFortranCallback(snes,_cb.test,(SNES*,PetscInt*,PetscReal*,PetscReal*,PetscReal*,SNESConvergedReason*,void*,PetscErrorCode*),(&snes,&it,&a,&d,&c,reason,_ctx,&ierr)); 92 } 93 94 #undef __FUNCT__ 95 #define __FUNCT__ "ourdestroy" 96 static PetscErrorCode ourdestroy(void *ctx) 97 { 98 PetscObjectUseFortranCallback(ctx,_cb.destroy,(void*,PetscErrorCode*),(_ctx,&ierr)); 99 } 100 101 #undef __FUNCT__ 102 #define __FUNCT__ "oursnesjacobian" 103 static PetscErrorCode oursnesjacobian(SNES snes,Vec x,Mat m,Mat p,void *ctx) 104 { 105 PetscObjectUseFortranCallback(snes,_cb.jacobian,(SNES*,Vec*,Mat*,Mat*,void*,PetscErrorCode*),(&snes,&x,&m,&p,_ctx,&ierr)); 106 } 107 108 #undef __FUNCT__ 109 #define __FUNCT__ "oursnesupdate" 110 static PetscErrorCode oursnesupdate(SNES snes,PetscInt i) 111 { 112 PetscObjectUseFortranCallback(snes,_cb.update,(SNES*,PetscInt *,PetscErrorCode*),(&snes,&i,&ierr)); 113 } 114 #undef __FUNCT__ 115 #define __FUNCT__ "oursnesngs" 116 static PetscErrorCode oursnesngs(SNES snes,Vec x,Vec b,void *ctx) 117 { 118 PetscObjectUseFortranCallback(snes,_cb.ngs,(SNES*,Vec*,Vec*,void*,PetscErrorCode*),(&snes,&x,&b,_ctx,&ierr)); 119 } 120 #undef __FUNCT__ 121 #define __FUNCT__ "oursnesmonitor" 122 static PetscErrorCode oursnesmonitor(SNES snes,PetscInt i,PetscReal d,void *ctx) 123 { 124 PetscObjectUseFortranCallback(snes,_cb.monitor,(SNES*,PetscInt*,PetscReal*,void*,PetscErrorCode*),(&snes,&i,&d,_ctx,&ierr)); 125 } 126 #undef __FUNCT__ 127 #define __FUNCT__ "ourmondestroy" 128 static PetscErrorCode ourmondestroy(void **ctx) 129 { 130 SNES snes = (SNES)*ctx; 131 PetscObjectUseFortranCallback(snes,_cb.mondestroy,(void*,PetscErrorCode*),(_ctx,&ierr)); 132 } 133 134 /* ---------------------------------------------------------*/ 135 /* 136 snescomputejacobiandefault() and snescomputejacobiandefaultcolor() 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 PETSC_EXTERN void matmffdcomputejacobian_(SNES *snes,Vec *x,Mat *m,Mat *p,void *ctx,PetscErrorCode *ierr) 143 { 144 *ierr = MatMFFDComputeJacobian(*snes,*x,*m,*p,ctx); 145 } 146 PETSC_EXTERN void snescomputejacobiandefault_(SNES *snes,Vec *x,Mat *m,Mat *p,void *ctx,PetscErrorCode *ierr) 147 { 148 *ierr = SNESComputeJacobianDefault(*snes,*x,*m,*p,ctx); 149 } 150 PETSC_EXTERN void snescomputejacobiandefaultcolor_(SNES *snes,Vec *x,Mat *m,Mat *p,void *ctx,PetscErrorCode *ierr) 151 { 152 *ierr = SNESComputeJacobianDefaultColor(*snes,*x,*m,*p,*(MatFDColoring*)ctx); 153 } 154 155 PETSC_EXTERN void PETSC_STDCALL snessetjacobian_(SNES *snes,Mat *A,Mat *B, 156 void (PETSC_STDCALL *func)(SNES*,Vec*,Mat*,Mat*,void*,PetscErrorCode*), 157 void *ctx,PetscErrorCode *ierr) 158 { 159 CHKFORTRANNULLFUNCTION(func); 160 if ((PetscVoidFunction)func == (PetscVoidFunction)snescomputejacobiandefault_) { 161 *ierr = SNESSetJacobian(*snes,*A,*B,SNESComputeJacobianDefault,ctx); 162 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snescomputejacobiandefaultcolor_) { 163 if (!ctx) { 164 *ierr = PETSC_ERR_ARG_NULL; 165 return; 166 } 167 *ierr = SNESSetJacobian(*snes,*A,*B,SNESComputeJacobianDefaultColor,*(MatFDColoring*)ctx); 168 } else if ((PetscVoidFunction)func == (PetscVoidFunction)matmffdcomputejacobian_) { 169 *ierr = SNESSetJacobian(*snes,*A,*B,MatMFFDComputeJacobian,ctx); 170 } else if (!func) { 171 *ierr = SNESSetJacobian(*snes,*A,*B,0,ctx); 172 } else { 173 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.jacobian,(PetscVoidFunction)func,ctx); 174 if (!*ierr) *ierr = SNESSetJacobian(*snes,*A,*B,oursnesjacobian,NULL); 175 } 176 } 177 /* -------------------------------------------------------------*/ 178 179 PETSC_EXTERN void PETSC_STDCALL snessolve_(SNES *snes,Vec *b,Vec *x, int *ierr) 180 { 181 CHKFORTRANNULLOBJECTDEREFERENCE(b); 182 CHKFORTRANNULLOBJECTDEREFERENCE(x); 183 *ierr = SNESSolve(*snes,*b,*x); 184 } 185 186 PETSC_EXTERN void PETSC_STDCALL snesgetoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len)) 187 { 188 const char *tname; 189 190 *ierr = SNESGetOptionsPrefix(*snes,&tname); 191 *ierr = PetscStrncpy(prefix,tname,len);if (*ierr) return; 192 FIXRETURNCHAR(PETSC_TRUE,prefix,len); 193 } 194 195 PETSC_EXTERN void PETSC_STDCALL snesgettype_(SNES *snes,CHAR name PETSC_MIXED_LEN(len), PetscErrorCode *ierr PETSC_END_LEN(len)) 196 { 197 const char *tname; 198 199 *ierr = SNESGetType(*snes,&tname); 200 *ierr = PetscStrncpy(name,tname,len);if (*ierr) return; 201 FIXRETURNCHAR(PETSC_TRUE,name,len); 202 } 203 204 /* ---------------------------------------------------------*/ 205 206 /* 207 These are not usually called from Fortran but allow Fortran users 208 to transparently set these monitors from .F code 209 210 functions, hence no STDCALL 211 */ 212 213 PETSC_EXTERN void PETSC_STDCALL snessetfunction_(SNES *snes,Vec *r,void (PETSC_STDCALL *func)(SNES*,Vec*,Vec*,void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(ptr)) 214 { 215 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.function,(PetscVoidFunction)func,ctx); 216 #if defined(PETSC_HAVE_F90_2PTR_ARG) 217 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.function_pgiptr,PETSC_NULL,ptr); 218 #endif 219 if (!*ierr) *ierr = SNESSetFunction(*snes,*r,oursnesfunction,NULL); 220 } 221 222 223 PETSC_EXTERN void PETSC_STDCALL snessetngs_(SNES *snes,void (PETSC_STDCALL *func)(SNES*,Vec*,Vec*,void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr) 224 { 225 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.ngs,(PetscVoidFunction)func,ctx); 226 if (!*ierr) *ierr = SNESSetNGS(*snes,oursnesngs,NULL); 227 } 228 PETSC_EXTERN void PETSC_STDCALL snessetupdate_(SNES *snes,void (PETSC_STDCALL *func)(SNES*,PetscInt*,PetscErrorCode*),PetscErrorCode *ierr) 229 { 230 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.update,(PetscVoidFunction)func,NULL); 231 if (!*ierr) *ierr = SNESSetUpdate(*snes,oursnesupdate); 232 } 233 /* ---------------------------------------------------------*/ 234 235 /* the func argument is ignored */ 236 PETSC_EXTERN void PETSC_STDCALL snesgetfunction_(SNES *snes,Vec *r,void *func,void **ctx,PetscErrorCode *ierr) 237 { 238 CHKFORTRANNULLOBJECT(r); 239 *ierr = SNESGetFunction(*snes,r,NULL,NULL); if (*ierr) return; 240 *ierr = PetscObjectGetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,_cb.function,NULL,ctx); 241 } 242 243 PETSC_EXTERN void PETSC_STDCALL snesgetngs_(SNES *snes,void *func,void **ctx,PetscErrorCode *ierr) 244 { 245 *ierr = PetscObjectGetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,_cb.ngs,NULL,ctx); 246 } 247 248 /*----------------------------------------------------------------------*/ 249 250 PETSC_EXTERN void snesconvergeddefault_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r, void *ct,PetscErrorCode *ierr) 251 { 252 *ierr = SNESConvergedDefault(*snes,*it,*a,*b,*c,r,ct); 253 } 254 255 PETSC_EXTERN void snesconvergedskip_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r,void *ct,PetscErrorCode *ierr) 256 { 257 *ierr = SNESConvergedSkip(*snes,*it,*a,*b,*c,r,ct); 258 } 259 260 PETSC_EXTERN 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) 261 { 262 CHKFORTRANNULLOBJECT(cctx); 263 CHKFORTRANNULLFUNCTION(destroy); 264 265 if ((PetscVoidFunction)func == (PetscVoidFunction)snesconvergeddefault_) { 266 *ierr = SNESSetConvergenceTest(*snes,SNESConvergedDefault,0,0); 267 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesconvergedskip_) { 268 *ierr = SNESSetConvergenceTest(*snes,SNESConvergedSkip,0,0); 269 } else { 270 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.test,(PetscVoidFunction)func,cctx); 271 if (*ierr) return; 272 if (!destroy) { 273 *ierr = SNESSetConvergenceTest(*snes,oursnestest,*snes,NULL); 274 } else { 275 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.destroy,(PetscVoidFunction)destroy,cctx); 276 if (!*ierr) *ierr = SNESSetConvergenceTest(*snes,oursnestest,*snes,ourdestroy); 277 } 278 } 279 } 280 /*----------------------------------------------------------------------*/ 281 282 PETSC_EXTERN void PETSC_STDCALL snesview_(SNES *snes,PetscViewer *viewer, PetscErrorCode *ierr) 283 { 284 PetscViewer v; 285 PetscPatchDefaultViewers_Fortran(viewer,v); 286 *ierr = SNESView(*snes,v); 287 } 288 289 /* func is currently ignored from Fortran */ 290 PETSC_EXTERN void PETSC_STDCALL snesgetjacobian_(SNES *snes,Mat *A,Mat *B,int *func,void **ctx,PetscErrorCode *ierr) 291 { 292 CHKFORTRANNULLINTEGER(ctx); 293 CHKFORTRANNULLOBJECT(A); 294 CHKFORTRANNULLOBJECT(B); 295 *ierr = SNESGetJacobian(*snes,A,B,0,NULL); if (*ierr) return; 296 *ierr = PetscObjectGetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,_cb.jacobian,NULL,ctx); 297 298 } 299 300 PETSC_EXTERN void PETSC_STDCALL snesgetconvergencehistory_(SNES *snes,PetscInt *na,PetscErrorCode *ierr) 301 { 302 *ierr = SNESGetConvergenceHistory(*snes,NULL,NULL,na); 303 } 304 305 PETSC_EXTERN void PETSC_STDCALL snessettype_(SNES *snes,CHAR type PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len)) 306 { 307 char *t; 308 309 FIXCHAR(type,len,t); 310 *ierr = SNESSetType(*snes,t); 311 FREECHAR(type,t); 312 } 313 314 PETSC_EXTERN void PETSC_STDCALL snesappendoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len)) 315 { 316 char *t; 317 318 FIXCHAR(prefix,len,t); 319 *ierr = SNESAppendOptionsPrefix(*snes,t); 320 FREECHAR(prefix,t); 321 } 322 323 PETSC_EXTERN void PETSC_STDCALL snessetoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len)) 324 { 325 char *t; 326 327 FIXCHAR(prefix,len,t); 328 *ierr = SNESSetOptionsPrefix(*snes,t); 329 FREECHAR(prefix,t); 330 } 331 332 /*----------------------------------------------------------------------*/ 333 /* functions, hence no STDCALL */ 334 335 PETSC_EXTERN void snesmonitorlgresidualnorm_(SNES *snes,PetscInt *its,PetscReal *fgnorm,PetscObject *dummy,PetscErrorCode *ierr) 336 { 337 *ierr = SNESMonitorLGResidualNorm(*snes,*its,*fgnorm,dummy); 338 } 339 340 PETSC_EXTERN void snesmonitordefault_(SNES *snes,PetscInt *its,PetscReal *fgnorm,PetscViewerAndFormat **dummy,PetscErrorCode *ierr) 341 { 342 *ierr = SNESMonitorDefault(*snes,*its,*fgnorm,*dummy); 343 } 344 345 PETSC_EXTERN void snesmonitorsolution_(SNES *snes,PetscInt *its,PetscReal *fgnorm,PetscViewerAndFormat **dummy,PetscErrorCode *ierr) 346 { 347 *ierr = SNESMonitorSolution(*snes,*its,*fgnorm,*dummy); 348 } 349 350 PETSC_EXTERN void snesmonitorsolutionupdate_(SNES *snes,PetscInt *its,PetscReal *fgnorm,PetscViewerAndFormat **dummy,PetscErrorCode *ierr) 351 { 352 *ierr = SNESMonitorSolutionUpdate(*snes,*its,*fgnorm,*dummy); 353 } 354 355 356 PETSC_EXTERN 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) 357 { 358 CHKFORTRANNULLOBJECT(mctx); 359 if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitordefault_) { 360 *ierr = SNESMonitorSet(*snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorDefault,*(PetscViewerAndFormat**)mctx,(PetscErrorCode (*)(void **))PetscViewerAndFormatDestroy); 361 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorsolution_) { 362 *ierr = SNESMonitorSet(*snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorSolution,*(PetscViewerAndFormat**)mctx,(PetscErrorCode (*)(void **))PetscViewerAndFormatDestroy); 363 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorsolutionupdate_) { 364 *ierr = SNESMonitorSet(*snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorSolutionUpdate,*(PetscViewerAndFormat**)mctx,(PetscErrorCode (*)(void **))PetscViewerAndFormatDestroy); 365 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorlgresidualnorm_) { 366 *ierr = SNESMonitorSet(*snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorLGResidualNorm,0,0); 367 } else { 368 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.monitor,(PetscVoidFunction)func,mctx); 369 if (*ierr) return; 370 if (FORTRANNULLFUNCTION(mondestroy)) { 371 *ierr = SNESMonitorSet(*snes,oursnesmonitor,*snes,NULL); 372 } else { 373 CHKFORTRANNULLFUNCTION(mondestroy); 374 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.mondestroy,(PetscVoidFunction)mondestroy,mctx); 375 if (!*ierr) *ierr = SNESMonitorSet(*snes,oursnesmonitor,*snes,ourmondestroy); 376 } 377 } 378 } 379 380