Actual source code: zda.c

  1: /*$Id: zda.c,v 1.49 2001/08/06 21:19:11 bsmith Exp $*/

 3:  #include src/dm/da/daimpl.h
 4:  #include src/fortran/custom/zpetsc.h
 5:  #include petscmat.h
 6:  #include petscda.h

  8: #ifdef PETSC_HAVE_FORTRAN_CAPS
  9: #define dasetlocalfunction_          DASETLOCALFUNCTION
 10: #define dasetLocaladiforfunction_    DASETLOCALADIFORFUNCTION
 11: #define dasetlocaladiformffunction_  DASETLOCALADIFORMFFUNCTION
 12: #define dasetlocaljacobian_          DASETLOCALJACOBIAN
 13: #define dagetlocalinfo_              DAGETLOCALINFO
 14: #define dagetinterpolation_          DAGETINTERPOLATION
 15: #define dacreate1d_                  DACREATE1D
 16: #define dacreate3d_                  DACREATE3D
 17: #define dacreate2d_                  DACREATE2D
 18: #define dadestroy_                   DADESTROY
 19: #define dacreateglobalvector_        DACREATEGLOBALVECTOR
 20: #define dacreatelocalvector_         DACREATELOCALVECTOR
 21: #define dagetlocalvector_            DAGETLOCALVECTOR
 22: #define darestorelocalvector_        DARESTORELOCALVECTOR
 23: #define dagetscatter_                DAGETSCATTER
 24: #define dagetglobalindices_          DAGETGLOBALINDICES
 25: #define daview_                      DAVIEW
 26: #define dagetinfo_                   DAGETINFO
 27: #define dagetcoloring_               DAGETCOLORING
 28: #define dagetmatrix_                 DAGETMATRIX
 29: #define dagetislocaltoglobalmapping_ DAGETISLOCALTOGLOBALMAPPING
 30: #define dagetislocaltoglobalmappingblck_ DAGETISLOCALTOGLOBALMAPPINGBLCK
 31: #define daload_                      DALOAD
 32: #define dasetfieldname_              DASETFIELDNAME
 33: #define dagetfieldname_              DAGETFIELDNAME
 34: #define darefine_                    DAREFINE
 35: #define dagetao_                     DAGETAO
 36: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 37: #define dagetlocalinfo_              dagetlocalinfo
 38: #define dagetlocalvector_            dagetlocalvector
 39: #define darestorelocalvector_        darestorelocalvector
 40: #define dagetinterpolation_          dagetinterpolation
 41: #define daload_                      daload
 42: #define dacreateglobalvector_        dacreateglobalvector
 43: #define dacreatelocalvector_         dacreatelocalvector
 44: #define daview_                      daview
 45: #define dacreate1d_                  dacreate1d
 46: #define dacreate3d_                  dacreate3d
 47: #define dacreate2d_                  dacreate2d
 48: #define dadestroy_                   dadestroy
 49: #define dagetscatter_                dagetscatter
 50: #define dagetglobalindices_          dagetglobalindices
 51: #define dagetinfo_                   dagetinfo
 52: #define dagetcoloring_               dagetcoloring
 53: #define dagetmatrix_                 dagetmatrix
 54: #define dagetislocaltoglobalmapping_ dagetislocaltoglobalmapping
 55: #define dagetislocaltoglobalmappingblck_ dagetislocaltoglobalmappingblck
 56: #define dasetfieldname_              dasetfieldname
 57: #define dagetfieldname_              dagetfieldname
 58: #define darefine_                    darefine
 59: #define dagetao_                     dagetao
 60: #define dasetlocalfunction_          dasetlocalfunction
 61: #define dasetlocaladiforfunction_       dasetlocaladiforfunction
 62: #define dasetlocaladiformffunction_       dasetlocaladiformffunction
 63: #define dasetlocaljacobian_          dasetlocaljacobian
 64: #endif



 68: EXTERN_C_BEGIN


 71: static void (PETSC_STDCALL *j1d)(DALocalInfo*,void*,void*,void*,int*);
 72: static int ourlj1d(DALocalInfo *info,PetscScalar *in,Mat m,void *ptr)
 73: {
 74:   int 0;
 75:   (*j1d)(info,&in[info->gxs],&m,ptr,&ierr);
 76:   return 0;
 77: }

 79: static void (PETSC_STDCALL *j2d)(DALocalInfo*,void*,void*,void*,int*);
 80: static int ourlj2d(DALocalInfo *info,PetscScalar **in,Mat m,void *ptr)
 81: {
 82:   int 0;
 83:   (*j2d)(info,&in[info->gys][info->gxs],&m,ptr,&ierr);
 84:   return 0;
 85: }

 87: static void (PETSC_STDCALL *j3d)(DALocalInfo*,void*,void*,void*,int*);
 88: static int ourlj3d(DALocalInfo *info,PetscScalar ***in,Mat m,void *ptr)
 89: {
 90:   int 0;
 91:   (*j3d)(info,&in[info->gzs][info->gys][info->gxs],&m,ptr,&ierr);
 92:   return 0;
 93: }

 95: static void (PETSC_STDCALL *f1d)(DALocalInfo*,void*,void*,void*,int*);
 96: static int ourlf1d(DALocalInfo *info,PetscScalar *in,PetscScalar *out,void *ptr)
 97: {
 98:   int 0;
 99:   (*f1d)(info,&in[info->gxs],&out[info->xs],ptr,&ierr);
100:   return 0;
101: }

