Actual source code: zis.c

  1: /*$Id: zis.c,v 1.41 2001/06/21 21:19:50 bsmith Exp $*/

 3:  #include src/fortran/custom/zpetsc.h
 4:  #include petscis.h
  5: #ifdef PETSC_HAVE_FORTRAN_CAPS
  6: #define ispartitioningcount_   ISPARTITIONINGCOUNT
  7: #define isdestroy_             ISDESTROY
  8: #define iscreatestride_        ISCREATESTRIDE
  9: #define iscreategeneral_       ISCREATEGENERAL
 10: #define isgetindices_          ISGETINDICES
 11: #define isrestoreindices_      ISRESTOREINDICES
 12: #define isblockgetindices_     ISBLOCKGETINDICES
 13: #define isblockrestoreindices_ ISBLOCKRESTOREINDICES
 14: #define iscreateblock_         ISCREATEBLOCK
 15: #define isblock_               ISBLOCK
 16: #define isstride_              ISSTRIDE
 17: #define ispermutation_         ISPERMUTATION
 18: #define isidentity_            ISIDENTITY
 19: #define issorted_              ISSORTED
 20: #define isequal_               ISEQUAL
 21: #define isinvertpermutation_   ISINVERTPERMUTATION
 22: #define isview_                ISVIEW
 23: #define iscoloringcreate_      ISCOLORINGCREATE
 24: #define islocaltoglobalmappingcreate_ ISLOCALTOGLOBALMAPPINGCREATE
 25: #define islocaltoglobalmappingblock_ ISLOCALTOGLOBALMAPPINGBLOCK
 26: #define isallgather_                  ISALLGATHER
 27: #define iscoloringdestroy_            ISCOLORINGDESTROY
 28: #define iscoloringview_               ISCOLORINGVIEW
 29: #define ispartitioningtonumbering_    ISPARTITIONINGTONUMBERING
 30: #define islocaltoglobalmappingapply_  ISLOCALTOGLOBALMAPPINGAPPLY
 31: #define islocaltoglobalmappingview_  ISLOCALTOGLOBALMAPPINGVIEW
 32: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 33: #define islocaltoglobalmappingview_   islocaltoglobalmappingview
 34: #define islocaltoglobalmappingapply_  islocaltoglobalmappingapply
 35: #define iscoloringview_        iscoloringview
 36: #define iscoloringdestroy_     iscoloringdestroy
 37: #define isview_                isview
 38: #define isinvertpermutation_   isinvertpermutation
 39: #define isdestroy_             isdestroy
 40: #define iscreatestride_        iscreatestride
 41: #define iscreategeneral_       iscreategeneral
 42: #define isgetindices_          isgetindices
 43: #define isrestoreindices_      isrestoreindices
 44: #define isblockgetindices_     isblockgetindices
 45: #define isblockrestoreindices_ isblockrestoreindices
 46: #define iscreateblock_         iscreateblock
 47: #define isblock_               isblock
 48: #define isstride_              isstride
 49: #define ispermutation_         ispermutation
 50: #define isidentity_            isidentity
 51: #define issorted_              issorted
 52: #define isequal_               isequal
 53: #define iscoloringcreate_      iscoloringcreate
 54: #define islocaltoglobalmappingcreate_ islocaltoglobalmappingcreate
 55: #define islocaltoglobalmappingblock_ islocaltoglobalmappingblock
 56: #define isallgather_                  isallgather
 57: #define ispartitioningcount_          ispartitioningcount
 58: #define ispartitioningtonumbering_    ispartitioningtonumbering
 59: #endif

 61: EXTERN_C_BEGIN

 63: void PETSC_STDCALL islocaltoglobalmappingview_(ISLocalToGlobalMapping *mapping,PetscViewer *viewer,int *ierr)
 64: {
 65:   PetscViewer v = PETSC_NULL;
 66:   CHKFORTRANNULLOBJECT(viewer) else v = *viewer;
 67:   *ISLocalToGlobalMappingView(*mapping,v);
 68: }

 70: /*
 71:    This is the same as the macro ISLocalToGlobalMappingApply() except it does not
 72:   return error codes.
 73: */
 74: void PETSC_STDCALL islocaltoglobalmappingapply_(ISLocalToGlobalMapping *mapping,int *N,int *in,int *out,int *ierr)
 75: {
 76:   int i,*idx = (*mapping)->indices,Nmax = (*mapping)->n;
 77:   for (i=0; i<(*N); i++) {
 78:     if (in[i] < 0) {out[i] = in[i]; continue;}
 79:     if (in[i] >= Nmax) {
 80:       *PetscError(__LINE__,"ISLocalToGlobalMappingApply_Fortran",__FILE__,__SDIR__,1,1,"Index out of range");
 81:       return;
 82:     }
 83:     out[i] = idx[in[i]];
 84:   }
 85: }

 87: void PETSC_STDCALL ispartitioningtonumbering_(IS *is,IS *isout,int *ierr)
 88: {
 89:   *ISPartitioningToNumbering(*is,isout);
 90: }

 92: void PETSC_STDCALL ispartitioningcount_(IS *is,int *count,int *ierr)
 93: {
 94:   *ISPartitioningCount(*is,count);
 95: }

 97: void PETSC_STDCALL iscoloringdestroy_(ISColoring *iscoloring,int *ierr)
 98: {
 99:   *ISColoringDestroy(*iscoloring);
100: }

