Actual source code: f90_cwrap.c

petsc-3.3-p0 2012-06-05
  1: #include <../src/sys/f90-src/f90impl.h>

  3: /*************************************************************************/

  5: #if defined(PETSC_HAVE_FORTRAN_CAPS)
  6: #define f90array1dcreatescalar_           F90ARRAY1DCREATESCALAR
  7: #define f90array1daccessscalar_           F90ARRAY1DACCESSSCALAR
  8: #define f90array1ddestroyscalar_          F90ARRAY1DDESTROYSCALAR
  9: #define f90array1dcreatereal_             F90ARRAY1DCREATEREAL
 10: #define f90array1daccessreal_             F90ARRAY1DACCESSREAL
 11: #define f90array1ddestroyreal_            F90ARRAY1DDESTROYREAL
 12: #define f90array1dcreateint_              F90ARRAY1DCREATEINT
 13: #define f90array1daccessint_              F90ARRAY1DACCESSINT
 14: #define f90array1ddestroyint_             F90ARRAY1DDESTROYINT
 15: #define f90array1dcreatefortranaddr_      F90ARRAY1DCREATEFORTRANADDR
 16: #define f90array1daccessfortranaddr_      F90ARRAY1DACCESSFORTRANADDR
 17: #define f90array1ddestroyfortranaddr_     F90ARRAY1DDESTROYFORTRANADDR
 18: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 19: #define f90array1dcreatescalar_           f90array1dcreatescalar
 20: #define f90array1daccessscalar_           f90array1daccessscalar
 21: #define f90array1ddestroyscalar_          f90array1ddestroyscalar
 22: #define f90array1dcreatereal_             f90array1dcreatereal
 23: #define f90array1daccessreal_             f90array1daccessreal
 24: #define f90array1ddestroyreal_            f90array1ddestroyreal
 25: #define f90array1dcreateint_              f90array1dcreateint
 26: #define f90array1daccessint_              f90array1daccessint
 27: #define f90array1ddestroyint_             f90array1ddestroyint
 28: #define f90array1dcreatefortranaddr_      f90array1dcreatefortranaddr
 29: #define f90array1daccessfortranaddr_      f90array1daccessfortranaddr
 30: #define f90array1ddestroyfortranaddr_     f90array1ddestroyfortranaddr
 31: #endif


 48: PetscErrorCode F90Array1dCreate(void *array,PetscDataType type,PetscInt start,PetscInt len,F90Array1d *ptr PETSC_F90_2PTR_PROTO(ptrd))
 49: {
 51:   if (type == PETSC_SCALAR) {
 52:     if (!len) array = PETSC_NULL_SCALAR_Fortran;
 53:     f90array1dcreatescalar_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
 54:   } else if (type == PETSC_REAL) {
 55:     if (!len) array = PETSC_NULL_REAL_Fortran;
 56:     f90array1dcreatereal_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
 57:   } else if (type == PETSC_INT) {
 58:     if (!len) array = PETSC_NULL_INTEGER_Fortran;
 59:     f90array1dcreateint_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
 60:   } else if (type == PETSC_FORTRANADDR) {
 61:     f90array1dcreatefortranaddr_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
 62:   } else {
 63:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported PetscDataType: %d",(PetscInt)type);
 64:   }
 65:   return(0);
 66: }

 70: PetscErrorCode  F90Array1dAccess(F90Array1d *ptr,PetscDataType type,void **array PETSC_F90_2PTR_PROTO(ptrd))
 71: {
 73:   if (type == PETSC_SCALAR) {
 74:     f90array1daccessscalar_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
 75:     if (*array == PETSC_NULL_SCALAR_Fortran) *array = 0;
 76:   } else if (type == PETSC_REAL) {
 77:     f90array1daccessreal_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
 78:     if (*array == PETSC_NULL_REAL_Fortran) *array = 0;
 79:   } else if (type == PETSC_INT) {
 80:     f90array1daccessint_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
 81:     if (*array == PETSC_NULL_INTEGER_Fortran) *array = 0;
 82:   } else if (type == PETSC_FORTRANADDR) {
 83:     f90array1daccessfortranaddr_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
 84:   } else {
 85:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported PetscDataType: %d",(PetscInt)type);
 86:   }
 87:   return(0);
 88: }

 92: PetscErrorCode  F90Array1dDestroy(F90Array1d *ptr,PetscDataType type PETSC_F90_2PTR_PROTO(ptrd))
 93: {
 95:   if (type == PETSC_SCALAR) {
 96:     f90array1ddestroyscalar_(ptr PETSC_F90_2PTR_PARAM(ptrd));
 97:   } else if (type == PETSC_REAL) {
 98:     f90array1ddestroyreal_(ptr PETSC_F90_2PTR_PARAM(ptrd));
 99:   } else if (type == PETSC_INT) {
100:     f90array1ddestroyint_(ptr PETSC_F90_2PTR_PARAM(ptrd));
101:   } else if (type == PETSC_FORTRANADDR) {
102:     f90array1ddestroyfortranaddr_(ptr PETSC_F90_2PTR_PARAM(ptrd));
103:   } else {
104:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported PetscDataType: %d",(PetscInt)type);
105:   }
106:   return(0);
107: }