103: static void (PETSC_STDCALL *f2d)(DALocalInfo*,void*,void*,void*,int*);
104: static int ourlf2d(DALocalInfo *info,PetscScalar **in,PetscScalar **out,void *ptr)
105: {
106:   int 0;
107:   (*f2d)(info,&in[info->gys][info->gxs],&out[info->ys][info->xs],ptr,&ierr);
108:   return 0;
109: }

111: static void (PETSC_STDCALL *f3d)(DALocalInfo*,void*,void*,void*,int*);
112: static int ourlf3d(DALocalInfo *info,PetscScalar ***in,PetscScalar ***out,void *ptr)
113: {
114:   int 0;
115:   (*f3d)(info,&in[info->gzs][info->gys][info->gxs],&out[info->zs][info->ys][info->xs],ptr,&ierr);
116:   return 0;
117: }

119: void PETSC_STDCALL dasetlocalfunction_(DA *da,void (PETSC_STDCALL *func)(DALocalInfo*,void*,void*,void*,int*),int *ierr)
120: {
121:   int dim;

123:   *DAGetInfo(*da,&dim,0,0,0,0,0,0,0,0,0,0); if (*ierr) return;
124:   if (dim == 2) {
125:      f2d    = (void (PETSC_STDCALL *)(DALocalInfo*,void*,void*,void*,int*))func;
126:     *DASetLocalFunction(*da,(DALocalFunction1)ourlf2d);
127:   } else if (dim == 3) {
128:      f3d    = (void (PETSC_STDCALL *)(DALocalInfo*,void*,void*,void*,int*))func;
129:     *DASetLocalFunction(*da,(DALocalFunction1)ourlf3d);
130:   } else if (dim == 1) {
131:      f1d    = (void (PETSC_STDCALL *)(DALocalInfo*,void*,void*,void*,int*))func;
132:     *DASetLocalFunction(*da,(DALocalFunction1)ourlf1d);
133:   } else *1;
134: }


137: void PETSC_STDCALL dasetlocaladiforfunction_(DA *da,
138: void (PETSC_STDCALL *jfunc)(int*,DALocalInfo*,void*,void*,int*,void*,void*,int*,void*,int*),int *ierr)
139: {
140:   (*da)->adifor_lf = (DALocalFunction1)jfunc;
141: }

143: void PETSC_STDCALL dasetlocaladiformffunction_(DA *da,
144: void (PETSC_STDCALL *jfunc)(DALocalInfo*,void*,void*,void*,void*,void*,int*),int *ierr)
145: {
146:   (*da)->adiformf_lf = (DALocalFunction1)jfunc;
147: }

