Actual source code: ztaolinesearchf.c

  1: #include <petsc/private/fortranimpl.h>
  2: #include <petsc/private/taolinesearchimpl.h>

  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)

 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

 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;

 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 PETSC_SUCCESS;
 33: }

 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 PETSC_SUCCESS;
 39: }

 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 PETSC_SUCCESS;
 45: }

 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 PETSC_SUCCESS;
 51: }

 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, NULL, ctx);
 58:   } else {
 59:     ((PetscObject)*ls)->fortran_func_pointers[OBJ] = (PetscVoidFn *)func;
 60:     *ierr                                          = TaoLineSearchSetObjectiveRoutine(*ls, ourtaolinesearchobjectiveroutine, ctx);
 61:   }
 62: }

 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, NULL, ctx);
 69:   } else {
 70:     ((PetscObject)*ls)->fortran_func_pointers[GRAD] = (PetscVoidFn *)func;
 71:     *ierr                                           = TaoLineSearchSetGradientRoutine(*ls, ourtaolinesearchgradientroutine, ctx);
 72:   }
 73: }

 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, NULL, ctx);
 80:   } else {
 81:     ((PetscObject)*ls)->fortran_func_pointers[OBJGRAD] = (PetscVoidFn *)func;
 82:     *ierr                                              = TaoLineSearchSetObjectiveAndGradientRoutine(*ls, ourtaolinesearchobjectiveandgradientroutine, ctx);
 83:   }
 84: }

 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, NULL, ctx);
 91:   } else {
 92:     ((PetscObject)*ls)->fortran_func_pointers[OBJGTS] = (PetscVoidFn *)func;
 93:     *ierr                                             = TaoLineSearchSetObjectiveAndGTSRoutine(*ls, ourtaolinesearchobjectiveandgtsroutine, ctx);
 94:   }
 95: }

 97: PETSC_EXTERN void taolinesearchsettype_(TaoLineSearch *ls, char *type_name, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)

 99: {
100:   char *t;

102:   FIXCHAR(type_name, len, t);
103:   *ierr = TaoLineSearchSetType(*ls, t);
104:   if (*ierr) return;
105:   FREECHAR(type_name, t);
106: }

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: }

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);
120:   if (*ierr) return;
121:   FIXRETURNCHAR(PETSC_TRUE, prefix, len);
122: }

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);
129:   if (*ierr) return;
130:   FREECHAR(prefix, name);
131: }

133: PETSC_EXTERN void taolinesearchsetoptionsprefix_(TaoLineSearch *ls, char *prefix, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)
134: {
135:   char *t;
136:   FIXCHAR(prefix, len, t);
137:   *ierr = TaoLineSearchSetOptionsPrefix(*ls, t);
138:   if (*ierr) return;
139:   FREECHAR(prefix, t);
140: }

142: PETSC_EXTERN void taolinesearchgettype_(TaoLineSearch *ls, char *name, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)
143: {
144:   const char *tname;
145:   *ierr = TaoLineSearchGetType(*ls, &tname);
146:   *ierr = PetscStrncpy(name, tname, len);
147:   if (*ierr) return;
148:   FIXRETURNCHAR(PETSC_TRUE, name, len);
149: }
150: PETSC_EXTERN void taolinesearchviewfromoptions_(TaoLineSearch *ao, PetscObject obj, char *type, PetscErrorCode *ierr, PETSC_FORTRAN_CHARLEN_T len)
151: {
152:   char *t;

154:   FIXCHAR(type, len, t);
155:   CHKFORTRANNULLOBJECT(obj);
156:   *ierr = TaoLineSearchViewFromOptions(*ao, obj, t);
157:   if (*ierr) return;
158:   FREECHAR(type, t);
159: }