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 dacreatenaturalvector_       DACREATENATURALVECTOR
 21: #define dacreatelocalvector_         DACREATELOCALVECTOR
 22: #define dagetlocalvector_            DAGETLOCALVECTOR
 23: #define darestorelocalvector_        DARESTORELOCALVECTOR
 24: #define dagetscatter_                DAGETSCATTER
 25: #define dagetglobalindices_          DAGETGLOBALINDICES
 26: #define daview_                      DAVIEW
 27: #define dagetinfo_                   DAGETINFO
 28: #define dagetcoloring_               DAGETCOLORING
 29: #define dagetmatrix_                 DAGETMATRIX
 30: #define dagetislocaltoglobalmapping_ DAGETISLOCALTOGLOBALMAPPING
 31: #define dagetislocaltoglobalmappingblck_ DAGETISLOCALTOGLOBALMAPPINGBLCK
 32: #define daload_                      DALOAD
 33: #define dasetfieldname_              DASETFIELDNAME
 34: #define dagetfieldname_              DAGETFIELDNAME
 35: #define darefine_                    DAREFINE
 36: #define dagetao_                     DAGETAO
 37: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 38: #define dagetlocalinfo_              dagetlocalinfo
 39: #define dagetlocalvector_            dagetlocalvector
 40: #define darestorelocalvector_        darestorelocalvector
 41: #define dagetinterpolation_          dagetinterpolation
 42: #define daload_                      daload
 43: #define dacreateglobalvector_        dacreateglobalvector
 44: #define dacreatenaturalvector_       dacreatenaturalvector
 45: #define dacreatelocalvector_         dacreatelocalvector
 46: #define daview_                      daview
 47: #define dacreate1d_                  dacreate1d
 48: #define dacreate3d_                  dacreate3d
 49: #define dacreate2d_                  dacreate2d
 50: #define dadestroy_                   dadestroy
 51: #define dagetscatter_                dagetscatter
 52: #define dagetglobalindices_          dagetglobalindices
 53: #define dagetinfo_                   dagetinfo
 54: #define dagetcoloring_               dagetcoloring
 55: #define dagetmatrix_                 dagetmatrix
 56: #define dagetislocaltoglobalmapping_ dagetislocaltoglobalmapping
 57: #define dagetislocaltoglobalmappingblck_ dagetislocaltoglobalmappingblck
 58: #define dasetfieldname_              dasetfieldname
 59: #define dagetfieldname_              dagetfieldname
 60: #define darefine_                    darefine
 61: #define dagetao_                     dagetao
 62: #define dasetlocalfunction_          dasetlocalfunction
 63: #define dasetlocaladiforfunction_       dasetlocaladiforfunction
 64: #define dasetlocaladiformffunction_       dasetlocaladiformffunction
 65: #define dasetlocaljacobian_          dasetlocaljacobian
 66: #endif



 70: EXTERN_C_BEGIN


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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

260: void PETSC_STDCALL dacreatenaturalvector_(DA *da,Vec* g,int *ierr)
261: {
262:   *DACreateNaturalVector(*da,g);
263: }

265: void PETSC_STDCALL dacreatelocalvector_(DA *da,Vec* l,int *ierr)
266: {
267:   *DACreateLocalVector(*da,l);
268: }

270: void PETSC_STDCALL dagetlocalvector_(DA *da,Vec* l,int *ierr)
271: {
272:   *DAGetLocalVector(*da,l);
273: }

275: void PETSC_STDCALL darestorelocalvector_(DA *da,Vec* l,int *ierr)
276: {
277:   *DARestoreLocalVector(*da,l);
278: }

280: void PETSC_STDCALL dagetscatter_(DA *da,VecScatter *ltog,VecScatter *gtol,VecScatter *ltol,int *ierr)
281: {
282:   CHKFORTRANNULLINTEGER(ltog);
283:   CHKFORTRANNULLINTEGER(gtol);
284:   CHKFORTRANNULLINTEGER(ltol);
285:   *DAGetScatter(*da,ltog,gtol,ltol);
286: }

288: void PETSC_STDCALL dadestroy_(DA *da,int *ierr)
289: {
290:   *DADestroy(*da);
291: }

293: void PETSC_STDCALL dacreate2d_(MPI_Comm *comm,DAPeriodicType *wrap,DAStencilType
294:                   *stencil_type,int *M,int *N,int *m,int *n,int *w,
295:                   int *s,int *lx,int *ly,DA *inra,int *ierr)
296: {
297:   CHKFORTRANNULLINTEGER(lx);
298:   CHKFORTRANNULLINTEGER(ly);
299:   *DACreate2d((MPI_Comm)PetscToPointerComm(*comm),*wrap,
300:                        *stencil_type,*M,*N,*m,*n,*w,*s,lx,ly,inra);
301: }

303: void PETSC_STDCALL dacreate1d_(MPI_Comm *comm,DAPeriodicType *wrap,int *M,int *w,int *s,
304:                  int *lc,DA *inra,int *ierr)
305: {
306:   CHKFORTRANNULLINTEGER(lc);
307:   *DACreate1d((MPI_Comm)PetscToPointerComm(*comm),*wrap,*M,*w,*s,lc,inra);
308: }

310: void PETSC_STDCALL dacreate3d_(MPI_Comm *comm,DAPeriodicType *wrap,DAStencilType 
311:                  *stencil_type,int *M,int *N,int *P,int *m,int *n,int *p,
312:                  int *w,int *s,int *lx,int *ly,int *lz,DA *inra,int *ierr)
313: {
314:   CHKFORTRANNULLINTEGER(lx);
315:   CHKFORTRANNULLINTEGER(ly);
316:   CHKFORTRANNULLINTEGER(lz);
317:   *DACreate3d((MPI_Comm)PetscToPointerComm(*comm),*wrap,*stencil_type,
318:                         *M,*N,*P,*m,*n,*p,*w,*s,lx,ly,lz,inra);
319: }

321: void PETSC_STDCALL dagetinfo_(DA *da,int *dim,int *M,int *N,int *P,int *m,int *n,int *p,int *w,int *s,
322:                 DAPeriodicType *wrap,DAStencilType *st,int *ierr)
323: {
324:   CHKFORTRANNULLINTEGER(dim);
325:   CHKFORTRANNULLINTEGER(M);
326:   CHKFORTRANNULLINTEGER(N);
327:   CHKFORTRANNULLINTEGER(P);
328:   CHKFORTRANNULLINTEGER(m);
329:   CHKFORTRANNULLINTEGER(n);
330:   CHKFORTRANNULLINTEGER(p);
331:   CHKFORTRANNULLINTEGER(w);
332:   CHKFORTRANNULLINTEGER(s);
333:   CHKFORTRANNULLINTEGER(wrap);
334:   CHKFORTRANNULLINTEGER(st);
335:   *DAGetInfo(*da,dim,M,N,P,m,n,p,w,s,wrap,st);
336: }

338: EXTERN_C_END