109: /*************************************************************************/

111: #if defined(PETSC_HAVE_FORTRAN_CAPS)
112: #define f90array2dcreatescalar_           F90ARRAY2DCREATESCALAR
113: #define f90array2daccessscalar_           F90ARRAY2DACCESSSCALAR
114: #define f90array2ddestroyscalar_          F90ARRAY2DDESTROYSCALAR
115: #define f90array2dcreatereal_             F90ARRAY2DCREATEREAL
116: #define f90array2daccessreal_             F90ARRAY2DACCESSREAL
117: #define f90array2ddestroyreal_            F90ARRAY2DDESTROYREAL
118: #define f90array2dcreateint_              F90ARRAY2DCREATEINT
119: #define f90array2daccessint_              F90ARRAY2DACCESSINT
120: #define f90array2ddestroyint_             F90ARRAY2DDESTROYINT
121: #define f90array2dcreatefortranaddr_      F90ARRAY2DCREATEFORTRANADDR
122: #define f90array2daccessfortranaddr_      F90ARRAY2DACCESSFORTRANADDR
123: #define f90array2ddestroyfortranaddr_     F90ARRAY2DDESTROYFORTRANADDR
124: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
125: #define f90array2dcreatescalar_           f90array2dcreatescalar
126: #define f90array2daccessscalar_           f90array2daccessscalar
127: #define f90array2ddestroyscalar_          f90array2ddestroyscalar
128: #define f90array2dcreatereal_             f90array2dcreatereal
129: #define f90array2daccessreal_             f90array2daccessreal
130: #define f90array2ddestroyreal_            f90array2ddestroyreal
131: #define f90array2dcreateint_              f90array2dcreateint
132: #define f90array2daccessint_              f90array2daccessint
133: #define f90array2ddestroyint_             f90array2ddestroyint
134: #define f90array2dcreatefortranaddr_      f90array2dcreatefortranaddr
135: #define f90array2daccessfortranaddr_      f90array2daccessfortranaddr
136: #define f90array2ddestroyfortranaddr_     f90array2ddestroyfortranaddr
137: #endif


154: PetscErrorCode F90Array2dCreate(void *array,PetscDataType type,PetscInt start1,PetscInt len1,PetscInt start2,PetscInt len2,F90Array2d *ptr PETSC_F90_2PTR_PROTO(ptrd))
155: {
157:   if (type == PETSC_SCALAR) {
158:     f90array2dcreatescalar_(array,&start1,&len1,&start2,&len2,ptr PETSC_F90_2PTR_PARAM(ptrd));
159:   } else if (type == PETSC_REAL) {
160:     f90array2dcreatereal_(array,&start1,&len1,&start2,&len2,ptr PETSC_F90_2PTR_PARAM(ptrd));
161:   } else if (type == PETSC_INT) {
162:     f90array2dcreateint_(array,&start1,&len1,&start2,&len2,ptr PETSC_F90_2PTR_PARAM(ptrd));
163:   } else if (type == PETSC_FORTRANADDR) {
164:     f90array2dcreatefortranaddr_(array,&start1,&len1,&start2,&len2,ptr PETSC_F90_2PTR_PARAM(ptrd));
165:   } else {
166:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported PetscDataType: %d",(PetscInt)type);
167:   }
168:   return(0);
169: }

