1 #include "private/zpetsc.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 EXTERN_C_BEGIN 61 static void (PETSC_STDCALL *f3)(SNES*,Vec*,Mat*,Mat*,MatStructure*,void*,PetscErrorCode*); 62 static void (PETSC_STDCALL *f2)(SNES*,Vec*,Vec*,void*,PetscErrorCode*); 63 static void (PETSC_STDCALL *f8)(SNES*,PetscInt *,PetscReal*,PetscReal*,PetscReal*,SNESConvergedReason*,void*,PetscErrorCode*); 64 static void (PETSC_STDCALL *f7)(SNES*,PetscInt*,PetscReal*,void*,PetscErrorCode*); 65 static void (PETSC_STDCALL *f71)(void*,PetscErrorCode*); 66 EXTERN_C_END 67 68 static PetscErrorCode oursnesfunction(SNES snes,Vec x,Vec f,void *ctx) 69 { 70 PetscErrorCode ierr = 0; 71 (*f2)(&snes,&x,&f,ctx,&ierr);CHKERRQ(ierr); 72 return 0; 73 } 74 static PetscErrorCode oursnestest(SNES snes,PetscInt it,PetscReal a,PetscReal d,PetscReal c,SNESConvergedReason*reason,void*ctx) 75 { 76 PetscErrorCode ierr = 0; 77 78 (*f8)(&snes,&it,&a,&d,&c,reason,ctx,&ierr);CHKERRQ(ierr); 79 return 0; 80 } 81 82 static PetscErrorCode oursnesjacobian(SNES snes,Vec x,Mat* m,Mat* p,MatStructure* type,void*ctx) 83 { 84 PetscErrorCode ierr = 0; 85 (*f3)(&snes,&x,m,p,type,ctx,&ierr);CHKERRQ(ierr); 86 return 0; 87 } 88 static PetscErrorCode oursnesmonitor(SNES snes,PetscInt i,PetscReal d,void*ctx) 89 { 90 PetscErrorCode ierr = 0; 91 92 (*f7)(&snes,&i,&d,ctx,&ierr);CHKERRQ(ierr); 93 return 0; 94 } 95 static PetscErrorCode ourmondestroy(void* ctx) 96 { 97 PetscErrorCode ierr = 0; 98 99 (*f71)(ctx,&ierr);CHKERRQ(ierr); 100 return 0; 101 } 102 103 EXTERN_C_BEGIN 104 /* ---------------------------------------------------------*/ 105 /* 106 snesdefaultcomputejacobian() and snesdefaultcomputejacobiancolor() 107 These can be used directly from Fortran but are mostly so that 108 Fortran SNESSetJacobian() will properly handle the defaults being passed in. 109 110 functions, hence no STDCALL 111 */ 112 void snesdefaultcomputejacobian_(SNES *snes,Vec *x,Mat *m,Mat *p,MatStructure* type,void *ctx,PetscErrorCode *ierr) 113 { 114 *ierr = SNESDefaultComputeJacobian(*snes,*x,m,p,type,ctx); 115 } 116 void snesdefaultcomputejacobiancolor_(SNES *snes,Vec *x,Mat *m,Mat *p,MatStructure* type,void *ctx,PetscErrorCode *ierr) 117 { 118 *ierr = SNESDefaultComputeJacobianColor(*snes,*x,m,p,type,*(MatFDColoring*)ctx); 119 } 120 121 void snesdacomputejacobianwithadifor_(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 snesdacomputejacobian_(SNES *snes,Vec *X,Mat *m,Mat *p,MatStructure* type,void *ctx,PetscErrorCode *ierr) 128 { 129 (*PetscErrorPrintf)("Cannot call this function from Fortran"); 130 *ierr = 1; 131 } 132 133 void PETSC_STDCALL snessetjacobian_(SNES *snes,Mat *A,Mat *B,void (PETSC_STDCALL *func)(SNES*,Vec*,Mat*,Mat*, 134 MatStructure*,void*,PetscErrorCode*),void *ctx,PetscErrorCode *ierr) 135 { 136 CHKFORTRANNULLOBJECT(ctx); 137 if ((PetscVoidFunction)func == (PetscVoidFunction)snesdefaultcomputejacobian_) { 138 *ierr = SNESSetJacobian(*snes,*A,*B,SNESDefaultComputeJacobian,ctx); 139 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesdefaultcomputejacobiancolor_) { 140 *ierr = SNESSetJacobian(*snes,*A,*B,SNESDefaultComputeJacobianColor,*(MatFDColoring*)ctx); 141 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesdacomputejacobianwithadifor_) { 142 *ierr = SNESSetJacobian(*snes,*A,*B,SNESDAComputeJacobianWithAdifor,ctx); 143 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesdacomputejacobian_) { 144 *ierr = SNESSetJacobian(*snes,*A,*B,SNESDAComputeJacobian,ctx); 145 } else { 146 f3 = func; 147 *ierr = SNESSetJacobian(*snes,*A,*B,oursnesjacobian,ctx); 148 } 149 } 150 /* -------------------------------------------------------------*/ 151 152 void PETSC_STDCALL snessolve_(SNES *snes,Vec *b,Vec *x, int *__ierr ) 153 { 154 Vec B = *b; 155 if (*b == PETSC_NULL_OBJECT_Fortran) B = PETSC_NULL; 156 *__ierr = SNESSolve(*snes,B,*x); 157 } 158 159 void PETSC_STDCALL snesgetoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len), 160 PetscErrorCode *ierr PETSC_END_LEN(len)) 161 { 162 const char *tname; 163 164 *ierr = SNESGetOptionsPrefix(*snes,&tname); 165 #if defined(PETSC_USES_CPTOFCD) 166 { 167 char *t = _fcdtocp(prefix); int len1 = _fcdlen(prefix); 168 *ierr = PetscStrncpy(t,tname,len1);if (*ierr) return; 169 } 170 #else 171 *ierr = PetscStrncpy(prefix,tname,len);if (*ierr) return; 172 #endif 173 } 174 175 void PETSC_STDCALL snesgettype_(SNES *snes,CHAR name PETSC_MIXED_LEN(len), 176 PetscErrorCode *ierr PETSC_END_LEN(len)) 177 { 178 const char *tname; 179 180 *ierr = SNESGetType(*snes,&tname); 181 #if defined(PETSC_USES_CPTOFCD) 182 { 183 char *t = _fcdtocp(name); int len1 = _fcdlen(name); 184 *ierr = PetscStrncpy(t,tname,len1);if (*ierr) return; 185 } 186 #else 187 *ierr = PetscStrncpy(name,tname,len);if (*ierr) return; 188 #endif 189 FIXRETURNCHAR(PETSC_TRUE,name,len); 190 } 191 192 void PETSC_STDCALL snesgetapplicationcontext_(SNES *snes,void **ctx,PetscErrorCode *ierr) 193 { 194 *ierr = SNESGetApplicationContext(*snes,ctx); 195 } 196 /* ---------------------------------------------------------*/ 197 198 /* 199 These are not usually called from Fortran but allow Fortran users 200 to transparently set these monitors from .F code 201 202 functions, hence no STDCALL 203 */ 204 void snesdaformfunction_(SNES *snes,Vec *X, Vec *F,void *ptr,PetscErrorCode *ierr) 205 { 206 *ierr = SNESDAFormFunction(*snes,*X,*F,ptr); 207 } 208 209 void PETSC_STDCALL snessetfunction_(SNES *snes,Vec *r,void (PETSC_STDCALL *func)(SNES*,Vec*,Vec*,void*,PetscErrorCode*), 210 void *ctx,PetscErrorCode *ierr) 211 { 212 CHKFORTRANNULLOBJECT(ctx); 213 f2 = func; 214 if ((PetscVoidFunction)func == (PetscVoidFunction)snesdaformfunction_) { 215 *ierr = SNESSetFunction(*snes,*r,SNESDAFormFunction,ctx); 216 } else { 217 *ierr = SNESSetFunction(*snes,*r,oursnesfunction,ctx); 218 } 219 } 220 /* ---------------------------------------------------------*/ 221 222 /* the func argument is ignored */ 223 void PETSC_STDCALL snesgetfunction_(SNES *snes,Vec *r,void *func,void **ctx,PetscErrorCode *ierr) 224 { 225 CHKFORTRANNULLINTEGER(ctx); 226 CHKFORTRANNULLOBJECT(r); 227 *ierr = SNESGetFunction(*snes,r,PETSC_NULL,ctx); 228 } 229 /*----------------------------------------------------------------------*/ 230 231 void snesdefaultconverged_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r, 232 void *ct,PetscErrorCode *ierr) 233 { 234 *ierr = SNESDefaultConverged(*snes,*it,*a,*b,*c,r,ct); 235 } 236 237 void snesskipconverged_(SNES *snes,PetscInt *it,PetscReal *a,PetscReal *b,PetscReal *c,SNESConvergedReason *r, 238 void *ct,PetscErrorCode *ierr) 239 { 240 *ierr = SNESSkipConverged(*snes,*it,*a,*b,*c,r,ct); 241 } 242 243 void PETSC_STDCALL snessetconvergencetest_(SNES *snes, 244 void (PETSC_STDCALL *func)(SNES*,PetscInt*,PetscReal*,PetscReal*,PetscReal*,SNESConvergedReason*,void*,PetscErrorCode*), 245 void *cctx,PetscErrorCode *ierr) 246 { 247 CHKFORTRANNULLOBJECT(cctx); 248 249 if ((PetscVoidFunction)func == (PetscVoidFunction)snesdefaultconverged_){ 250 *ierr = SNESSetConvergenceTest(*snes,SNESDefaultConverged,0); 251 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesskipconverged_){ 252 *ierr = SNESSetConvergenceTest(*snes,SNESSkipConverged,0); 253 } else { 254 f8 = func; 255 *ierr = SNESSetConvergenceTest(*snes,oursnestest,cctx); 256 } 257 } 258 /*----------------------------------------------------------------------*/ 259 260 void PETSC_STDCALL snesview_(SNES *snes,PetscViewer *viewer, PetscErrorCode *ierr) 261 { 262 PetscViewer v; 263 PetscPatchDefaultViewers_Fortran(viewer,v); 264 *ierr = SNESView(*snes,v); 265 } 266 267 /* func is currently ignored from Fortran */ 268 void PETSC_STDCALL snesgetjacobian_(SNES *snes,Mat *A,Mat *B,int *func,void **ctx,PetscErrorCode *ierr) 269 { 270 CHKFORTRANNULLINTEGER(ctx); 271 CHKFORTRANNULLOBJECT(A); 272 CHKFORTRANNULLOBJECT(B); 273 *ierr = SNESGetJacobian(*snes,A,B,0,ctx); 274 } 275 276 void PETSC_STDCALL snesgetconvergencehistory_(SNES *snes,PetscInt *na,PetscErrorCode *ierr) 277 { 278 *ierr = SNESGetConvergenceHistory(*snes,PETSC_NULL,PETSC_NULL,na); 279 } 280 281 void PETSC_STDCALL snessettype_(SNES *snes,CHAR type PETSC_MIXED_LEN(len), 282 PetscErrorCode *ierr PETSC_END_LEN(len)) 283 { 284 char *t; 285 286 FIXCHAR(type,len,t); 287 *ierr = SNESSetType(*snes,t); 288 FREECHAR(type,t); 289 } 290 291 void PETSC_STDCALL snesappendoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len), 292 PetscErrorCode *ierr PETSC_END_LEN(len)) 293 { 294 char *t; 295 296 FIXCHAR(prefix,len,t); 297 *ierr = SNESAppendOptionsPrefix(*snes,t); 298 FREECHAR(prefix,t); 299 } 300 301 void PETSC_STDCALL snessetoptionsprefix_(SNES *snes,CHAR prefix PETSC_MIXED_LEN(len), 302 PetscErrorCode *ierr PETSC_END_LEN(len)) 303 { 304 char *t; 305 306 FIXCHAR(prefix,len,t); 307 *ierr = SNESSetOptionsPrefix(*snes,t); 308 FREECHAR(prefix,t); 309 } 310 311 /*----------------------------------------------------------------------*/ 312 /* functions, hence no STDCALL */ 313 314 void snesmonitorlg_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr) 315 { 316 *ierr = SNESMonitorLG(*snes,*its,*fgnorm,dummy); 317 } 318 319 void snesmonitordefault_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr) 320 { 321 *ierr = SNESMonitorDefault(*snes,*its,*fgnorm,dummy); 322 } 323 324 void snesmonitorsolution_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr) 325 { 326 *ierr = SNESMonitorSolution(*snes,*its,*fgnorm,dummy); 327 } 328 329 void snesmonitorsolutionupdate_(SNES *snes,PetscInt *its,PetscReal *fgnorm,void *dummy,PetscErrorCode *ierr) 330 { 331 *ierr = SNESMonitorSolutionUpdate(*snes,*its,*fgnorm,dummy); 332 } 333 334 335 void PETSC_STDCALL snesmonitorset_(SNES *snes,void (PETSC_STDCALL *func)(SNES*,PetscInt*,PetscReal*,void*,PetscErrorCode*), 336 void *mctx,void (PETSC_STDCALL *mondestroy)(void*,PetscErrorCode*),PetscErrorCode *ierr) 337 { 338 CHKFORTRANNULLOBJECT(mctx); 339 if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitordefault_) { 340 *ierr = SNESMonitorSet(*snes,SNESMonitorDefault,0,0); 341 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorsolution_) { 342 *ierr = SNESMonitorSet(*snes,SNESMonitorSolution,0,0); 343 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorsolutionupdate_) { 344 *ierr = SNESMonitorSet(*snes,SNESMonitorSolutionUpdate,0,0); 345 } else if ((PetscVoidFunction)func == (PetscVoidFunction)snesmonitorlg_) { 346 *ierr = SNESMonitorSet(*snes,SNESMonitorLG,0,0); 347 } else { 348 f7 = func; 349 if (FORTRANNULLFUNCTION(mondestroy)){ 350 *ierr = SNESMonitorSet(*snes,oursnesmonitor,mctx,0); 351 } else { 352 f71 = mondestroy; 353 *ierr = SNESMonitorSet(*snes,oursnesmonitor,mctx,ourmondestroy); 354 } 355 } 356 } 357 358 359 360 EXTERN_C_END 361