xref: /petsc/src/dm/impls/da/ftn-custom/zdaf.c (revision 1fd49c2520df80b4b23f7bac1a7aae848a1a6ec1)
1af0996ceSBarry Smith #include <petsc/private/fortranimpl.h>
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>
32219e2e3SSatish Balay 
42219e2e3SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
52219e2e3SSatish Balay #define dmdagetownershipranges_        DMDAGETOWNERSHIPRANGES
62219e2e3SSatish Balay #define dmdagetneighbors_              DMDAGETNEIGHBORS
7*1fd49c25SBarry Smith #define dmdagetrefinementfactor_       DMDAGETREFINEMENTFACTOR
82219e2e3SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
92219e2e3SSatish Balay #define dmdagetownershipranges_        dmdagetownershipranges
102219e2e3SSatish Balay #define dmdagetneighbors_              dmdagetneighbors
11*1fd49c25SBarry Smith #define dmdagetrefinementfactor_       dmdagetrefinementfactor
122219e2e3SSatish Balay #endif
132219e2e3SSatish Balay 
148cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmdagetneighbors_(DM *da,PetscMPIInt *ranks,PetscErrorCode *ierr)
152219e2e3SSatish Balay {
162219e2e3SSatish Balay   const PetscMPIInt *r;
17c73cfb54SMatthew G. Knepley   PetscInt          n,dim;
182219e2e3SSatish Balay 
192219e2e3SSatish Balay   *ierr = DMDAGetNeighbors(*da,&r);if (*ierr) return;
20c73cfb54SMatthew G. Knepley   *ierr = DMGetDimension(*da,&dim);if (*ierr) return;
21c73cfb54SMatthew G. Knepley   if (dim == 2) n = 9;
228865f1eaSKarl Rupp   else n = 27;
232219e2e3SSatish Balay   *ierr = PetscMemcpy(ranks,r,n*sizeof(PetscMPIInt));
242219e2e3SSatish Balay }
252219e2e3SSatish Balay 
268cc058d9SJed Brown PETSC_EXTERN void PETSC_STDCALL dmdagetownershipranges_(DM *da,PetscInt lx[],PetscInt ly[],PetscInt lz[],PetscErrorCode *ierr)
272219e2e3SSatish Balay {
282219e2e3SSatish Balay   const PetscInt *gx,*gy,*gz;
292219e2e3SSatish Balay   PetscInt       M,N,P,i;
302219e2e3SSatish Balay 
312219e2e3SSatish Balay   CHKFORTRANNULLINTEGER(lx);
322219e2e3SSatish Balay   CHKFORTRANNULLINTEGER(ly);
332219e2e3SSatish Balay   CHKFORTRANNULLINTEGER(lz);
342219e2e3SSatish Balay   *ierr = DMDAGetInfo(*da,0,0,0,0,&M,&N,&P,0,0,0,0,0,0);if (*ierr) return;
352219e2e3SSatish Balay   *ierr = DMDAGetOwnershipRanges(*da,&gx,&gy,&gz);if (*ierr) return;
368865f1eaSKarl Rupp   if (lx) {
378865f1eaSKarl Rupp     for (i=0; i<M; i++) lx[i] = gx[i];
388865f1eaSKarl Rupp   }
398865f1eaSKarl Rupp   if (ly) {
408865f1eaSKarl Rupp     for (i=0; i<N; i++) ly[i] = gy[i];
418865f1eaSKarl Rupp   }
428865f1eaSKarl Rupp   if (lz) {
438865f1eaSKarl Rupp     for (i=0; i<P; i++) lz[i] = gz[i];
448865f1eaSKarl Rupp   }
452219e2e3SSatish Balay }
462219e2e3SSatish Balay 
47*1fd49c25SBarry Smith PETSC_EXTERN void PETSC_STDCALL dmdagetrefinementfactor_(DM *da,PetscInt *lx,PetscInt *ly,PetscInt *lz,PetscErrorCode *ierr)
48*1fd49c25SBarry Smith {
49*1fd49c25SBarry Smith   CHKFORTRANNULLINTEGER(lx);
50*1fd49c25SBarry Smith   CHKFORTRANNULLINTEGER(ly);
51*1fd49c25SBarry Smith   CHKFORTRANNULLINTEGER(lz);
52*1fd49c25SBarry Smith   *ierr = DMDAGetRefinementFactor(*da,lx,ly,lz);
53*1fd49c25SBarry Smith }
54*1fd49c25SBarry Smith 
55