173: PetscErrorCode  F90Array2dAccess(F90Array2d *ptr,PetscDataType type,void **array PETSC_F90_2PTR_PROTO(ptrd))
174: {
176:   if (type == PETSC_SCALAR) {
177:     f90array2daccessscalar_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
178:   } else if (type == PETSC_REAL) {
179:     f90array2daccessreal_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
180:   } else if (type == PETSC_INT) {
181:     f90array2daccessint_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
182:   } else if (type == PETSC_FORTRANADDR) {
183:     f90array2daccessfortranaddr_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
184:   } else {
185:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported PetscDataType: %d",(PetscInt)type);
186:   }
187:   return(0);
188: }

192: PetscErrorCode  F90Array2dDestroy(F90Array2d *ptr,PetscDataType type PETSC_F90_2PTR_PROTO(ptrd))
193: {
195:   if (type == PETSC_SCALAR) {
196:     f90array2ddestroyscalar_(ptr PETSC_F90_2PTR_PARAM(ptrd));
197:   } else if (type == PETSC_REAL) {
198:     f90array2ddestroyreal_(ptr PETSC_F90_2PTR_PARAM(ptrd));
199:   } else if (type == PETSC_INT) {
200:     f90array2ddestroyint_(ptr PETSC_F90_2PTR_PARAM(ptrd));
201:   } else if (type == PETSC_FORTRANADDR) {
202:     f90array2ddestroyfortranaddr_(ptr PETSC_F90_2PTR_PARAM(ptrd));
203:   } else {
204:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported PetscDataType: %d",(PetscInt)type);
205:   }
206:   return(0);
207: }

209: /*************************************************************************/

211: #if defined(PETSC_HAVE_FORTRAN_CAPS)
212: #define f90array3dcreatescalar_           F90ARRAY3DCREATESCALAR
213: #define f90array3daccessscalar_           F90ARRAY3DACCESSSCALAR
214: #define f90array3ddestroyscalar_          F90ARRAY3DDESTROYSCALAR
215: #define f90array3dcreatereal_             F90ARRAY3DCREATEREAL
216: #define f90array3daccessreal_             F90ARRAY3DACCESSREAL
217: #define f90array3ddestroyreal_            F90ARRAY3DDESTROYREAL
218: #define f90array3dcreateint_              F90ARRAY3DCREATEINT
219: #define f90array3daccessint_              F90ARRAY3DACCESSINT
220: #define f90array3ddestroyint_             F90ARRAY3DDESTROYINT
221: #define f90array3dcreatefortranaddr_      F90ARRAY3DCREATEFORTRANADDR
222: #define f90array3daccessfortranaddr_      F90ARRAY3DACCESSFORTRANADDR
223: #define f90array3ddestroyfortranaddr_     F90ARRAY3DDESTROYFORTRANADDR
224: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
225: #define f90array3dcreatescalar_           f90array3dcreatescalar
226: #define f90array3daccessscalar_           f90array3daccessscalar
227: #define f90array3ddestroyscalar_          f90array3ddestroyscalar
228: #define f90array3dcreatereal_             f90array3dcreatereal
229: #define f90array3daccessreal_             f90array3daccessreal
230: #define f90array3ddestroyreal_            f90array3ddestroyreal
231: #define f90array3dcreateint_              f90array3dcreateint
232: #define f90array3daccessint_              f90array3daccessint
233: #define f90array3ddestroyint_             f90array3ddestroyint
234: #define f90array3dcreatefortranaddr_      f90array3dcreatefortranaddr
235: #define f90array3daccessfortranaddr_      f90array3daccessfortranaddr
236: #define f90array3ddestroyfortranaddr_     f90array3ddestroyfortranaddr
237: #endif


