1 #include "zpetsc.h" 2 #include "petscsnes.h" 3 4 #ifdef PETSC_HAVE_FORTRAN_UNDERSCORE_UNDERSCORE 5 #define snesconverged_tr_ snesconverged_tr__ 6 #define snesconverged_ls_ snesconverged_ls__ 7 #endif 8 9 #if defined(PETSC_HAVE_FORTRAN_CAPS) 10 #define snessolve_ SNESSOLVE 11 #define snesdefaultcomputejacobian_ SNESDEFAULTCOMPUTEJACOBIAN 12 #define snesdefaultcomputejacobiancolor_ SNESDEFAULTCOMPUTEJACOBIANCOLOR 13 #define snesdacomputejacobian_ SNESDACOMPUTEJACOBIAN 14 #define snesdacomputejacobianwithadifor_ SNESDACOMPUTEJACOBIANWITHADIFOR 15 #define snessetjacobian_ SNESSETJACOBIAN 16 #define snesgetoptionsprefix_ SNESGETOPTIONSPREFIX 17 #define snesgettype_ SNESGETTYPE 18 #define snesdaformfunction_ SNESDAFORMFUNCTION 19 #define snessetfunction_ SNESSETFUNCTION 20 #define snesgetfunction_ SNESGETFUNCTION 21 #define snessetconvergencetest_ SNESSETCONVERGENCETEST 22 #define snesdefaultconverged_ SNESDEFAULTCONVERGED 23 #define snesskipconverged_ SNESSKIPCONVERGED 24 #define snesconverged_tr_ SNESCONVERGED_TR 25 #define snesconverged_ls_ SNESCONVERGED_LS 26 #define snesview_ SNESVIEW 27 #define snesgetconvergencehistory_ SNESGETCONVERGENCEHISTORY 28 #define snesgetjacobian_ SNESGETJACOBIAN 29 #define snessettype_ SNESSETTYPE 30 #define snesappendoptionsprefix_ SNESAPPENDOPTIONSPREFIX 31 #define snessetoptionsprefix_ SNESSETOPTIONSPREFIX 32 #define snesmonitordefault_ SNESMONITORDEFAULT 33 #define snesmonitorsolution_ SNESMONITORSOLUTION 34 #define snesmonitorlg_ SNESMONITORLG 35 #define snesmonitorsolutionupdate_ SNESMONITORSOLUTIONUPDATE 36 #define snesmonitorset_ SNESMONITORSET 37 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 38 #define snessolve_ snessolve 39 #define snesdefaultcomputejacobian_ snesdefaultcomputejacobian 40 #define snesdefaultcomputejacobiancolor_ snesdefaultcomputejacobiancolor 41 #define snesdacomputejacobian_ snesdacomputejacobian 42 #define snesdacomputejacobianwithadifor_ snesdacomputejacobianwithadifor 43 #define snessetjacobian_ snessetjacobian 44 #define snesgetoptionsprefix_ snesgetoptionsprefix 45 #define snesgettype_ snesgettype 46 #define snesdaformfunction_ snesdaformfunction 47 #define snessetfunction_ snessetfunction 48 #define snesgetfunction_ snesgetfunction 49 #define snessetconvergencetest_ snessetconvergencetest 50 #define snesdefaultconverged_ snesdefaultconverged 51 #define snesskipconverged_ snesskipconverged 52 #define snesconverged_tr_ snesconverged_tr 53 #define snesconverged_ls_ snesconverged_ls 54 #define snesview_ snesview 55 #define snesgetjacobian_ snesgetjacobian 56 #define snesgetconvergencehistory_ snesgetconvergencehistory 57 #define snessettype_ snessettype 58 #define snesappendoptionsprefix_ snesappendoptionsprefix 59 #define snessetoptionsprefix_ snessetoptionsprefix 60 #define snesmonitorlg_ snesmonitorlg 61 #define snesmonitordefault_ snesmonitordefault 62 #define snesmonitorsolution_ snesmonitorsolution 63 #define snesmonitorsolutionupdate_ snesmonitorsolutionupdate 64 #define snesmonitorset_ snesmonitorset 65 #endif 66 67 EXTERN_C_BEGIN 68 static void (PETSC_STDCALL *f3)(SNES*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*); 69 static void (PETSC_STDCALL *f2)(SNES*,Vec*,Vec*,void*,PetscErrorCode*); 70 static void (PETSC_STDCALL *f8)(SNES*,PetscInt *,PetscReal*,PetscReal*,PetscReal*,SNESConvergedReason*,void*,PetscErrorCode*); 71 static void (PETSC_STDCALL *f7)(SNES*,PetscInt*,PetscReal*,void*,PetscErrorCode*); 72 static void (PETSC_STDCALL *f71)(void*,PetscErrorCode*); 73 EXTERN_C_END 74 75 static PetscErrorCode oursnesfunction(SNES snes,Vec x,Vec f,void *ctx) 76 { 77 PetscErrorCode ierr = 0; 78 (*f2)(&snes,&x,&f,ctx,&ierr);CHKERRQ(ierr); 79 return 0; 80 } 81 static PetscErrorCode oursnestest(SNES snes,PetscInt it,PetscReal a,PetscReal d,PetscReal c,SNESConvergedReason*reason,void*ctx) 82 { 83 PetscErrorCode ierr = 0; 84 85 (*f8)(&snes,&it,&a,&d,&c,reason,ctx,&ierr);CHKERRQ(ierr); 86 return 0; 87 } 88 89 static PetscErrorCode oursnesjacobian(SNES snes,Vec x,Mat* m,Mat* p,MatStructure* type,void*ctx) 90 { 91 PetscErrorCode ierr = 0; 92 (*f3)(&snes,&x,m,p,type,ctx,&ierr);CHKERRQ(ierr); 93 return 0; 94 } 95 static PetscErrorCode oursnesmonitor(SNES snes,PetscInt i,PetscReal d,void*ctx) 96 { 97 PetscErrorCode ierr = 0; 98 99 (*f7)(&snes,&i,&d,ctx,&ierr);CHKERRQ(ierr); 100 return 0; 101 } 102 static PetscErrorCode ourmondestroy(void* ctx) 103 { 104 PetscErrorCode ierr = 0; 105 106 (*f71)(ctx,&ierr);CHKERRQ(ierr); 107 return 0; 108 } 109 110 EXTERN_C_BEGIN 111 /* ---------------------------------------------------------*/ 112 /* 113 snesdefaultcomputejacobian() and snesdefaultcomputejacobiancolor() 114 These can be used directly from Fortran but are mostly so that 115 Fortran SNESSetJacobian() will properly handle the defaults being passed in. 116 117 functions, hence no STDCALL 118 */ 119 void snesdefaultcomputejacobian_(SNES *snes,Vec *x,Mat *m,Mat *p,MatStructure* type,void *ctx,PetscErrorCode *ierr) 120 { 121 *ierr = SNESDefaultComputeJacobian(*snes,*x,m,p,type,ctx); 122 } 123 void snesdefaultcomputejacobiancolor_(SNES *snes,Vec *x,Mat *m,Mat *p,MatStructure* type,void *ctx,PetscErrorCode *ierr) 124 { 125 *ierr = SNESDefaultComputeJacobianColor(*snes,*x,m,p,type,*(MatFDColoring*)ctx); 126 } 127 128 void snesdacomputejacobianwithadifor_(SNES *snes,Vec *X,Mat *m,Mat *p,MatStructure* type,void *ctx,PetscErrorCode *ierr) 129 { 130 (*PetscErrorPrintf)("Cannot call this function from Fortran"); 131 *ierr = 1; 132 } 133 134 void snesdacomputejacobian_(SNES *snes,Vec *X,Mat *m,Mat *p,MatStructure* type,void *ctx,PetscErrorCode *ierr) 135 { 136 (*PetscErrorPrintf)("Cannot call this function from Fortran"); 137 *ierr = 1; 138 } 139 140 void PETSC_STDCALL snessetjacobian_(SNES *snes,Mat *A,Mat *B,void (PETSC_STDCALL *func)(SNES*,Vec*,Mat*,Mat*, 141 MatStructure*,void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr) 142 { 143 CHKFORTRANNULLOBJECT(ctx); 144 if ((PetscVoidFunction)func == (PetscVoidFunction)snesdefaultcomputejacobian_) { 145 *ierr = SNESSetJacobian(*snes,*A,*B,SNESDefaultComputeJacobian,ctx); 146 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesdefaultcomputejacobiancolor_) { 147 *ierr = SNESSetJacobian(*snes,*A,*B,SNESDefaultComputeJacobianColor,*(MatFDColoring*)ctx); 148 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesdacomputejacobianwithadifor_) { 149 *ierr = SNESSetJacobian(*snes,*A,*B,SNESDAComputeJacobianWithAdifor,ctx); 150 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesdacomputejacobian_) { 151 *ierr = SNESSetJacobian(*snes,*A,*B,SNESDAComputeJacobian,ctx); 152 } else { 153 f3 = func; 154 *ierr = SNESSetJacobian(*snes,*A,*B,oursnesjacobian,ctx); 155 } 156 } 157 /* -------------------------------------------------------------*/ 158 159 void PETSC_STDCALL snessolve_(SNES *snes,Vec *b,Vec *x, int *__ierr ) 160 { 161 Vec B = *b; 162 if (*b == PETSC_NULL_OBJECT_Fortran) B = PETSC_NULL; 163 *__ierr = SNESSolve(*snes,B,*x); 164 } 165 166 void PETSC_STDCALL snesgetoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len), 167 PetscErrorCode *ierr PETSC_END_LEN(len)) 168 { 169 const char *tname; 170 171 *ierr = SNESGetOptionsPrefix(*snes,&tname); 172 #if defined(PETSC_USES_CPTOFCD) 173 { 174 char *t = _fcdtocp(prefix); int len1 = _fcdlen(prefix); 175 *ierr = PetscStrncpy(t,tname,len1);if (*ierr) return; 176 } 177 #else 178 *ierr = PetscStrncpy(prefix,tname,len);if (*ierr) return; 179 #endif 180 } 181 182 void PETSC_STDCALL snesgettype_(SNES *snes,CHAR name PETSC_MIXED_LEN(len), 183 PetscErrorCode *ierr PETSC_END_LEN(len)) 184 { 185 const char *tname; 186 187 *ierr = SNESGetType(*snes,&tname); 188 #if defined(PETSC_USES_CPTOFCD) 189 { 190 char *t = _fcdtocp(name); int len1 = _fcdlen(name); 191 *ierr = PetscStrncpy(t,tname,len1);if (*ierr) return; 192 } 193 #else 194 *ierr = PetscStrncpy(name,tname,len);if (*ierr) return; 195 #endif 196 FIXRETURNCHAR(name,len); 197 } 198 /* ---------------------------------------------------------*/ 199 200 /* 201 These are not usually called from Fortran but allow Fortran users 202 to transparently set these monitors from .F code 203 204 functions, hence no STDCALL 205 */ 206 void snesdaformfunction_(SNES *snes,Vec *X, Vec *F,void *ptr,PetscErrorCode *ierr) 207 { 208 *ierr = SNESDAFormFunction(*snes,*X,*F,ptr); 209 } 210 211 void PETSC_STDCALL snessetfunction_(SNES *snes,Vec *r,void (PETSC_STDCALL *func)(SNES*,Vec*,Vec*,void*,PetscErrorCode*), 212 void *ctx,PetscErrorCode *ierr) 213 { 214 CHKFORTRANNULLOBJECT(ctx); 215 f2 = func; 216 if ((PetscVoidFunction)func == (PetscVoidFunction)snesdaformfunction_) { 217 *ierr = SNESSetFunction(*snes,*r,SNESDAFormFunction,ctx); 218 } else { 219 *ierr = SNESSetFunction(*snes,*r,oursnesfunction,ctx); 220 } 221 } 222 /* ---------------------------------------------------------*/ 223 224 /* the func argument is ignored */ 225 void PETSC_STDCALL snesgetfunction_(SNES *snes,Vec *r,void *func,void **ctx,PetscErrorCode *ierr) 226 { 227 CHKFORTRANNULLINTEGER(ctx); 228 CHKFORTRANNULLOBJECT(r); 229 *ierr = SNESGetFunction(*snes,r,PETSC_NULL,ctx); 230 } 231 /*----------------------------------------------------------------------*/ 232 233 void snesdefaultconverged_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r, 234 void *ct,PetscErrorCode *ierr) 235 { 236 *ierr = SNESDefaultConverged(*snes,*it,*a,*b,*c,r,ct); 237 } 238 239 void snesskipconverged_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r, 240 void *ct,PetscErrorCode *ierr) 241 { 242 *ierr = SNESSkipConverged(*snes,*it,*a,*b,*c,r,ct); 243 } 244 245 void snesconverged_tr_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r, 246 void *ct,PetscErrorCode *ierr) 247 { 248 *ierr = SNESConverged_TR(*snes,*it,*a,*b,*c,r,ct); 249 } 250 251 void snesconverged_ls_(SNES *snes,PetscInt *it, PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r, 252 void *ct,PetscErrorCode *ierr) 253 { 254 *ierr = SNESConverged_LS(*snes,*it,*a,*b,*c,r,ct); 255 } 256 257 258 void PETSC_STDCALL snessetconvergencetest_(SNES *snes, 259 void (PETSC_STDCALL *func)(SNES*,PetscInt*,PetscReal*,PetscReal*,PetscReal*,SNESConvergedReason*,void*,PetscErrorCode*), 260 void *cctx,PetscErrorCode *ierr) 261 { 262 CHKFORTRANNULLOBJECT(cctx); 263 264 if ((PetscVoidFunction)func == (PetscVoidFunction)snesdefaultconverged_){ 265 *ierr = SNESSetConvergenceTest(*snes,SNESDefaultConverged,0); 266 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesskipconverged_){ 267 *ierr = SNESSetConvergenceTest(*snes,SNESSkipConverged,0); 268 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesconverged_ls_){ 269 *ierr = SNESSetConvergenceTest(*snes,SNESConverged_LS,0); 270 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesconverged_tr_){ 271 *ierr = SNESSetConvergenceTest(*snes,SNESConverged_TR,0); 272 } else { 273 f8 = func; 274 *ierr = SNESSetConvergenceTest(*snes,oursnestest,cctx); 275 } 276 } 277 /*----------------------------------------------------------------------*/ 278 279 void PETSC_STDCALL snesview_(SNES *snes,PetscViewer *viewer, PetscErrorCode *ierr) 280 { 281 PetscViewer v; 282 PetscPatchDefaultViewers_Fortran(viewer,v); 283 *ierr = SNESView(*snes,v); 284 } 285 286 /* func is currently ignored from Fortran */ 287 void PETSC_STDCALL snesgetjacobian_(SNES *snes,Mat *A,Mat *B,int *func,void **ctx,PetscErrorCode *ierr) 288 { 289 CHKFORTRANNULLINTEGER(ctx); 290 CHKFORTRANNULLOBJECT(A); 291 CHKFORTRANNULLOBJECT(B); 292 *ierr = SNESGetJacobian(*snes,A,B,0,ctx); 293 } 294 295 void PETSC_STDCALL snesgetconvergencehistory_(SNES *snes,PetscInt *na,PetscErrorCode *ierr) 296 { 297 *ierr = SNESGetConvergenceHistory(*snes,PETSC_NULL,PETSC_NULL,na); 298 } 299 300 void PETSC_STDCALL snessettype_(SNES *snes,CHAR type PETSC_MIXED_LEN(len), 301 PetscErrorCode *ierr PETSC_END_LEN(len)) 302 { 303 char *t; 304 305 FIXCHAR(type,len,t); 306 *ierr = SNESSetType(*snes,t); 307 FREECHAR(type,t); 308 } 309 310 void PETSC_STDCALL snesappendoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len), 311 PetscErrorCode *ierr PETSC_END_LEN(len)) 312 { 313 char *t; 314 315 FIXCHAR(prefix,len,t); 316 *ierr = SNESAppendOptionsPrefix(*snes,t); 317 FREECHAR(prefix,t); 318 } 319 320 void PETSC_STDCALL snessetoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len), 321 PetscErrorCode *ierr PETSC_END_LEN(len)) 322 { 323 char *t; 324 325 FIXCHAR(prefix,len,t); 326 *ierr = SNESSetOptionsPrefix(*snes,t); 327 FREECHAR(prefix,t); 328 } 329 330 /*----------------------------------------------------------------------*/ 331 /* functions, hence no STDCALL */ 332 333 void snesmonitorlg_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr) 334 { 335 *ierr = SNESMonitorLG(*snes,*its,*fgnorm,dummy); 336 } 337 338 void snesmonitordefault_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr) 339 { 340 *ierr = SNESMonitorDefault(*snes,*its,*fgnorm,dummy); 341 } 342 343 void snesmonitorsolution_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr) 344 { 345 *ierr = SNESMonitorSolution(*snes,*its,*fgnorm,dummy); 346 } 347 348 void snesmonitorsolutionupdate_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr) 349 { 350 *ierr = SNESMonitorSolutionUpdate(*snes,*its,*fgnorm,dummy); 351 } 352 353 354 void PETSC_STDCALL snesmonitorset_(SNES *snes,void (PETSC_STDCALL *func)(SNES*,PetscInt*,PetscReal*,void*,PetscErrorCode*), 355 void *mctx,void (PETSC_STDCALL *mondestroy)(void*,PetscErrorCode*),PetscErrorCode *ierr) 356 { 357 CHKFORTRANNULLOBJECT(mctx); 358 if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitordefault_) { 359 *ierr = SNESMonitorSet(*snes,SNESMonitorDefault,0,0); 360 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorsolution_) { 361 *ierr = SNESMonitorSet(*snes,SNESMonitorSolution,0,0); 362 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorsolutionupdate_) { 363 *ierr = SNESMonitorSet(*snes,SNESMonitorSolutionUpdate,0,0); 364 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorlg_) { 365 *ierr = SNESMonitorSet(*snes,SNESMonitorLG,0,0); 366 } else { 367 f7 = func; 368 if (FORTRANNULLFUNCTION(mondestroy)){ 369 *ierr = SNESMonitorSet(*snes,oursnesmonitor,mctx,0); 370 } else { 371 f71 = mondestroy; 372 *ierr = SNESMonitorSet(*snes,oursnesmonitor,mctx,ourmondestroy); 373 } 374 } 375 } 376 377 378 379 EXTERN_C_END 380