149: void PETSC_STDCALL dasetlocaljacobian_(DA *da,void (PETSC_STDCALL *jac)(DALocalInfo*,void*,void*,void*,int*),int *ierr)
150: {
151:   int dim;

153:   *DAGetInfo(*da,&dim,0,0,0,0,0,0,0,0,0,0); if (*ierr) return;
154:   if (dim == 2) {
155:      j2d    = (void (PETSC_STDCALL *)(DALocalInfo*,void*,void*,void*,int*))jac;
156:     *DASetLocalJacobian(*da,(DALocalFunction1)ourlj2d);
157:   } else if (dim == 3) {
158:      j3d    = (void (PETSC_STDCALL *)(DALocalInfo*,void*,void*,void*,int*))jac;
159:     *DASetLocalJacobian(*da,(DALocalFunction1)ourlj3d);
160:   } else if (dim == 1) {
161:      j1d    = (void (PETSC_STDCALL *)(DALocalInfo*,void*,void*,void*,int*))jac;
162:     *DASetLocalJacobian(*da,(DALocalFunction1)ourlj1d);
163:   } else *1;
164: }

166: void PETSC_STDCALL dagetlocalinfo_(DA *da,DALocalInfo *ao,int *ierr)
167: {
168:   *DAGetLocalInfo(*da,ao);
169: }

171: void PETSC_STDCALL dagetao_(DA *da,AO *ao,int *ierr)
172: {
173:   *DAGetAO(*da,ao);
174: }

176: void PETSC_STDCALL darefine_(DA *da,MPI_Comm *comm,DA *daref, int *ierr)
177: {
178:   *DARefine(*da,(MPI_Comm)PetscToPointerComm(*comm),daref);
179: }

181: void PETSC_STDCALL dagetinterpolation_(DA *dac,DA *daf,Mat *A,Vec *scale,int *ierr)
182: {
183:   CHKFORTRANNULLOBJECT(scale);
184:   *DAGetInterpolation(*dac,*daf,A,scale);
185: }

187: void PETSC_STDCALL dasetfieldname_(DA *da,int *nf,CHAR name PETSC_MIXED_LEN(len),int *ierr PETSC_END_LEN(len))
188: {
189:   char *t;
190:   FIXCHAR(name,len,t);
191:   *DASetFieldName(*da,*nf,t);
192:   FREECHAR(name,t);
193: }
194: void PETSC_STDCALL dagetfieldname(DA *da,int *nf,CHAR name PETSC_MIXED_LEN(len),int *ierr PETSC_END_LEN(len))
195: {
196:   char *tname;

198:   *DAGetFieldName(*da,*nf,&tname);
199: #if defined(PETSC_USES_CPTOFCD)
200:   {
201:     char *t = _fcdtocp(name); int len1 = _fcdlen(name);
202:     *PetscStrncpy(t,tname,len1);
203:   }
204: #else
205:   *PetscStrncpy(name,tname,len);
206: #endif
207: }

209: void PETSC_STDCALL daload_(PetscViewer *viewer,int *M,int *N,int *P,DA *da,int *ierr)
210: {
211:   PetscViewer v;
212:   PetscPatchDefaultViewers_Fortran(viewer,v);
213:   *DALoad(v,*M,*N,*P,da);
214: }

216: void PETSC_STDCALL dagetislocaltoglobalmapping_(DA *da,ISLocalToGlobalMapping *map,int *ierr)
217: {
218:   *DAGetISLocalToGlobalMapping(*da,map);
219: }

221: void PETSC_STDCALL dagetislocaltoglobalmappingblck_(DA *da,ISLocalToGlobalMapping *map,int *ierr)
222: {
223:   *DAGetISLocalToGlobalMappingBlck(*da,map);
224: }

226: void PETSC_STDCALL dagetcoloring_(DA *da,ISColoringType *ctype,ISColoring *coloring,int *ierr PETSC_END_LEN(len))
227: {
228:   *DAGetColoring(*da,*ctype,coloring);
229: }

231: void PETSC_STDCALL dagetmatrix_(DA *da,CHAR mat_type PETSC_MIXED_LEN(len),Mat *J,int *ierr PETSC_END_LEN(len))
232: {
233:   char *t;
234:   FIXCHAR(mat_type,len,t);
235:   *DAGetMatrix(*da,t,J);
236:   FREECHAR(mat_type,t);
237: }

239: void PETSC_STDCALL daview_(DA *da,PetscViewer *vin,int *ierr)
240: {
241:   PetscViewer v;
242:   PetscPatchDefaultViewers_Fortran(vin,v);
243:   *DAView(*da,v);
244: }

