Actual source code: finit.c

  1: #include <petsc/private/petscimpl.h>

  3: #if defined(PETSC_USE_FORTRAN_BINDINGS)
  4:   #if defined(PETSC_HAVE_FORTRAN_CAPS)
  5:     #define petscinitializefortran_     PETSCINITIALIZEFORTRAN
  6:     #define petscsetmoduleblock_        PETSCSETMODULEBLOCK
  7:     #define petscsetmoduleblockmpi_     PETSCSETMODULEBLOCKMPI
  8:     #define petscsetmoduleblocknumeric_ PETSCSETMODULEBLOCKNUMERIC
  9:     #define petscsetcomm_               PETSCSETCOMM
 10:   #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 11:     #define petscinitializefortran_     petscinitializefortran
 12:     #define petscsetmoduleblock_        petscsetmoduleblock
 13:     #define petscsetmoduleblockmpi_     petscsetmoduleblockmpi
 14:     #define petscsetmoduleblocknumeric_ petscsetmoduleblocknumeric
 15:     #define petscsetcomm_               petscsetcomm
 16:   #endif

 18: PETSC_EXTERN void petscsetmoduleblock_(void);
 19: PETSC_EXTERN void petscsetmoduleblockmpi_(MPI_Fint *, MPI_Fint *, MPI_Fint *, MPI_Fint *);
 20: PETSC_EXTERN void petscsetmoduleblocknumeric_(PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscReal *);
 21: PETSC_EXTERN void petscsetcomm_(MPI_Fint *, MPI_Fint *);
 22: #endif

 24: /*@C
 25:   PetscInitializeFortran - Routine that should be called soon AFTER
 26:   the call to `PetscInitialize()` if one is using a C main program
 27:   that calls Fortran routines that in turn call PETSc routines.

 29:   Collective on `PETSC_COMM_WORLD`

 31:   Level: beginner

 33:   Notes:
 34:   `PetscInitializeFortran()` initializes some of the default viewers,
 35:   communicators, etc. for use in the Fortran if a user's main program is
 36:   written in C.  `PetscInitializeFortran()` is NOT needed if a user's main
 37:   program is written in Fortran; in this case, just calling
 38:   `PetscInitialize()` in the main (Fortran) program is sufficient.

 40:   This function exists and can be called even if PETSc has been configured
 41:   with --with-fortran-bindings=0 or --with-fc=0. It just does nothing
 42:   in that case.

 44: .seealso: `PetscInitialize()`
 45: @*/
 46: PetscErrorCode PetscInitializeFortran(void)
 47: {
 48: #if defined(PETSC_USE_FORTRAN_BINDINGS)
 49:   MPI_Fint c1 = 0, c2 = 0;

 51:   if (PETSC_COMM_WORLD) c1 = MPI_Comm_c2f(PETSC_COMM_WORLD);
 52:   c2 = MPI_Comm_c2f(PETSC_COMM_SELF);
 53:   petscsetmoduleblock_();
 54:   petscsetcomm_(&c1, &c2);

 56:   {
 57:     MPI_Fint freal, fscalar, fsum, fint;
 58:     freal   = MPI_Type_c2f(MPIU_REAL);
 59:     fscalar = MPI_Type_c2f(MPIU_SCALAR);
 60:     fsum    = MPI_Op_c2f(MPIU_SUM);
 61:     fint    = MPI_Type_c2f(MPIU_INT);
 62:     petscsetmoduleblockmpi_(&freal, &fscalar, &fsum, &fint);
 63:   }

 65:   {
 66:     PetscReal pi      = PETSC_PI;
 67:     PetscReal maxreal = PETSC_MAX_REAL;
 68:     PetscReal minreal = PETSC_MIN_REAL;
 69:     PetscReal eps     = PETSC_MACHINE_EPSILON;
 70:     PetscReal seps    = PETSC_SQRT_MACHINE_EPSILON;
 71:     PetscReal small   = PETSC_SMALL;
 72:     PetscReal pinf    = PETSC_INFINITY;
 73:     PetscReal pninf   = PETSC_NINFINITY;
 74:     petscsetmoduleblocknumeric_(&pi, &maxreal, &minreal, &eps, &seps, &small, &pinf, &pninf);
 75:   }
 76: #endif
 77:   return PETSC_SUCCESS;
 78: }