Actual source code: zshellpcf.c

  1: #include <petsc/private/ftnimpl.h>
  2: #include <petscpc.h>
  3: #include <petscksp.h>

  5: #if defined(PETSC_HAVE_FORTRAN_CAPS)
  6:   #define pcshellsetapply_               PCSHELLSETAPPLY
  7:   #define pcshellsetapplysymmetricleft_  PCSHELLSETAPPLYSYMMETRICLEFT
  8:   #define pcshellsetapplysymmetricright_ PCSHELLSETAPPLYSYMMETRICRIGHT
  9:   #define pcshellsetapplyba_             PCSHELLSETAPPLYBA
 10:   #define pcshellsetapplyrichardson_     PCSHELLSETAPPLYRICHARDSON
 11:   #define pcshellsetapplytranspose_      PCSHELLSETAPPLYTRANSPOSE
 12:   #define pcshellsetsetup_               PCSHELLSETSETUP
 13:   #define pcshellsetdestroy_             PCSHELLSETDESTROY
 14:   #define pcshellsetpresolve_            PCSHELLSETPRESOLVE
 15:   #define pcshellsetpostsolve_           PCSHELLSETPOSTSOLVE
 16:   #define pcshellsetview_                PCSHELLSETVIEW
 17: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 18:   #define pcshellsetapply_               pcshellsetapply
 19:   #define pcshellsetapplysymmetricleft_  pcshellsetapplysymmetricleft
 20:   #define pcshellsetapplysymmetricright_ pcshellsetapplysymmetricright
 21:   #define pcshellsetapplyba_             pcshellsetapplyba
 22:   #define pcshellsetapplyrichardson_     pcshellsetapplyrichardson
 23:   #define pcshellsetapplytranspose_      pcshellsetapplytranspose
 24:   #define pcshellsetsetup_               pcshellsetsetup
 25:   #define pcshellsetdestroy_             pcshellsetdestroy
 26:   #define pcshellsetpresolve_            pcshellsetpresolve
 27:   #define pcshellsetpostsolve_           pcshellsetpostsolve
 28:   #define pcshellsetview_                pcshellsetview
 29: #endif

 31: /* These are not extern C because they are passed into non-extern C user level functions */
 32: static PetscErrorCode ourshellapply(PC pc, Vec x, Vec y)
 33: {
 34:   PetscCallFortranVoidFunction((*(void (*)(PC *, Vec *, Vec *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[0]))(&pc, &x, &y, &ierr));
 35:   return PETSC_SUCCESS;
 36: }

 38: static PetscErrorCode ourshellapplysymmetricleft(PC pc, Vec x, Vec y)
 39: {
 40:   PetscCallFortranVoidFunction((*(void (*)(PC *, Vec *, Vec *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[9]))(&pc, &x, &y, &ierr));
 41:   return PETSC_SUCCESS;
 42: }

 44: static PetscErrorCode ourshellapplysymmetricright(PC pc, Vec x, Vec y)
 45: {
 46:   PetscCallFortranVoidFunction((*(void (*)(PC *, Vec *, Vec *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[10]))(&pc, &x, &y, &ierr));
 47:   return PETSC_SUCCESS;
 48: }

 50: static PetscErrorCode ourshellapplyctx(PC pc, Vec x, Vec y)
 51: {
 52:   void *ctx;
 53:   PetscCall(PCShellGetContext(pc, &ctx));
 54:   PetscCallFortranVoidFunction((*(void (*)(PC *, void *, Vec *, Vec *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[0]))(&pc, ctx, &x, &y, &ierr));
 55:   return PETSC_SUCCESS;
 56: }

 58: static PetscErrorCode ourshellapplyba(PC pc, PCSide side, Vec x, Vec y, Vec work)
 59: {
 60:   PetscCallFortranVoidFunction((*(void (*)(PC *, PCSide *, Vec *, Vec *, Vec *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[1]))(&pc, &side, &x, &y, &work, &ierr));
 61:   return PETSC_SUCCESS;
 62: }

 64: static PetscErrorCode ourapplyrichardson(PC pc, Vec x, Vec y, Vec w, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt m, PetscBool guesszero, PetscInt *outits, PCRichardsonConvergedReason *reason)
 65: {
 66:   PetscCallFortranVoidFunction((*(void (*)(PC *, Vec *, Vec *, Vec *, PetscReal *, PetscReal *, PetscReal *, PetscInt *, PetscBool *, PetscInt *, PCRichardsonConvergedReason *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[2]))(&pc, &x, &y, &w, &rtol, &abstol, &dtol, &m, &guesszero, outits, reason, &ierr));
 67:   return PETSC_SUCCESS;
 68: }

 70: static PetscErrorCode ourshellapplytranspose(PC pc, Vec x, Vec y)
 71: {
 72:   PetscCallFortranVoidFunction((*(void (*)(void *, Vec *, Vec *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[3]))(&pc, &x, &y, &ierr));
 73:   return PETSC_SUCCESS;
 74: }

 76: static PetscErrorCode ourshellsetup(PC pc)
 77: {
 78:   PetscCallFortranVoidFunction((*(void (*)(PC *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[4]))(&pc, &ierr));
 79:   return PETSC_SUCCESS;
 80: }

 82: static PetscErrorCode ourshellsetupctx(PC pc)
 83: {
 84:   void *ctx;
 85:   PetscCall(PCShellGetContext(pc, &ctx));
 86:   PetscCallFortranVoidFunction((*(void (*)(PC *, void *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[4]))(&pc, ctx, &ierr));
 87:   return PETSC_SUCCESS;
 88: }

 90: static PetscErrorCode ourshelldestroy(PC pc)
 91: {
 92:   PetscCallFortranVoidFunction((*(void (*)(void *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[5]))(&pc, &ierr));
 93:   return PETSC_SUCCESS;
 94: }

 96: static PetscErrorCode ourshellpresolve(PC pc, KSP ksp, Vec x, Vec y)
 97: {
 98:   PetscCallFortranVoidFunction((*(void (*)(PC *, KSP *, Vec *, Vec *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[6]))(&pc, &ksp, &x, &y, &ierr));
 99:   return PETSC_SUCCESS;
100: }

102: static PetscErrorCode ourshellpostsolve(PC pc, KSP ksp, Vec x, Vec y)
103: {
104:   PetscCallFortranVoidFunction((*(void (*)(PC *, KSP *, Vec *, Vec *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[7]))(&pc, &ksp, &x, &y, &ierr));
105:   return PETSC_SUCCESS;
106: }

108: static PetscErrorCode ourshellview(PC pc, PetscViewer view)
109: {
110:   PetscCallFortranVoidFunction((*(void (*)(PC *, PetscViewer *, PetscErrorCode *))(((PetscObject)pc)->fortran_func_pointers[8]))(&pc, &view, &ierr));
111:   return PETSC_SUCCESS;
112: }

114: PETSC_EXTERN void pcshellsetapply_(PC *pc, void (*apply)(void *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
115: {
116:   PetscObjectAllocateFortranPointers(*pc, 11);
117:   ((PetscObject)*pc)->fortran_func_pointers[0] = (PetscVoidFn *)apply;

119:   *ierr = PCShellSetApply(*pc, ourshellapply);
120: }

122: PETSC_EXTERN void pcshellsetapplysymmetricleft_(PC *pc, void (*apply)(void *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
123: {
124:   PetscObjectAllocateFortranPointers(*pc, 11);
125:   ((PetscObject)*pc)->fortran_func_pointers[9] = (PetscVoidFn *)apply;

127:   *ierr = PCShellSetApplySymmetricLeft(*pc, ourshellapplysymmetricleft);
128: }

130: PETSC_EXTERN void pcshellsetapplysymmetricright_(PC *pc, void (*apply)(void *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
131: {
132:   PetscObjectAllocateFortranPointers(*pc, 11);
133:   ((PetscObject)*pc)->fortran_func_pointers[10] = (PetscVoidFn *)apply;

135:   *ierr = PCShellSetApplySymmetricRight(*pc, ourshellapplysymmetricright);
136: }

138: PETSC_EXTERN void pcshellsetapplyctx_(PC *pc, void (*apply)(void *, void *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
139: {
140:   PetscObjectAllocateFortranPointers(*pc, 11);
141:   ((PetscObject)*pc)->fortran_func_pointers[0] = (PetscVoidFn *)apply;

143:   *ierr = PCShellSetApply(*pc, ourshellapplyctx);
144: }

146: PETSC_EXTERN void pcshellsetapplyba_(PC *pc, void (*apply)(void *, PCSide *, Vec *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
147: {
148:   PetscObjectAllocateFortranPointers(*pc, 11);
149:   ((PetscObject)*pc)->fortran_func_pointers[1] = (PetscVoidFn *)apply;

151:   *ierr = PCShellSetApplyBA(*pc, ourshellapplyba);
152: }

154: PETSC_EXTERN void pcshellsetapplyrichardson_(PC *pc, void (*apply)(void *, Vec *, Vec *, Vec *, PetscReal *, PetscReal *, PetscReal *, PetscInt *, PetscBool *, PetscInt *, PCRichardsonConvergedReason *, PetscErrorCode *), PetscErrorCode *ierr)
155: {
156:   PetscObjectAllocateFortranPointers(*pc, 11);
157:   ((PetscObject)*pc)->fortran_func_pointers[2] = (PetscVoidFn *)apply;
158:   *ierr                                        = PCShellSetApplyRichardson(*pc, ourapplyrichardson);
159: }

161: PETSC_EXTERN void pcshellsetapplytranspose_(PC *pc, void (*applytranspose)(void *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
162: {
163:   PetscObjectAllocateFortranPointers(*pc, 11);
164:   ((PetscObject)*pc)->fortran_func_pointers[3] = (PetscVoidFn *)applytranspose;

166:   *ierr = PCShellSetApplyTranspose(*pc, ourshellapplytranspose);
167: }

169: PETSC_EXTERN void pcshellsetsetupctx_(PC *pc, void (*setup)(void *, void *, PetscErrorCode *), PetscErrorCode *ierr)
170: {
171:   PetscObjectAllocateFortranPointers(*pc, 11);
172:   ((PetscObject)*pc)->fortran_func_pointers[4] = (PetscVoidFn *)setup;

174:   *ierr = PCShellSetSetUp(*pc, ourshellsetupctx);
175: }

177: PETSC_EXTERN void pcshellsetsetup_(PC *pc, void (*setup)(void *, PetscErrorCode *), PetscErrorCode *ierr)
178: {
179:   PetscObjectAllocateFortranPointers(*pc, 11);
180:   ((PetscObject)*pc)->fortran_func_pointers[4] = (PetscVoidFn *)setup;

182:   *ierr = PCShellSetSetUp(*pc, ourshellsetup);
183: }

185: PETSC_EXTERN void pcshellsetdestroy_(PC *pc, void (*setup)(void *, PetscErrorCode *), PetscErrorCode *ierr)
186: {
187:   PetscObjectAllocateFortranPointers(*pc, 11);
188:   ((PetscObject)*pc)->fortran_func_pointers[5] = (PetscVoidFn *)setup;

190:   *ierr = PCShellSetDestroy(*pc, ourshelldestroy);
191: }

193: PETSC_EXTERN void pcshellsetpresolve_(PC *pc, void (*presolve)(void *, void *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
194: {
195:   PetscObjectAllocateFortranPointers(*pc, 11);
196:   ((PetscObject)*pc)->fortran_func_pointers[6] = (PetscVoidFn *)presolve;

198:   *ierr = PCShellSetPreSolve(*pc, ourshellpresolve);
199: }

201: PETSC_EXTERN void pcshellsetpostsolve_(PC *pc, void (*postsolve)(void *, void *, Vec *, Vec *, PetscErrorCode *), PetscErrorCode *ierr)
202: {
203:   PetscObjectAllocateFortranPointers(*pc, 11);
204:   ((PetscObject)*pc)->fortran_func_pointers[7] = (PetscVoidFn *)postsolve;

206:   *ierr = PCShellSetPostSolve(*pc, ourshellpostsolve);
207: }

209: PETSC_EXTERN void pcshellsetview_(PC *pc, void (*view)(void *, PetscViewer *, PetscErrorCode *), PetscErrorCode *ierr)
210: {
211:   PetscObjectAllocateFortranPointers(*pc, 11);
212:   ((PetscObject)*pc)->fortran_func_pointers[8] = (PetscVoidFn *)view;

214:   *ierr = PCShellSetView(*pc, ourshellview);
215: }