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 CHKFORTRANNULLOBJECT(ctx); 160 CHKFORTRANNULLFUNCTION(func); 161 if ((PetscVoidFunction)func == (PetscVoidFunction)snescomputejacobiandefault_) { 162 *ierr = SNESSetJacobian(*snes,*A,*B,SNESComputeJacobianDefault,ctx); 163 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snescomputejacobiandefaultcolor_) { 164 *ierr = SNESSetJacobian(*snes,*A,*B,SNESComputeJacobianDefaultColor,*(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 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.jacobian,(PetscVoidFunction)func,ctx); 171 if (!*ierr) *ierr = SNESSetJacobian(*snes,*A,*B,oursnesjacobian,NULL); 172 } 173 } 174 /* -------------------------------------------------------------*/ 175 176 PETSC_EXTERN void PETSC_STDCALL snessolve_(SNES *snes,Vec *b,Vec *x, int *__ierr) 177 { 178 Vec B = *b,X = *x; 179 if (FORTRANNULLOBJECT(b)) B = NULL; 180 if (FORTRANNULLOBJECT(x)) X = NULL; 181 *__ierr = SNESSolve(*snes,B,X); 182 } 183 184 PETSC_EXTERN 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 PETSC_EXTERN 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 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)) 211 { 212 CHKFORTRANNULLOBJECT(ctx); 213 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.function,(PetscVoidFunction)func,ctx); 214 #if defined(PETSC_HAVE_F90_2PTR_ARG) 215 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.function_pgiptr,PETSC_NULL,ptr); 216 #endif 217 if (!*ierr) *ierr = SNESSetFunction(*snes,*r,oursnesfunction,NULL); 218 } 219 220 221 PETSC_EXTERN void PETSC_STDCALL snessetngs_(SNES *snes,void (PETSC_STDCALL *func)(SNES*,Vec*,Vec*,void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr) 222 { 223 CHKFORTRANNULLOBJECT(ctx); 224 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.ngs,(PetscVoidFunction)func,ctx); 225 if (!*ierr) *ierr = SNESSetNGS(*snes,oursnesngs,NULL); 226 } 227 PETSC_EXTERN void PETSC_STDCALL snessetupdate_(SNES *snes,void (PETSC_STDCALL *func)(SNES*,PetscInt*,PetscErrorCode*),PetscErrorCode *ierr) 228 { 229 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.update,(PetscVoidFunction)func,NULL); 230 if (!*ierr) *ierr = SNESSetUpdate(*snes,oursnesupdate); 231 } 232 /* ---------------------------------------------------------*/ 233 234 /* the func argument is ignored */ 235 PETSC_EXTERN void PETSC_STDCALL snesgetfunction_(SNES *snes,Vec *r,void *func,void **ctx,PetscErrorCode *ierr) 236 { 237 CHKFORTRANNULLINTEGER(ctx); 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 CHKFORTRANNULLINTEGER(ctx); 246 *ierr = PetscObjectGetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,_cb.ngs,NULL,ctx); 247 } 248 249 /*----------------------------------------------------------------------*/ 250 251 PETSC_EXTERN void snesconvergeddefault_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r, void *ct,PetscErrorCode *ierr) 252 { 253 *ierr = SNESConvergedDefault(*snes,*it,*a,*b,*c,r,ct); 254 } 255 256 PETSC_EXTERN void snesconvergedskip_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r,void *ct,PetscErrorCode *ierr) 257 { 258 *ierr = SNESConvergedSkip(*snes,*it,*a,*b,*c,r,ct); 259 } 260 261 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) 262 { 263 CHKFORTRANNULLOBJECT(cctx); 264 CHKFORTRANNULLFUNCTION(destroy); 265 266 if ((PetscVoidFunction)func == (PetscVoidFunction)snesconvergeddefault_) { 267 *ierr = SNESSetConvergenceTest(*snes,SNESConvergedDefault,0,0); 268 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesconvergedskip_) { 269 *ierr = SNESSetConvergenceTest(*snes,SNESConvergedSkip,0,0); 270 } else { 271 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.test,(PetscVoidFunction)func,cctx); 272 if (*ierr) return; 273 if (!destroy) { 274 *ierr = SNESSetConvergenceTest(*snes,oursnestest,*snes,NULL); 275 } else { 276 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.destroy,(PetscVoidFunction)destroy,cctx); 277 if (!*ierr) *ierr = SNESSetConvergenceTest(*snes,oursnestest,*snes,ourdestroy); 278 } 279 } 280 } 281 /*----------------------------------------------------------------------*/ 282 283 PETSC_EXTERN void PETSC_STDCALL snesview_(SNES *snes,PetscViewer *viewer, PetscErrorCode *ierr) 284 { 285 PetscViewer v; 286 PetscPatchDefaultViewers_Fortran(viewer,v); 287 *ierr = SNESView(*snes,v); 288 } 289 290 /* func is currently ignored from Fortran */ 291 PETSC_EXTERN void PETSC_STDCALL snesgetjacobian_(SNES *snes,Mat *A,Mat *B,int *func,void **ctx,PetscErrorCode *ierr) 292 { 293 CHKFORTRANNULLINTEGER(ctx); 294 CHKFORTRANNULLOBJECT(A); 295 CHKFORTRANNULLOBJECT(B); 296 *ierr = SNESGetJacobian(*snes,A,B,0,NULL); if (*ierr) return; 297 *ierr = PetscObjectGetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,_cb.jacobian,NULL,ctx); 298 299 } 300 301 PETSC_EXTERN void PETSC_STDCALL snesgetconvergencehistory_(SNES *snes,PetscInt *na,PetscErrorCode *ierr) 302 { 303 *ierr = SNESGetConvergenceHistory(*snes,NULL,NULL,na); 304 } 305 306 PETSC_EXTERN void PETSC_STDCALL snessettype_(SNES *snes,CHAR type PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len)) 307 { 308 char *t; 309 310 FIXCHAR(type,len,t); 311 *ierr = SNESSetType(*snes,t); 312 FREECHAR(type,t); 313 } 314 315 PETSC_EXTERN void PETSC_STDCALL snesappendoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len)) 316 { 317 char *t; 318 319 FIXCHAR(prefix,len,t); 320 *ierr = SNESAppendOptionsPrefix(*snes,t); 321 FREECHAR(prefix,t); 322 } 323 324 PETSC_EXTERN void PETSC_STDCALL snessetoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len)) 325 { 326 char *t; 327 328 FIXCHAR(prefix,len,t); 329 *ierr = SNESSetOptionsPrefix(*snes,t); 330 FREECHAR(prefix,t); 331 } 332 333 /*----------------------------------------------------------------------*/ 334 /* functions, hence no STDCALL */ 335 336 PETSC_EXTERN void snesmonitorlgresidualnorm_(SNES *snes,PetscInt *its,PetscReal *fgnorm,PetscObject *dummy,PetscErrorCode *ierr) 337 { 338 *ierr = SNESMonitorLGResidualNorm(*snes,*its,*fgnorm,dummy); 339 } 340 341 PETSC_EXTERN void snesmonitordefault_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr) 342 { 343 *ierr = SNESMonitorDefault(*snes,*its,*fgnorm,dummy); 344 } 345 346 PETSC_EXTERN void snesmonitorsolution_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr) 347 { 348 *ierr = SNESMonitorSolution(*snes,*its,*fgnorm,dummy); 349 } 350 351 PETSC_EXTERN void snesmonitorsolutionupdate_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr) 352 { 353 *ierr = SNESMonitorSolutionUpdate(*snes,*its,*fgnorm,dummy); 354 } 355 356 357 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) 358 { 359 CHKFORTRANNULLOBJECT(mctx); 360 if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitordefault_) { 361 *ierr = SNESMonitorSet(*snes,SNESMonitorDefault,0,0); 362 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorsolution_) { 363 *ierr = SNESMonitorSet(*snes,SNESMonitorSolution,0,0); 364 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorsolutionupdate_) { 365 *ierr = SNESMonitorSet(*snes,SNESMonitorSolutionUpdate,0,0); 366 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorlgresidualnorm_) { 367 *ierr = SNESMonitorSet(*snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorLGResidualNorm,0,0); 368 } else { 369 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.monitor,(PetscVoidFunction)func,mctx); 370 if (*ierr) return; 371 if (FORTRANNULLFUNCTION(mondestroy)) { 372 *ierr = SNESMonitorSet(*snes,oursnesmonitor,*snes,NULL); 373 } else { 374 CHKFORTRANNULLFUNCTION(mondestroy); 375 *ierr = PetscObjectSetFortranCallback((PetscObject)*snes,PETSC_FORTRAN_CALLBACK_CLASS,&_cb.mondestroy,(PetscVoidFunction)mondestroy,mctx); 376 if (!*ierr) *ierr = SNESMonitorSet(*snes,oursnesmonitor,*snes,ourmondestroy); 377 } 378 } 379 } 380 381