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