Actual source code: ex7.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  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

 20: PETSC_EXTERN void PETSC_STDCALL ex7f_(Vec*,int*);

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

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

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

 36:   PetscInitializeFortran();

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

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

 49:   ex7f_(&vec,&fcomm);

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

 59: PETSC_EXTERN void PETSC_STDCALL ex7c_(Vec *fvec,int *fcomm,PetscErrorCode *ierr)
 60: {
 61:   MPI_Comm comm;
 62:   PetscInt vsize;

 64:   /*
 65:     Translate Fortran integer pointer back to C and
 66:     Fortran Communicator back to C communicator
 67:   */
 68:   comm = MPI_Comm_f2c(*fcomm);

 70:   /* Some PETSc/MPI operations on Vec/Communicator objects */
 71:   *VecGetSize(*fvec,&vsize);
 72:   *MPI_Barrier(comm);

 74: }