254: PetscErrorCode F90Array3dCreate(void *array,PetscDataType type,PetscInt start1,PetscInt len1,PetscInt start2,PetscInt len2,PetscInt start3,PetscInt len3,F90Array3d *ptr PETSC_F90_2PTR_PROTO(ptrd))
255: {
257:   if (type == PETSC_SCALAR) {
258:     f90array3dcreatescalar_(array,&start1,&len1,&start2,&len2,&start3,&len3,ptr PETSC_F90_2PTR_PARAM(ptrd));
259:   } else if (type == PETSC_REAL) {
260:     f90array3dcreatereal_(array,&start1,&len1,&start2,&len2,&start3,&len3,ptr PETSC_F90_2PTR_PARAM(ptrd));
261:   } else if (type == PETSC_INT) {
262:     f90array3dcreateint_(array,&start1,&len1,&start2,&len2,&start3,&len3,ptr PETSC_F90_2PTR_PARAM(ptrd));
263:   } else if (type == PETSC_FORTRANADDR) {
264:     f90array3dcreatefortranaddr_(array,&start1,&len1,&start2,&len2,&start3,&len3,ptr PETSC_F90_2PTR_PARAM(ptrd));
265:   } else {
266:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported PetscDataType: %d",(PetscInt)type);
267:   }
268:   return(0);
269: }

273: PetscErrorCode  F90Array3dAccess(F90Array3d *ptr,PetscDataType type,void **array PETSC_F90_2PTR_PROTO(ptrd))
274: {
276:   if (type == PETSC_SCALAR) {
277:     f90array3daccessscalar_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
278:   } else if (type == PETSC_REAL) {
279:     f90array3daccessreal_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
280:   } else if (type == PETSC_INT) {
281:     f90array3daccessint_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
282:   } else if (type == PETSC_FORTRANADDR) {
283:     f90array3daccessfortranaddr_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
284:   } else {
285:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported PetscDataType: %d",(PetscInt)type);
286:   }
287:   return(0);
288: }

292: PetscErrorCode  F90Array3dDestroy(F90Array3d *ptr,PetscDataType type PETSC_F90_2PTR_PROTO(ptrd))
293: {
295:   if (type == PETSC_SCALAR) {
296:     f90array3ddestroyscalar_(ptr PETSC_F90_2PTR_PARAM(ptrd));
297:   } else if (type == PETSC_REAL) {
298:     f90array3ddestroyreal_(ptr PETSC_F90_2PTR_PARAM(ptrd));
299:   } else if (type == PETSC_INT) {
300:     f90array3ddestroyint_(ptr PETSC_F90_2PTR_PARAM(ptrd));
301:   } else if (type == PETSC_FORTRANADDR) {
302:     f90array3ddestroyfortranaddr_(ptr PETSC_F90_2PTR_PARAM(ptrd));
303:   } else {
304:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported PetscDataType: %d",(PetscInt)type);
305:   }
306:   return(0);
307: }

309: /*************************************************************************/

311: #if defined(PETSC_HAVE_FORTRAN_CAPS)
312: #define f90array4dcreatescalar_           F90ARRAY4DCREATESCALAR
313: #define f90array4ddestroyscalar_          F90ARRAY4DDESTROYSCALAR
314: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
315: #define f90array4dcreatescalar_           f90array4dcreatescalar
316: #define f90array4ddestroyscalar_          f90array4ddestroyscalar
317: #endif


324: PetscErrorCode F90Array4dCreate(void *array,PetscDataType type,PetscInt start1,PetscInt len1,PetscInt start2,PetscInt len2,PetscInt start3,PetscInt len3,PetscInt start4,PetscInt len4,F90Array4d *ptr PETSC_F90_2PTR_PROTO(ptrd))
325: {
327:   if (type == PETSC_SCALAR) {
328:     f90array4dcreatescalar_(array,&start1,&len1,&start2,&len2,&start3,&len3,&start4,&len4,ptr PETSC_F90_2PTR_PARAM(ptrd));
329:   } else {
330:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported PetscDataType: %d",(PetscInt)type);
331:   }
332:   return(0);
333: }

