1 #include <petsc/private/fortranimpl.h> 2 #include <petsc/private/dmdaimpl.h> 3 4 #if defined(PETSC_HAVE_FORTRAN_CAPS) 5 #define dmdagetownershipranges_ DMDAGETOWNERSHIPRANGES 6 #define dmdagetneighbors_ DMDAGETNEIGHBORS 7 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 8 #define dmdagetownershipranges_ dmdagetownershipranges 9 #define dmdagetneighbors_ dmdagetneighbors 10 #endif 11 12 PETSC_EXTERN void dmdagetneighbors_(DM *da, PetscMPIInt *ranks, PetscErrorCode *ierr) 13 { 14 const PetscMPIInt *r; 15 PetscInt n, dim; 16 17 *ierr = DMDAGetNeighbors(*da, &r); 18 if (*ierr) return; 19 *ierr = DMGetDimension(*da, &dim); 20 if (*ierr) return; 21 if (dim == 2) n = 9; 22 else n = 27; 23 *ierr = PetscArraycpy(ranks, r, n); 24 } 25 26 PETSC_EXTERN void dmdagetownershipranges_(DM *da, PetscInt lx[], PetscInt ly[], PetscInt lz[], PetscErrorCode *ierr) 27 { 28 const PetscInt *gx, *gy, *gz; 29 PetscInt M, N, P, i; 30 31 CHKFORTRANNULLINTEGER(lx); 32 CHKFORTRANNULLINTEGER(ly); 33 CHKFORTRANNULLINTEGER(lz); 34 *ierr = DMDAGetInfo(*da, NULL, NULL, NULL, NULL, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL); 35 if (*ierr) return; 36 *ierr = DMDAGetOwnershipRanges(*da, &gx, &gy, &gz); 37 if (*ierr) return; 38 if (lx) { 39 for (i = 0; i < M; i++) lx[i] = gx[i]; 40 } 41 if (ly) { 42 for (i = 0; i < N; i++) ly[i] = gy[i]; 43 } 44 if (lz) { 45 for (i = 0; i < P; i++) lz[i] = gz[i]; 46 } 47 } 48