102: void PETSC_STDCALL iscoloringview_(ISColoring *iscoloring,PetscViewer *viewer,int *ierr)
103: {
104:   PetscViewer v;
105:   PetscPatchDefaultViewers_Fortran(viewer,v);
106:   *ISColoringView(*iscoloring,v);
107: }

109: void PETSC_STDCALL isview_(IS *is,PetscViewer *vin,int *ierr)
110: {
111:   PetscViewer v;
112:   PetscPatchDefaultViewers_Fortran(vin,v);
113:   *ISView(*is,v);
114: }

116: void PETSC_STDCALL isequal_(IS *is1,IS *is2,PetscTruth *flg,int *ierr)
117: {
118:   *ISEqual(*is1,*is2,flg);
119: }

121: void PETSC_STDCALL isidentity_(IS *is,PetscTruth *ident,int *ierr)
122: {
123:   *ISIdentity(*is,ident);
124: }

126: void PETSC_STDCALL issorted_(IS *is,PetscTruth *flg,int *ierr)
127: {
128:   *ISSorted(*is,flg);
129: }

131: void PETSC_STDCALL ispermutation_(IS *is,PetscTruth *perm,int *ierr){
132:   *ISPermutation(*is,perm);
133: }

135: void PETSC_STDCALL isstride_(IS *is,PetscTruth *flag,int *ierr)
136: {
137:   *ISStride(*is,flag);
138: }

140: void PETSC_STDCALL isblockgetindices_(IS *x,int *fa,long *ia,int *ierr)
141: {
142:   int   *lx;

144:   *ISGetIndices(*x,&lx); if (*ierr) return;
145:   *ia      = PetscIntAddressToFortran(fa,lx);
146: }

148: void PETSC_STDCALL isblockrestoreindices_(IS *x,int *fa,long *ia,int *ierr)
149: {
150:   int *lx = PetscIntAddressFromFortran(fa,*ia);

152:   *ISRestoreIndices(*x,&lx);
153: }

155: void PETSC_STDCALL isblock_(IS *is,PetscTruth *flag,int *ierr)
156: {
157:   *ISBlock(*is,flag);
158: }

160: void PETSC_STDCALL isgetindices_(IS *x,int *fa,long *ia,int *ierr)
161: {
162:   int   *lx;

164:   *ISGetIndices(*x,&lx); if (*ierr) return;
165:   *ia      = PetscIntAddressToFortran(fa,lx);
166: }

168: void PETSC_STDCALL isrestoreindices_(IS *x,int *fa,long *ia,int *ierr)
169: {
170:   int *lx = PetscIntAddressFromFortran(fa,*ia);

172:   *ISRestoreIndices(*x,&lx);
173: }

175: void PETSC_STDCALL iscreategeneral_(MPI_Comm *comm,int *n,int *idx,IS *is,int *ierr)
176: {
177:   *ISCreateGeneral((MPI_Comm)PetscToPointerComm(*comm),*n,idx,is);
178: }

180: void PETSC_STDCALL isinvertpermutation_(IS *is,int *nlocal,IS *isout,int *ierr)
181: {
182:   *ISInvertPermutation(*is,*nlocal,isout);
183: }

185: void PETSC_STDCALL iscreateblock_(MPI_Comm *comm,int *bs,int *n,int *idx,IS *is,int *ierr)
186: {
187:   *ISCreateBlock((MPI_Comm)PetscToPointerComm(*comm),*bs,*n,idx,is);
188: }

190: void PETSC_STDCALL iscreatestride_(MPI_Comm *comm,int *n,int *first,int *step,
191:                                IS *is,int *ierr)
192: {
193:   *ISCreateStride((MPI_Comm)PetscToPointerComm(*comm),*n,*first,*step,is);
194: }

196: void PETSC_STDCALL isdestroy_(IS *is,int *ierr)
197: {
198:   *ISDestroy(*is);
199: }

201: void PETSC_STDCALL iscoloringcreate_(MPI_Comm *comm,int *n,int *colors,ISColoring *iscoloring,int *ierr)
202: {
203:   int *color;
204:   /* copies the colors[] array since that is kept by the ISColoring that is created */
205:   *PetscMalloc((*n+1)*sizeof(int),&color);if (*ierr) return;
206:   *PetscMemcpy(color,colors,(*n)*sizeof(int));if (*ierr) return;
207:   *ISColoringCreate((MPI_Comm)PetscToPointerComm(*comm),*n,color,iscoloring);
208: }

210: void PETSC_STDCALL islocaltoglobalmappingcreate_(MPI_Comm *comm,int *n,int *indices,ISLocalToGlobalMapping *mapping,int *ierr)
211: {
212:   *ISLocalToGlobalMappingCreate((MPI_Comm)PetscToPointerComm(*comm),*n,indices,mapping);
213: }

215: void PETSC_STDCALL islocaltoglobalmappingblock_(ISLocalToGlobalMapping *inmap,int bs,ISLocalToGlobalMapping *outmap,int *ierr)
216: {
217:   *ISLocalToGlobalMappingBlock(*inmap,bs,outmap);
218: }

220: void PETSC_STDCALL isallgather_(IS *is,IS *isout,int *ierr)
221: {
222:   *ISAllGather(*is,isout);

224: }

226: EXTERN_C_END