1 #include <petsc/private/fortranimpl.h> 2 #include <petsc/private/taolinesearchimpl.h> 3 4 #if defined(PETSC_HAVE_FORTRAN_CAPS) 5 #define taolinesearchsetobjectiveroutine_ TAOLINESEARCHSETOBJECTIVEROUTINE 6 #define taolinesearchsetgradientroutine_ TAOLINESEARCHSETGRADIENTROUTINE 7 #define taolinesearchsetobjectiveandgradientroutine_ TAOLINESEARCHSETOBJECTIVEANDGRADIENTROUTINE 8 #define taolinesearchsetobjectiveandgtsroutine_ TAOLINESEARCHSETOBJECTIVEANDGTSROUTINE 9 #define taolinesearchview_ TAOLINESEARCHVIEW 10 #define taolinesearchsettype_ TAOLINESEARCHSETTYPE 11 #define taolinesearchviewfromoptions_ TAOLINESEARCHVIEWFROMOPTIONS 12 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 13 14 #define taolinesearchsetobjectiveroutine_ taolinesearchsetobjectiveroutine 15 #define taolinesearchsetgradientroutine_ taolinesearchsetgradientroutine 16 #define taolinesearchsetobjectiveandgradientroutine_ taolinesearchsetobjectiveandgradientroutine 17 #define taolinesearchsetobjectiveandgtsroutine_ taolinesearchsetobjectiveandgtsroutine 18 #define taolinesearchview_ taolinesearchview 19 #define taolinesearchsettype_ taolinesearchsettype 20 #define taolinesearchviewfromoptions_ taolinesearchviewfromoptions 21 #endif 22 23 static int OBJ=0; 24 static int GRAD=1; 25 static int OBJGRAD=2; 26 static int OBJGTS=3; 27 static size_t NFUNCS=4; 28 29 static PetscErrorCode ourtaolinesearchobjectiveroutine(TaoLineSearch ls, Vec x, PetscReal *f, void *ctx) 30 { 31 PetscCallFortranVoidFunction((*(void (*)(TaoLineSearch*,Vec*,PetscReal*,void*,PetscErrorCode*)) (((PetscObject)ls)->fortran_func_pointers[OBJ]))(&ls,&x,f,ctx,&ierr)); 32 return 0; 33 } 34 35 static PetscErrorCode ourtaolinesearchgradientroutine(TaoLineSearch ls, Vec x, Vec g, void *ctx) 36 { 37 PetscCallFortranVoidFunction((*(void (*)(TaoLineSearch*,Vec*,Vec*,void*,PetscErrorCode*)) (((PetscObject)ls)->fortran_func_pointers[GRAD]))(&ls,&x,&g,ctx,&ierr)); 38 return 0; 39 } 40 41 static PetscErrorCode ourtaolinesearchobjectiveandgradientroutine(TaoLineSearch ls, Vec x, PetscReal *f, Vec g, void* ctx) 42 { 43 PetscCallFortranVoidFunction((*(void (*)(TaoLineSearch*,Vec*,PetscReal*,Vec*,void*,PetscErrorCode*))(((PetscObject)ls)->fortran_func_pointers[OBJGRAD]))(&ls,&x,f,&g,ctx,&ierr)); 44 return 0; 45 } 46 47 static PetscErrorCode ourtaolinesearchobjectiveandgtsroutine(TaoLineSearch ls, Vec x, Vec s, PetscReal *f, PetscReal *gts, void* ctx) 48 { 49 PetscCallFortranVoidFunction((*(void (*)(TaoLineSearch*,Vec*,Vec*,PetscReal*,PetscReal*,void*,PetscErrorCode*)) (((PetscObject)ls)->fortran_func_pointers[OBJGTS]))(&ls,&x,&s,f,gts,ctx,&ierr)); 50 return 0; 51 } 52 53 PETSC_EXTERN void taolinesearchsetobjectiveroutine_(TaoLineSearch *ls, void (*func)(TaoLineSearch*, Vec *, PetscReal *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 54 { 55 PetscObjectAllocateFortranPointers(*ls,NFUNCS); 56 if (!func) { 57 *ierr = TaoLineSearchSetObjectiveRoutine(*ls,0,ctx); 58 } else { 59 ((PetscObject)*ls)->fortran_func_pointers[OBJ] = (PetscVoidFunction)func; 60 *ierr = TaoLineSearchSetObjectiveRoutine(*ls, ourtaolinesearchobjectiveroutine,ctx); 61 } 62 } 63 64 PETSC_EXTERN void taolinesearchsetgradientroutine_(TaoLineSearch *ls, void (*func)(TaoLineSearch*, Vec *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 65 { 66 PetscObjectAllocateFortranPointers(*ls,NFUNCS); 67 if (!func) { 68 *ierr = TaoLineSearchSetGradientRoutine(*ls,0,ctx); 69 } else { 70 ((PetscObject)*ls)->fortran_func_pointers[GRAD] = (PetscVoidFunction)func; 71 *ierr = TaoLineSearchSetGradientRoutine(*ls, ourtaolinesearchgradientroutine,ctx); 72 } 73 } 74 75 PETSC_EXTERN void taolinesearchsetobjectiveandgradientroutine_(TaoLineSearch *ls, void (*func)(TaoLineSearch*, Vec *, PetscReal *, Vec *, void *, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 76 { 77 PetscObjectAllocateFortranPointers(*ls,NFUNCS); 78 if (!func) { 79 *ierr = TaoLineSearchSetObjectiveAndGradientRoutine(*ls,0,ctx); 80 } else { 81 ((PetscObject)*ls)->fortran_func_pointers[OBJGRAD] = (PetscVoidFunction)func; 82 *ierr = TaoLineSearchSetObjectiveAndGradientRoutine(*ls, ourtaolinesearchobjectiveandgradientroutine,ctx); 83 } 84 } 85 86 PETSC_EXTERN void taolinesearchsetobjectiveandgtsroutine_(TaoLineSearch *ls, void (*func)(TaoLineSearch*, Vec *, Vec *, PetscReal*, PetscReal*,void*, PetscErrorCode *), void *ctx, PetscErrorCode *ierr) 87 { 88 PetscObjectAllocateFortranPointers(*ls,NFUNCS); 89 if (!func) { 90 *ierr = TaoLineSearchSetObjectiveAndGTSRoutine(*ls,0,ctx); 91 } else { 92 ((PetscObject)*ls)->fortran_func_pointers[OBJGTS] = (PetscVoidFunction)func; 93 *ierr = TaoLineSearchSetObjectiveAndGTSRoutine(*ls, ourtaolinesearchobjectiveandgtsroutine,ctx); 94 } 95 } 96 97 PETSC_EXTERN void taolinesearchsettype_(TaoLineSearch *ls, char* type_name, PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len) 98 99 { 100 char *t; 101 102 FIXCHAR(type_name,len,t); 103 *ierr = TaoLineSearchSetType(*ls,t);if (*ierr) return; 104 FREECHAR(type_name,t); 105 106 } 107 108 PETSC_EXTERN void taolinesearchview_(TaoLineSearch *ls, PetscViewer *viewer, PetscErrorCode *ierr) 109 { 110 PetscViewer v; 111 PetscPatchDefaultViewers_Fortran(viewer,v); 112 *ierr = TaoLineSearchView(*ls,v); 113 } 114 115 PETSC_EXTERN void taolinesearchgetoptionsprefix_(TaoLineSearch *ls, char* prefix, PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len) 116 { 117 const char *name; 118 *ierr = TaoLineSearchGetOptionsPrefix(*ls,&name); 119 *ierr = PetscStrncpy(prefix,name,len); if (*ierr) return; 120 FIXRETURNCHAR(PETSC_TRUE,prefix,len); 121 122 } 123 124 PETSC_EXTERN void taolinesearchappendoptionsprefix_(TaoLineSearch *ls, char* prefix,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len) 125 { 126 char *name; 127 FIXCHAR(prefix,len,name); 128 *ierr = TaoLineSearchAppendOptionsPrefix(*ls,name);if (*ierr) return; 129 FREECHAR(prefix,name); 130 } 131 132 PETSC_EXTERN void taolinesearchsetoptionsprefix_(TaoLineSearch *ls, char* prefix,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len) 133 { 134 char *t; 135 FIXCHAR(prefix,len,t); 136 *ierr = TaoLineSearchSetOptionsPrefix(*ls,t);if (*ierr) return; 137 FREECHAR(prefix,t); 138 } 139 140 PETSC_EXTERN void taolinesearchgettype_(TaoLineSearch *ls, char* name, PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len) 141 { 142 const char *tname; 143 *ierr = TaoLineSearchGetType(*ls,&tname); 144 *ierr = PetscStrncpy(name,tname,len); if (*ierr) return; 145 FIXRETURNCHAR(PETSC_TRUE,name,len); 146 147 } 148 PETSC_EXTERN void taolinesearchviewfromoptions_(TaoLineSearch *ao,PetscObject obj,char* type,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len) 149 { 150 char *t; 151 152 FIXCHAR(type,len,t); 153 CHKFORTRANNULLOBJECT(obj); 154 *ierr = TaoLineSearchViewFromOptions(*ao,obj,t);if (*ierr) return; 155 FREECHAR(type,t); 156 } 157