xref: /petsc/src/vec/is/sf/interface/ftn-custom/zsf.c (revision 030f984af8d8bb4c203755d35bded3c05b3d83ce)
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 #define petscsfdestroy_       PETSCSFDESTROY
12 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
13 #define petscsfgetgraph_      petscsfgetgraph
14 #define petscsfview_          petscsfview
15 #define petscsfbcastbegin_    petscsfbcastbegin
16 #define petscsfbcastend_      petscsfbcastend
17 #define f90arraysfnodecreate_ f90arraysfnodecreate
18 #define petscsfviewfromoptions_ petscsfviewfromoptions
19 #define petscsfdestroy_       petscsfdestroy
20 #endif
21 
22 PETSC_EXTERN void f90arraysfnodecreate_(const PetscInt *,PetscInt *,void * PETSC_F90_2PTR_PROTO_NOVAR);
23 
24 PETSC_EXTERN void petscsfview_(PetscSF *sf, PetscViewer *vin, PetscErrorCode *ierr)
25 {
26   PetscViewer v;
27 
28   PetscPatchDefaultViewers_Fortran(vin, v);
29   *ierr = PetscSFView(*sf, v);
30 }
31 
32 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))
33 {
34   const PetscInt    *ilocal;
35   const PetscSFNode *iremote;
36 
37   *ierr = PetscSFGetGraph(*sf,nroots,nleaves,&ilocal,&iremote);if (*ierr) return;
38   *ierr = F90Array1dCreate((void*)ilocal,MPIU_INT,1,*nleaves, ailocal PETSC_F90_2PTR_PARAM(pilocal));
39   /* this creates a memory leak */
40   f90arraysfnodecreate_((PetscInt*)iremote,nleaves, airemote PETSC_F90_2PTR_PARAM(piremote));
41 }
42 
43 #if defined(PETSC_HAVE_F90_ASSUMED_TYPE_NOT_PTR)
44 PETSC_EXTERN void petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, MPI_Fint *op, PetscErrorCode *ierr)
45 {
46   MPI_Datatype dtype;
47   MPI_Op       cop = MPI_Op_f2c(*op);
48 
49   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
50   *ierr = PetscSFBcastBegin(*sf, dtype, rptr, lptr, cop);
51 }
52 
53 PETSC_EXTERN void petscsfbcastend_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, MPI_Fint *op, PetscErrorCode *ierr)
54 {
55   MPI_Datatype dtype;
56   MPI_Op       cop = MPI_Op_f2c(*op);
57 
58   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
59   *ierr = PetscSFBcastEnd(*sf, dtype, rptr, lptr, cop);
60 }
61 
62 #else
63 
64 PETSC_EXTERN void petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit,F90Array1d *rptr, F90Array1d *lptr, MPI_Fint *op, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(rptrd) PETSC_F90_2PTR_PROTO(lptrd))
65 {
66   MPI_Datatype dtype;
67   const void   *rootdata;
68   void         *leafdata;
69   MPI_Op       cop = MPI_Op_f2c(*op);
70 
71   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
72   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
73   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
74   *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata, cop);
75 }
76 
77 PETSC_EXTERN void petscsfbcastend_(PetscSF *sf, MPI_Fint *unit,F90Array1d *rptr, F90Array1d *lptr, MPI_Fint *op, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(rptrd) PETSC_F90_2PTR_PROTO(lptrd))
78 {
79   MPI_Datatype dtype;
80   const void   *rootdata;
81   void         *leafdata;
82   MPI_Op       cop = MPI_Op_f2c(*op);
83 
84   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
85   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
86   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
87   *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata, cop);
88 }
89 PETSC_EXTERN void petscsfviewfromoptions_(PetscSF *ao,PetscObject obj,char* type,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
90 {
91   char *t;
92 
93   FIXCHAR(type,len,t);
94   CHKFORTRANNULLOBJECT(obj);
95   *ierr = PetscSFViewFromOptions(*ao,obj,t);if (*ierr) return;
96   FREECHAR(type,t);
97 }
98 
99 PETSC_EXTERN void petscsfdestroy_(PetscSF *x,int *ierr)
100 {
101   PETSC_FORTRAN_OBJECT_F_DESTROYED_TO_C_NULL(x);
102   *ierr = PetscSFDestroy(x); if (*ierr) return;
103   PETSC_FORTRAN_OBJECT_C_NULL_TO_F_DESTROYED(x);
104 }
105 
106 #endif
107