337: PetscErrorCode  F90Array4dDestroy(F90Array4d *ptr,PetscDataType type PETSC_F90_2PTR_PROTO(ptrd))
338: {
340:   if (type == PETSC_SCALAR) {
341:     f90array4ddestroyscalar_(ptr PETSC_F90_2PTR_PARAM(ptrd));
342:   } else {
343:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unsupported PetscDataType: %d",(PetscInt)type);
344:   }
345:   return(0);
346: }

348: /*************************************************************************/
349: #if defined(PETSC_HAVE_FORTRAN_CAPS)
350: #define f90array1dgetaddrscalar_            F90ARRAY1DGETADDRSCALAR
351: #define f90array1dgetaddrreal_              F90ARRAY1DGETADDRREAL
352: #define f90array1dgetaddrint_               F90ARRAY1DGETADDRINT
353: #define f90array1dgetaddrfortranaddr_       F90ARRAY1DGETADDRFORTRANADDR
354: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
355: #define f90array1dgetaddrscalar_            f90array1dgetaddrscalar
356: #define f90array1dgetaddrreal_              f90array1dgetaddrreal
357: #define f90array1dgetaddrint_               f90array1dgetaddrint
358: #define f90array1dgetaddrfortranaddr_       f90array1dgetaddrfortranaddr
359: #endif

362: {
363:   *address = (PetscFortranAddr)array;
364: }
366: {
367:   *address = (PetscFortranAddr)array;
368: }
370: {
371:   *address = (PetscFortranAddr)array;
372: }
374: {
375:   *address = (PetscFortranAddr)array;
376: }

378: /*************************************************************************/
379: #if defined(PETSC_HAVE_FORTRAN_CAPS)
380: #define f90array2dgetaddrscalar_            F90ARRAY2DGETADDRSCALAR
381: #define f90array2dgetaddrreal_              F90ARRAY2DGETADDRREAL
382: #define f90array2dgetaddrint_               F90ARRAY2DGETADDRINT
383: #define f90array2dgetaddrfortranaddr_       F90ARRAY2DGETADDRFORTRANADDR
384: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
385: #define f90array2dgetaddrscalar_            f90array2dgetaddrscalar
386: #define f90array2dgetaddrreal_              f90array2dgetaddrreal
387: #define f90array2dgetaddrint_               f90array2dgetaddrint
388: #define f90array2dgetaddrfortranaddr_       f90array2dgetaddrfortranaddr
389: #endif

392: {
393:   *address = (PetscFortranAddr)array;
394: }
396: {
397:   *address = (PetscFortranAddr)array;
398: }
400: {
401:   *address = (PetscFortranAddr)array;
402: }
404: {
405:   *address = (PetscFortranAddr)array;
406: }

408: /*************************************************************************/
409: #if defined(PETSC_HAVE_FORTRAN_CAPS)
410: #define f90array3dgetaddrscalar_            F90ARRAY3DGETADDRSCALAR
411: #define f90array3dgetaddrreal_              F90ARRAY3DGETADDRREAL
412: #define f90array3dgetaddrint_               F90ARRAY3DGETADDRINT
413: #define f90array3dgetaddrfortranaddr_       F90ARRAY3DGETADDRFORTRANADDR
414: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
415: #define f90array3dgetaddrscalar_            f90array3dgetaddrscalar
416: #define f90array3dgetaddrreal_              f90array3dgetaddrreal
417: #define f90array3dgetaddrint_               f90array3dgetaddrint
418: #define f90array3dgetaddrfortranaddr_       f90array3dgetaddrfortranaddr
419: #endif

422: {
423:   *address = (PetscFortranAddr)array;
424: }
426: {
427:   *address = (PetscFortranAddr)array;
428: }
430: {
431:   *address = (PetscFortranAddr)array;
432: }
434: {
435:   *address = (PetscFortranAddr)array;
436: }

438: /*************************************************************************/