1 #include <petsc/private/f90impl.h> 2 #include <petsc/private/sfimpl.h> 3 4 #if defined(PETSC_HAVE_FORTRAN_CAPS) 5 #define petscsfview_ PETSCSFVIEW 6 #define petscsfgetgraph_ PETSCSFGETGRAPH 7 #define petscsfbcastbegin_ PETSCSFBCASTBEGIN 8 #define petscsfbcastend_ PETSCSFBCASTEND 9 #define f90arraysfnodecreate_ F90ARRAYSFNODECREATE 10 #define petscsfviewfromoptions_ PETSCSFVIEWFROMOPTIONS 11 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 12 #define petscsfgetgraph_ petscsfgetgraph 13 #define petscsfview_ petscsfview 14 #define petscsfbcastbegin_ petscsfbcastbegin 15 #define petscsfbcastend_ petscsfbcastend 16 #define f90arraysfnodecreate_ f90arraysfnodecreate 17 #define petscsfviewfromoptions_ petscsfviewfromoptions 18 #endif 19 20 PETSC_EXTERN void f90arraysfnodecreate_(const PetscInt *,PetscInt *,void * PETSC_F90_2PTR_PROTO_NOVAR); 21 22 PETSC_EXTERN void petscsfview_(PetscSF *sf, PetscViewer *vin, PetscErrorCode *ierr) 23 { 24 PetscViewer v; 25 26 PetscPatchDefaultViewers_Fortran(vin, v); 27 *ierr = PetscSFView(*sf, v); 28 } 29 30 31 PETSC_EXTERN void petscsfgetgraph_(PetscSF *sf,PetscInt *nroots,PetscInt *nleaves, F90Array1d *ailocal, F90Array1d *airemote, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(pilocal) PETSC_F90_2PTR_PROTO(piremote)) 32 { 33 const PetscInt *ilocal; 34 const PetscSFNode *iremote; 35 36 *ierr = PetscSFGetGraph(*sf,nroots,nleaves,&ilocal,&iremote);if (*ierr) return; 37 *ierr = F90Array1dCreate((void*)ilocal,MPIU_INT,1,*nleaves, ailocal PETSC_F90_2PTR_PARAM(pilocal)); 38 /* this creates a memory leak */ 39 f90arraysfnodecreate_((PetscInt*)iremote,nleaves, airemote PETSC_F90_2PTR_PARAM(piremote)); 40 } 41 42 #if defined(PETSC_HAVE_F90_ASSUMED_TYPE_NOT_PTR) 43 PETSC_EXTERN void petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, PetscErrorCode *ierr) 44 { 45 MPI_Datatype dtype; 46 47 *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return; 48 *ierr = PetscSFBcastBegin(*sf, dtype, rptr, lptr); 49 } 50 51 52 PETSC_EXTERN void petscsfbcastend_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, PetscErrorCode *ierr) 53 { 54 MPI_Datatype dtype; 55 56 *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return; 57 *ierr = PetscSFBcastEnd(*sf, dtype, rptr, lptr); 58 } 59 60 #else 61 62 PETSC_EXTERN void petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit,F90Array1d *rptr, F90Array1d *lptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(rptrd) PETSC_F90_2PTR_PROTO(lptrd)) 63 { 64 MPI_Datatype dtype; 65 const void *rootdata; 66 void *leafdata; 67 68 69 *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return; 70 *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return; 71 *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return; 72 *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata); 73 } 74 75 PETSC_EXTERN void petscsfbcastend_(PetscSF *sf, MPI_Fint *unit,F90Array1d *rptr, F90Array1d *lptr, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(rptrd) PETSC_F90_2PTR_PROTO(lptrd)) 76 { 77 MPI_Datatype dtype; 78 const void *rootdata; 79 void *leafdata; 80 81 *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return; 82 *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return; 83 *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return; 84 *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata); 85 } 86 PETSC_EXTERN void petscsfviewfromoptions_(PetscSF *ao,PetscObject obj,char* type,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len) 87 { 88 char *t; 89 90 FIXCHAR(type,len,t); 91 *ierr = PetscSFViewFromOptions(*ao,obj,t);if (*ierr) return; 92 FREECHAR(type,t); 93 } 94 95 #endif 96