xref: /petsc/src/dm/impls/da/ftn-custom/zdaf.c (revision 8865f1ea4ea560cd84ab8db62e98b7095cdff96f)
12219e2e3SSatish Balay #include <petsc-private/fortranimpl.h>
22219e2e3SSatish Balay #include <petsc-private/daimpl.h>
32219e2e3SSatish Balay 
42219e2e3SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
52219e2e3SSatish Balay #define dmdagetownershipranges_        DMDAGETOWNERSHIPRANGES
62219e2e3SSatish Balay #define dmdagetneighbors_              DMDAGETNEIGHBORS
72219e2e3SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
82219e2e3SSatish Balay #define dmdagetownershipranges_        dmdagetownershipranges
92219e2e3SSatish Balay #define dmdagetneighbors_              dmdagetneighbors
102219e2e3SSatish Balay #endif
112219e2e3SSatish Balay 
122219e2e3SSatish Balay EXTERN_C_BEGIN
132219e2e3SSatish Balay void PETSC_STDCALL dmdagetneighbors_(DM *da,PetscMPIInt *ranks,PetscErrorCode *ierr)
142219e2e3SSatish Balay {
152219e2e3SSatish Balay   const PetscMPIInt *r;
162219e2e3SSatish Balay   PetscInt          n;
172219e2e3SSatish Balay   DM_DA             *dd = (DM_DA*)(*da)->data;
182219e2e3SSatish Balay 
192219e2e3SSatish Balay   *ierr = DMDAGetNeighbors(*da,&r);if (*ierr) return;
20*8865f1eaSKarl Rupp   if (dd->dim == 2) n = 9;
21*8865f1eaSKarl Rupp   else n = 27;
222219e2e3SSatish Balay   *ierr = PetscMemcpy(ranks,r,n*sizeof(PetscMPIInt));
232219e2e3SSatish Balay }
242219e2e3SSatish Balay 
252219e2e3SSatish Balay void PETSC_STDCALL dmdagetownershipranges_(DM *da,PetscInt lx[],PetscInt ly[],PetscInt lz[],PetscErrorCode *ierr)
262219e2e3SSatish Balay {
272219e2e3SSatish Balay   const PetscInt *gx,*gy,*gz;
282219e2e3SSatish Balay   PetscInt       M,N,P,i;
292219e2e3SSatish Balay 
302219e2e3SSatish Balay   CHKFORTRANNULLINTEGER(lx);
312219e2e3SSatish Balay   CHKFORTRANNULLINTEGER(ly);
322219e2e3SSatish Balay   CHKFORTRANNULLINTEGER(lz);
332219e2e3SSatish Balay   *ierr = DMDAGetInfo(*da,0,0,0,0,&M,&N,&P,0,0,0,0,0,0);if (*ierr) return;
342219e2e3SSatish Balay   *ierr = DMDAGetOwnershipRanges(*da,&gx,&gy,&gz);if (*ierr) return;
35*8865f1eaSKarl Rupp   if (lx) {
36*8865f1eaSKarl Rupp     for (i=0; i<M; i++) lx[i] = gx[i];
37*8865f1eaSKarl Rupp   }
38*8865f1eaSKarl Rupp   if (ly) {
39*8865f1eaSKarl Rupp     for (i=0; i<N; i++) ly[i] = gy[i];
40*8865f1eaSKarl Rupp   }
41*8865f1eaSKarl Rupp   if (lz) {
42*8865f1eaSKarl Rupp     for (i=0; i<P; i++) lz[i] = gz[i];
43*8865f1eaSKarl Rupp   }
442219e2e3SSatish Balay }
452219e2e3SSatish Balay 
462219e2e3SSatish Balay EXTERN_C_END
47