246: void PETSC_STDCALL dagetglobalindices_(DA *da,int *n,int *indices,long *ia,int *ierr)
247: {
248:   int *idx;
249:   *DAGetGlobalIndices(*da,n,&idx);
250:   *ia   = PetscIntAddressToFortran(indices,idx);
251: }

253: void PETSC_STDCALL dacreateglobalvector_(DA *da,Vec* g,int *ierr)
254: {
255:   *DACreateGlobalVector(*da,g);
256: }

258: void PETSC_STDCALL dacreatelocalvector_(DA *da,Vec* l,int *ierr)
259: {
260:   *DACreateLocalVector(*da,l);
261: }

263: void PETSC_STDCALL dagetlocalvector_(DA *da,Vec* l,int *ierr)
264: {
265:   *DAGetLocalVector(*da,l);
266: }

268: void PETSC_STDCALL darestorelocalvector_(DA *da,Vec* l,int *ierr)
269: {
270:   *DARestoreLocalVector(*da,l);
271: }

273: void PETSC_STDCALL dagetscatter_(DA *da,VecScatter *ltog,VecScatter *gtol,VecScatter *ltol,int *ierr)
274: {
275:   CHKFORTRANNULLINTEGER(ltog);
276:   CHKFORTRANNULLINTEGER(gtol);
277:   CHKFORTRANNULLINTEGER(ltol);
278:   *DAGetScatter(*da,ltog,gtol,ltol);
279: }

281: void PETSC_STDCALL dadestroy_(DA *da,int *ierr)
282: {
283:   *DADestroy(*da);
284: }

286: void PETSC_STDCALL dacreate2d_(MPI_Comm *comm,DAPeriodicType *wrap,DAStencilType
287:                   *stencil_type,int *M,int *N,int *m,int *n,int *w,
288:                   int *s,int *lx,int *ly,DA *inra,int *ierr)
289: {
290:   CHKFORTRANNULLINTEGER(lx);
291:   CHKFORTRANNULLINTEGER(ly);
292:   *DACreate2d((MPI_Comm)PetscToPointerComm(*comm),*wrap,
293:                        *stencil_type,*M,*N,*m,*n,*w,*s,lx,ly,inra);
294: }

296: void PETSC_STDCALL dacreate1d_(MPI_Comm *comm,DAPeriodicType *wrap,int *M,int *w,int *s,
297:                  int *lc,DA *inra,int *ierr)
298: {
299:   CHKFORTRANNULLINTEGER(lc);
300:   *DACreate1d((MPI_Comm)PetscToPointerComm(*comm),*wrap,*M,*w,*s,lc,inra);
301: }

303: void PETSC_STDCALL dacreate3d_(MPI_Comm *comm,DAPeriodicType *wrap,DAStencilType 
304:                  *stencil_type,int *M,int *N,int *P,int *m,int *n,int *p,
305:                  int *w,int *s,int *lx,int *ly,int *lz,DA *inra,int *ierr)
306: {
307:   CHKFORTRANNULLINTEGER(lx);
308:   CHKFORTRANNULLINTEGER(ly);
309:   CHKFORTRANNULLINTEGER(lz);
310:   *DACreate3d((MPI_Comm)PetscToPointerComm(*comm),*wrap,*stencil_type,
311:                         *M,*N,*P,*m,*n,*p,*w,*s,lx,ly,lz,inra);
312: }

314: void PETSC_STDCALL dagetinfo_(DA *da,int *dim,int *M,int *N,int *P,int *m,int *n,int *p,int *w,int *s,
315:                 DAPeriodicType *wrap,DAStencilType *st,int *ierr)
316: {
317:   CHKFORTRANNULLINTEGER(dim);
318:   CHKFORTRANNULLINTEGER(M);
319:   CHKFORTRANNULLINTEGER(N);
320:   CHKFORTRANNULLINTEGER(P);
321:   CHKFORTRANNULLINTEGER(m);
322:   CHKFORTRANNULLINTEGER(n);
323:   CHKFORTRANNULLINTEGER(p);
324:   CHKFORTRANNULLINTEGER(w);
325:   CHKFORTRANNULLINTEGER(s);
326:   CHKFORTRANNULLINTEGER(wrap);
327:   CHKFORTRANNULLINTEGER(st);
328:   *DAGetInfo(*da,dim,M,N,P,m,n,p,w,s,wrap,st);
329: }

331: EXTERN_C_END