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