xref: /petsc/src/vec/is/sf/interface/ftn-custom/zsf.c (revision 2fb1d9daa8d289b4e375fe89ef72d82069342f85)
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 
33 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))
34 {
35   const PetscInt    *ilocal;
36   const PetscSFNode *iremote;
37 
38   *ierr = PetscSFGetGraph(*sf,nroots,nleaves,&ilocal,&iremote);if (*ierr) return;
39   *ierr = F90Array1dCreate((void*)ilocal,MPIU_INT,1,*nleaves, ailocal PETSC_F90_2PTR_PARAM(pilocal));
40   /* this creates a memory leak */
41   f90arraysfnodecreate_((PetscInt*)iremote,nleaves, airemote PETSC_F90_2PTR_PARAM(piremote));
42 }
43 
44 #if defined(PETSC_HAVE_F90_ASSUMED_TYPE_NOT_PTR)
45 PETSC_EXTERN void petscsfbcastbegin_(PetscSF *sf, MPI_Fint *unit, const void *rptr, void *lptr, MPI_Fint *op, PetscErrorCode *ierr)
46 {
47   MPI_Datatype dtype;
48   MPI_Op       cop = MPI_Op_f2c(*op);
49 
50   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
51   *ierr = PetscSFBcastBegin(*sf, dtype, rptr, lptr, cop);
52 }
53 
54 
55 PETSC_EXTERN void petscsfbcastend_(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 = PetscSFBcastEnd(*sf, dtype, rptr, lptr, cop);
62 }
63 
64 #else
65 
66 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))
67 {
68   MPI_Datatype dtype;
69   const void   *rootdata;
70   void         *leafdata;
71   MPI_Op       cop = MPI_Op_f2c(*op);
72 
73 
74   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
75   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
76   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
77   *ierr = PetscSFBcastBegin(*sf, dtype, rootdata, leafdata, cop);
78 }
79 
80 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))
81 {
82   MPI_Datatype dtype;
83   const void   *rootdata;
84   void         *leafdata;
85   MPI_Op       cop = MPI_Op_f2c(*op);
86 
87   *ierr = PetscMPIFortranDatatypeToC(*unit,&dtype);if (*ierr) return;
88   *ierr = F90Array1dAccess(rptr, dtype, (void**) &rootdata PETSC_F90_2PTR_PARAM(rptrd));if (*ierr) return;
89   *ierr = F90Array1dAccess(lptr, dtype, (void**) &leafdata PETSC_F90_2PTR_PARAM(lptrd));if (*ierr) return;
90   *ierr = PetscSFBcastEnd(*sf, dtype, rootdata, leafdata, cop);
91 }
92 PETSC_EXTERN void petscsfviewfromoptions_(PetscSF *ao,PetscObject obj,char* type,PetscErrorCode *ierr,PETSC_FORTRAN_CHARLEN_T len)
93 {
94   char *t;
95 
96   FIXCHAR(type,len,t);
97   CHKFORTRANNULLOBJECT(obj);
98   *ierr = PetscSFViewFromOptions(*ao,obj,t);if (*ierr) return;
99   FREECHAR(type,t);
100 }
101 
102 PETSC_EXTERN void petscsfdestroy_(PetscSF *x,int *ierr)
103 {
104   PETSC_FORTRAN_OBJECT_F_DESTROYED_TO_C_NULL(x);
105   *ierr = PetscSFDestroy(x); if (*ierr) return;
106   PETSC_FORTRAN_OBJECT_C_NULL_TO_F_DESTROYED(x);
107 }
108 
109 #endif
110