Actual source code: ex7.c

petsc-3.3-p7 2013-05-11
  2: static char help[] = "Demonstrates calling a Fortran computational routine from C.\n\
  3: Also demonstrates passing  PETSc objects, MPI Communicators from C to Fortran\n\
  4: and from Fortran to C\n\n";

  6: #include <petscvec.h>
  7: /*
  8:   Ugly stuff to insure the function names match between Fortran 
  9:   and C. Sorry, but this is out of our PETSc hands to cleanup.
 10: */
 11: #include <petsc-private/fortranimpl.h>
 12: #if defined(PETSC_HAVE_FORTRAN_CAPS)
 13: #define ex7f_ EX7F
 14: #define ex7c_ EX7C
 15: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 16: #define ex7f_ ex7f
 17: #define ex7c_ ex7c
 18: #endif
 19: EXTERN_C_BEGIN
 20: extern void PETSC_STDCALL ex7f_(Vec *,int*);
 21: EXTERN_C_END

 25: int main(int argc,char **args)
 26: {
 27:   PetscErrorCode   ierr;
 28:   PetscInt         m = 10;
 29:   int              fcomm;
 30:   Vec              vec;

 32:   PetscInitialize(&argc,&args,(char *)0,help);

 34:   /* This function should be called to be able to use PETSc routines
 35:      from the FORTRAN subroutines needed by this program */

 37:   PetscInitializeFortran();

 39:   VecCreate(PETSC_COMM_WORLD,&vec);
 40:   VecSetSizes(vec,PETSC_DECIDE,m);
 41:   VecSetFromOptions(vec);

 43:   /* 
 44:      Call Fortran routine - the use of MPI_Comm_c2f() allows
 45:      translation of the MPI_Comm from C so that it can be properly 
 46:      interpreted from Fortran.
 47:   */
 48:   fcomm = MPI_Comm_c2f(PETSC_COMM_WORLD);

 50:   ex7f_(&vec,&fcomm);

 52:   VecView(vec,PETSC_VIEWER_STDOUT_WORLD);
 53:   VecDestroy(&vec);
 54:   PetscFinalize();
 55:   return 0;
 56: }

 58: EXTERN_C_BEGIN
 61: void PETSC_STDCALL ex7c_(Vec *fvec,int *fcomm,PetscErrorCode* ierr)
 62: {
 63:   MPI_Comm comm;
 64:   PetscInt vsize;

 66:   /*
 67:     Translate Fortran integer pointer back to C and
 68:     Fortran Communicator back to C communicator
 69:   */
 70:   comm = MPI_Comm_f2c(*fcomm);
 71: 
 72:   /* Some PETSc/MPI operations on Vec/Communicator objects */
 73:   *VecGetSize(*fvec,&vsize);
 74:   *MPI_Barrier(comm);

 76: }
 77: EXTERN_C_END