Actual source code: f90_cwrap.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: #include <petsc/private/f90impl.h>

  3: /*@C

  5:    Converts a MPI_Fint that contains a Fortran MPI_Datatype to its C MPI_Datatype equivalent

  7:    Not Collective

  9:    Input Parameter:
 10: .  unit - The Fortran MPI_Datatype

 12:    Output Parameter:
 13: .  dtype - the corresponding C MPI_Datatype

 15:    Level: developer

 17:    Developer Notes:
 18:     The MPI documentation in multiple places says that one can never us
 19:    Fortran MPI_Datatypes in C (or vis-versa) but this is problematic since users could never
 20:    call C routines from Fortran that have MPI_Datatype arguments. Jed states that the Fortran
 21:    MPI_Datatypes will always be available in C if the MPI was built to support Fortran. This function
 22:    relys on this.
 23: @*/
 24: PetscErrorCode PetscMPIFortranDatatypeToC(MPI_Fint unit,MPI_Datatype *dtype)
 25: {
 26:   MPI_Datatype ftype;

 29:   ftype = MPI_Type_f2c(unit);
 30:   if (ftype == MPI_INTEGER) *dtype = MPI_INT;
 31:   else if (ftype == MPI_INTEGER8) *dtype = MPIU_INT64;
 32:   else if (ftype == MPI_DOUBLE_PRECISION) *dtype = MPI_DOUBLE;
 33: #if defined(PETSC_HAVE_COMPLEX)
 34:   else if (ftype == MPI_COMPLEX16) *dtype = MPIU_C_DOUBLE_COMPLEX;
 35: #endif
 36:   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown Fortran MPI_Datatype");
 37:   return(0);
 38: }

 40: /*************************************************************************/

 42: #if defined(PETSC_HAVE_FORTRAN_CAPS)
 43: #define f90array1dcreatescalar_           F90ARRAY1DCREATESCALAR
 44: #define f90array1daccessscalar_           F90ARRAY1DACCESSSCALAR
 45: #define f90array1ddestroyscalar_          F90ARRAY1DDESTROYSCALAR
 46: #define f90array1dcreatereal_             F90ARRAY1DCREATEREAL
 47: #define f90array1daccessreal_             F90ARRAY1DACCESSREAL
 48: #define f90array1ddestroyreal_            F90ARRAY1DDESTROYREAL
 49: #define f90array1dcreateint_              F90ARRAY1DCREATEINT
 50: #define f90array1daccessint_              F90ARRAY1DACCESSINT
 51: #define f90array1ddestroyint_             F90ARRAY1DDESTROYINT
 52: #define f90array1dcreatefortranaddr_      F90ARRAY1DCREATEFORTRANADDR
 53: #define f90array1daccessfortranaddr_      F90ARRAY1DACCESSFORTRANADDR
 54: #define f90array1ddestroyfortranaddr_     F90ARRAY1DDESTROYFORTRANADDR
 55: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
 56: #define f90array1dcreatescalar_           f90array1dcreatescalar
 57: #define f90array1daccessscalar_           f90array1daccessscalar
 58: #define f90array1ddestroyscalar_          f90array1ddestroyscalar
 59: #define f90array1dcreatereal_             f90array1dcreatereal
 60: #define f90array1daccessreal_             f90array1daccessreal
 61: #define f90array1ddestroyreal_            f90array1ddestroyreal
 62: #define f90array1dcreateint_              f90array1dcreateint
 63: #define f90array1daccessint_              f90array1daccessint
 64: #define f90array1ddestroyint_             f90array1ddestroyint
 65: #define f90array1dcreatefortranaddr_      f90array1dcreatefortranaddr
 66: #define f90array1daccessfortranaddr_      f90array1daccessfortranaddr
 67: #define f90array1ddestroyfortranaddr_     f90array1ddestroyfortranaddr
 68: #endif

 70: PETSC_EXTERN void PETSC_STDCALL f90array1dcreatescalar_(void *,PetscInt *,PetscInt *,F90Array1d * PETSC_F90_2PTR_PROTO_NOVAR);
 71: PETSC_EXTERN void PETSC_STDCALL f90array1daccessscalar_(F90Array1d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
 72: PETSC_EXTERN void PETSC_STDCALL f90array1ddestroyscalar_(F90Array1d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
 73: PETSC_EXTERN void PETSC_STDCALL f90array1dcreatereal_(void *,PetscInt *,PetscInt *,F90Array1d * PETSC_F90_2PTR_PROTO_NOVAR);
 74: PETSC_EXTERN void PETSC_STDCALL f90array1daccessreal_(F90Array1d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
 75: PETSC_EXTERN void PETSC_STDCALL f90array1ddestroyreal_(F90Array1d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
 76: PETSC_EXTERN void PETSC_STDCALL f90array1dcreateint_(void *,PetscInt *,PetscInt *,F90Array1d * PETSC_F90_2PTR_PROTO_NOVAR);
 77: PETSC_EXTERN void PETSC_STDCALL f90array1daccessint_(F90Array1d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
 78: PETSC_EXTERN void PETSC_STDCALL f90array1ddestroyint_(F90Array1d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
 79: PETSC_EXTERN void PETSC_STDCALL f90array1dcreatefortranaddr_(void *,PetscInt *,PetscInt *,F90Array1d * PETSC_F90_2PTR_PROTO_NOVAR);
 80: PETSC_EXTERN void PETSC_STDCALL f90array1daccessfortranaddr_(F90Array1d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
 81: PETSC_EXTERN void PETSC_STDCALL f90array1ddestroyfortranaddr_(F90Array1d *ptr PETSC_F90_2PTR_PROTO_NOVAR);

 83: PetscErrorCode F90Array1dCreate(void *array,MPI_Datatype type,PetscInt start,PetscInt len,F90Array1d *ptr PETSC_F90_2PTR_PROTO(ptrd))
 84: {
 86:   if (type == MPIU_SCALAR) {
 87:     if (!len) array = PETSC_NULL_SCALAR_Fortran;
 88:     f90array1dcreatescalar_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
 89:   } else if (type == MPIU_REAL) {
 90:     if (!len) array = PETSC_NULL_REAL_Fortran;
 91:     f90array1dcreatereal_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
 92:   } else if (type == MPIU_INT) {
 93:     if (!len) array = PETSC_NULL_INTEGER_Fortran;
 94:     f90array1dcreateint_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
 95:   } else if (type == MPIU_FORTRANADDR) {
 96:     f90array1dcreatefortranaddr_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
 97:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
 98:   return(0);
 99: }

101: PetscErrorCode  F90Array1dAccess(F90Array1d *ptr,MPI_Datatype type,void **array PETSC_F90_2PTR_PROTO(ptrd))
102: {
104:   if (type == MPIU_SCALAR) {
105:     f90array1daccessscalar_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
106:     if (*array == PETSC_NULL_SCALAR_Fortran) *array = 0;
107:   } else if (type == MPIU_REAL) {
108:     f90array1daccessreal_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
109:     if (*array == PETSC_NULL_REAL_Fortran) *array = 0;
110:   } else if (type == MPIU_INT) {
111:     f90array1daccessint_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
112:     if (*array == PETSC_NULL_INTEGER_Fortran) *array = 0;
113:   } else if (type == MPIU_FORTRANADDR) {
114:     f90array1daccessfortranaddr_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
115:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
116:   return(0);
117: }

119: PetscErrorCode  F90Array1dDestroy(F90Array1d *ptr,MPI_Datatype type PETSC_F90_2PTR_PROTO(ptrd))
120: {
122:   if (type == MPIU_SCALAR) {
123:     f90array1ddestroyscalar_(ptr PETSC_F90_2PTR_PARAM(ptrd));
124:   } else if (type == MPIU_REAL) {
125:     f90array1ddestroyreal_(ptr PETSC_F90_2PTR_PARAM(ptrd));
126:   } else if (type == MPIU_INT) {
127:     f90array1ddestroyint_(ptr PETSC_F90_2PTR_PARAM(ptrd));
128:   } else if (type == MPIU_FORTRANADDR) {
129:     f90array1ddestroyfortranaddr_(ptr PETSC_F90_2PTR_PARAM(ptrd));
130:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
131:   return(0);
132: }

134: /*************************************************************************/

136: #if defined(PETSC_HAVE_FORTRAN_CAPS)
137: #define f90array2dcreatescalar_           F90ARRAY2DCREATESCALAR
138: #define f90array2daccessscalar_           F90ARRAY2DACCESSSCALAR
139: #define f90array2ddestroyscalar_          F90ARRAY2DDESTROYSCALAR
140: #define f90array2dcreatereal_             F90ARRAY2DCREATEREAL
141: #define f90array2daccessreal_             F90ARRAY2DACCESSREAL
142: #define f90array2ddestroyreal_            F90ARRAY2DDESTROYREAL
143: #define f90array2dcreateint_              F90ARRAY2DCREATEINT
144: #define f90array2daccessint_              F90ARRAY2DACCESSINT
145: #define f90array2ddestroyint_             F90ARRAY2DDESTROYINT
146: #define f90array2dcreatefortranaddr_      F90ARRAY2DCREATEFORTRANADDR
147: #define f90array2daccessfortranaddr_      F90ARRAY2DACCESSFORTRANADDR
148: #define f90array2ddestroyfortranaddr_     F90ARRAY2DDESTROYFORTRANADDR
149: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
150: #define f90array2dcreatescalar_           f90array2dcreatescalar
151: #define f90array2daccessscalar_           f90array2daccessscalar
152: #define f90array2ddestroyscalar_          f90array2ddestroyscalar
153: #define f90array2dcreatereal_             f90array2dcreatereal
154: #define f90array2daccessreal_             f90array2daccessreal
155: #define f90array2ddestroyreal_            f90array2ddestroyreal
156: #define f90array2dcreateint_              f90array2dcreateint
157: #define f90array2daccessint_              f90array2daccessint
158: #define f90array2ddestroyint_             f90array2ddestroyint
159: #define f90array2dcreatefortranaddr_      f90array2dcreatefortranaddr
160: #define f90array2daccessfortranaddr_      f90array2daccessfortranaddr
161: #define f90array2ddestroyfortranaddr_     f90array2ddestroyfortranaddr
162: #endif

164: PETSC_EXTERN void PETSC_STDCALL f90array2dcreatescalar_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array2d * PETSC_F90_2PTR_PROTO_NOVAR);
165: PETSC_EXTERN void PETSC_STDCALL f90array2daccessscalar_(F90Array2d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
166: PETSC_EXTERN void PETSC_STDCALL f90array2ddestroyscalar_(F90Array2d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
167: PETSC_EXTERN void PETSC_STDCALL f90array2dcreatereal_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array2d * PETSC_F90_2PTR_PROTO_NOVAR);
168: PETSC_EXTERN void PETSC_STDCALL f90array2daccessreal_(F90Array2d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
169: PETSC_EXTERN void PETSC_STDCALL f90array2ddestroyreal_(F90Array2d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
170: PETSC_EXTERN void PETSC_STDCALL f90array2dcreateint_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array2d * PETSC_F90_2PTR_PROTO_NOVAR);
171: PETSC_EXTERN void PETSC_STDCALL f90array2daccessint_(F90Array2d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
172: PETSC_EXTERN void PETSC_STDCALL f90array2ddestroyint_(F90Array2d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
173: PETSC_EXTERN void PETSC_STDCALL f90array2dcreatefortranaddr_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array2d * PETSC_F90_2PTR_PROTO_NOVAR);
174: PETSC_EXTERN void PETSC_STDCALL f90array2daccessfortranaddr_(F90Array2d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
175: PETSC_EXTERN void PETSC_STDCALL f90array2ddestroyfortranaddr_(F90Array2d *ptr PETSC_F90_2PTR_PROTO_NOVAR);

177: PetscErrorCode F90Array2dCreate(void *array,MPI_Datatype type,PetscInt start1,PetscInt len1,PetscInt start2,PetscInt len2,F90Array2d *ptr PETSC_F90_2PTR_PROTO(ptrd))
178: {
180:   if (type == MPIU_SCALAR) {
181:     f90array2dcreatescalar_(array,&start1,&len1,&start2,&len2,ptr PETSC_F90_2PTR_PARAM(ptrd));
182:   } else if (type == MPIU_REAL) {
183:     f90array2dcreatereal_(array,&start1,&len1,&start2,&len2,ptr PETSC_F90_2PTR_PARAM(ptrd));
184:   } else if (type == MPIU_INT) {
185:     f90array2dcreateint_(array,&start1,&len1,&start2,&len2,ptr PETSC_F90_2PTR_PARAM(ptrd));
186:   } else if (type == MPIU_FORTRANADDR) {
187:     f90array2dcreatefortranaddr_(array,&start1,&len1,&start2,&len2,ptr PETSC_F90_2PTR_PARAM(ptrd));
188:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
189:   return(0);
190: }

192: PetscErrorCode  F90Array2dAccess(F90Array2d *ptr,MPI_Datatype type,void **array PETSC_F90_2PTR_PROTO(ptrd))
193: {
195:   if (type == MPIU_SCALAR) {
196:     f90array2daccessscalar_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
197:   } else if (type == MPIU_REAL) {
198:     f90array2daccessreal_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
199:   } else if (type == MPIU_INT) {
200:     f90array2daccessint_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
201:   } else if (type == MPIU_FORTRANADDR) {
202:     f90array2daccessfortranaddr_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
203:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
204:   return(0);
205: }

207: PetscErrorCode  F90Array2dDestroy(F90Array2d *ptr,MPI_Datatype type PETSC_F90_2PTR_PROTO(ptrd))
208: {
210:   if (type == MPIU_SCALAR) {
211:     f90array2ddestroyscalar_(ptr PETSC_F90_2PTR_PARAM(ptrd));
212:   } else if (type == MPIU_REAL) {
213:     f90array2ddestroyreal_(ptr PETSC_F90_2PTR_PARAM(ptrd));
214:   } else if (type == MPIU_INT) {
215:     f90array2ddestroyint_(ptr PETSC_F90_2PTR_PARAM(ptrd));
216:   } else if (type == MPIU_FORTRANADDR) {
217:     f90array2ddestroyfortranaddr_(ptr PETSC_F90_2PTR_PARAM(ptrd));
218:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
219:   return(0);
220: }

222: /*************************************************************************/

224: #if defined(PETSC_HAVE_FORTRAN_CAPS)
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: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
238: #define f90array3dcreatescalar_           f90array3dcreatescalar
239: #define f90array3daccessscalar_           f90array3daccessscalar
240: #define f90array3ddestroyscalar_          f90array3ddestroyscalar
241: #define f90array3dcreatereal_             f90array3dcreatereal
242: #define f90array3daccessreal_             f90array3daccessreal
243: #define f90array3ddestroyreal_            f90array3ddestroyreal
244: #define f90array3dcreateint_              f90array3dcreateint
245: #define f90array3daccessint_              f90array3daccessint
246: #define f90array3ddestroyint_             f90array3ddestroyint
247: #define f90array3dcreatefortranaddr_      f90array3dcreatefortranaddr
248: #define f90array3daccessfortranaddr_      f90array3daccessfortranaddr
249: #define f90array3ddestroyfortranaddr_     f90array3ddestroyfortranaddr
250: #endif

252: PETSC_EXTERN void PETSC_STDCALL f90array3dcreatescalar_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array3d * PETSC_F90_2PTR_PROTO_NOVAR);
253: PETSC_EXTERN void PETSC_STDCALL f90array3daccessscalar_(F90Array3d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
254: PETSC_EXTERN void PETSC_STDCALL f90array3ddestroyscalar_(F90Array3d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
255: PETSC_EXTERN void PETSC_STDCALL f90array3dcreatereal_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array3d * PETSC_F90_2PTR_PROTO_NOVAR);
256: PETSC_EXTERN void PETSC_STDCALL f90array3daccessreal_(F90Array3d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
257: PETSC_EXTERN void PETSC_STDCALL f90array3ddestroyreal_(F90Array3d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
258: PETSC_EXTERN void PETSC_STDCALL f90array3dcreateint_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array3d * PETSC_F90_2PTR_PROTO_NOVAR);
259: PETSC_EXTERN void PETSC_STDCALL f90array3daccessint_(F90Array3d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
260: PETSC_EXTERN void PETSC_STDCALL f90array3ddestroyint_(F90Array3d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
261: PETSC_EXTERN void PETSC_STDCALL f90array3dcreatefortranaddr_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array3d * PETSC_F90_2PTR_PROTO_NOVAR);
262: PETSC_EXTERN void PETSC_STDCALL f90array3daccessfortranaddr_(F90Array3d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
263: PETSC_EXTERN void PETSC_STDCALL f90array3ddestroyfortranaddr_(F90Array3d *ptr PETSC_F90_2PTR_PROTO_NOVAR);

265: PetscErrorCode F90Array3dCreate(void *array,MPI_Datatype type,PetscInt start1,PetscInt len1,PetscInt start2,PetscInt len2,PetscInt start3,PetscInt len3,F90Array3d *ptr PETSC_F90_2PTR_PROTO(ptrd))
266: {
268:   if (type == MPIU_SCALAR) {
269:     f90array3dcreatescalar_(array,&start1,&len1,&start2,&len2,&start3,&len3,ptr PETSC_F90_2PTR_PARAM(ptrd));
270:   } else if (type == MPIU_REAL) {
271:     f90array3dcreatereal_(array,&start1,&len1,&start2,&len2,&start3,&len3,ptr PETSC_F90_2PTR_PARAM(ptrd));
272:   } else if (type == MPIU_INT) {
273:     f90array3dcreateint_(array,&start1,&len1,&start2,&len2,&start3,&len3,ptr PETSC_F90_2PTR_PARAM(ptrd));
274:   } else if (type == MPIU_FORTRANADDR) {
275:     f90array3dcreatefortranaddr_(array,&start1,&len1,&start2,&len2,&start3,&len3,ptr PETSC_F90_2PTR_PARAM(ptrd));
276:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
277:   return(0);
278: }

280: PetscErrorCode  F90Array3dAccess(F90Array3d *ptr,MPI_Datatype type,void **array PETSC_F90_2PTR_PROTO(ptrd))
281: {
283:   if (type == MPIU_SCALAR) {
284:     f90array3daccessscalar_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
285:   } else if (type == MPIU_REAL) {
286:     f90array3daccessreal_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
287:   } else if (type == MPIU_INT) {
288:     f90array3daccessint_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
289:   } else if (type == MPIU_FORTRANADDR) {
290:     f90array3daccessfortranaddr_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
291:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
292:   return(0);
293: }

295: PetscErrorCode  F90Array3dDestroy(F90Array3d *ptr,MPI_Datatype type PETSC_F90_2PTR_PROTO(ptrd))
296: {
298:   if (type == MPIU_SCALAR) {
299:     f90array3ddestroyscalar_(ptr PETSC_F90_2PTR_PARAM(ptrd));
300:   } else if (type == MPIU_REAL) {
301:     f90array3ddestroyreal_(ptr PETSC_F90_2PTR_PARAM(ptrd));
302:   } else if (type == MPIU_INT) {
303:     f90array3ddestroyint_(ptr PETSC_F90_2PTR_PARAM(ptrd));
304:   } else if (type == MPIU_FORTRANADDR) {
305:     f90array3ddestroyfortranaddr_(ptr PETSC_F90_2PTR_PARAM(ptrd));
306:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
307:   return(0);
308: }

310: /*************************************************************************/
311: #if defined(PETSC_HAVE_FORTRAN_CAPS)
312: #define f90array4dcreatescalar_           F90ARRAY4DCREATESCALAR
313: #define f90array4daccessscalar_           F90ARRAY4DACCESSSCALAR
314: #define f90array4ddestroyscalar_          F90ARRAY4DDESTROYSCALAR
315: #define f90array4dcreatereal_             F90ARRAY4DCREATEREAL
316: #define f90array4daccessreal_             F90ARRAY4DACCESSREAL
317: #define f90array4ddestroyreal_            F90ARRAY4DDESTROYREAL
318: #define f90array4dcreateint_              F90ARRAY4DCREATEINT
319: #define f90array4daccessint_              F90ARRAY4DACCESSINT
320: #define f90array4ddestroyint_             F90ARRAY4DDESTROYINT
321: #define f90array4dcreatefortranaddr_      F90ARRAY4DCREATEFORTRANADDR
322: #define f90array4daccessfortranaddr_      F90ARRAY4DACCESSFORTRANADDR
323: #define f90array4ddestroyfortranaddr_     F90ARRAY4DDESTROYFORTRANADDR
324: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
325: #define f90array4dcreatescalar_           f90array4dcreatescalar
326: #define f90array4daccessscalar_           f90array4daccessscalar
327: #define f90array4ddestroyscalar_          f90array4ddestroyscalar
328: #define f90array4dcreatereal_             f90array4dcreatereal
329: #define f90array4daccessreal_             f90array4daccessreal
330: #define f90array4ddestroyreal_            f90array4ddestroyreal
331: #define f90array4dcreateint_              f90array4dcreateint
332: #define f90array4daccessint_              f90array4daccessint
333: #define f90array4ddestroyint_             f90array4ddestroyint
334: #define f90array4dcreatefortranaddr_      f90array4dcreatefortranaddr
335: #define f90array4daccessfortranaddr_      f90array4daccessfortranaddr
336: #define f90array4ddestroyfortranaddr_     f90array4ddestroyfortranaddr
337: #endif

339: PETSC_EXTERN void PETSC_STDCALL f90array4dcreatescalar_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt*,PetscInt*,F90Array4d * PETSC_F90_2PTR_PROTO_NOVAR);
340: PETSC_EXTERN void PETSC_STDCALL f90array4daccessscalar_(F90Array4d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
341: PETSC_EXTERN void PETSC_STDCALL f90array4ddestroyscalar_(F90Array4d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
342: PETSC_EXTERN void PETSC_STDCALL f90array4dcreatereal_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt*,PetscInt*,F90Array4d * PETSC_F90_2PTR_PROTO_NOVAR);
343: PETSC_EXTERN void PETSC_STDCALL f90array4daccessreal_(F90Array4d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
344: PETSC_EXTERN void PETSC_STDCALL f90array4ddestroyreal_(F90Array4d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
345: PETSC_EXTERN void PETSC_STDCALL f90array4dcreateint_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt*,PetscInt*,F90Array4d * PETSC_F90_2PTR_PROTO_NOVAR);
346: PETSC_EXTERN void PETSC_STDCALL f90array4daccessint_(F90Array4d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
347: PETSC_EXTERN void PETSC_STDCALL f90array4ddestroyint_(F90Array4d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
348: PETSC_EXTERN void PETSC_STDCALL f90array4dcreatefortranaddr_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt*,PetscInt*,F90Array4d * PETSC_F90_2PTR_PROTO_NOVAR);
349: PETSC_EXTERN void PETSC_STDCALL f90array4daccessfortranaddr_(F90Array4d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
350: PETSC_EXTERN void PETSC_STDCALL f90array4ddestroyfortranaddr_(F90Array4d *ptr PETSC_F90_2PTR_PROTO_NOVAR);

352: PetscErrorCode F90Array4dCreate(void *array,MPI_Datatype type,PetscInt start1,PetscInt len1,PetscInt start2,PetscInt len2,PetscInt start3,PetscInt len3,PetscInt start4,PetscInt len4,F90Array4d *ptr PETSC_F90_2PTR_PROTO(ptrd))
353: {
355:   if (type == MPIU_SCALAR) {
356:     f90array4dcreatescalar_(array,&start1,&len1,&start2,&len2,&start3,&len3,&start4,&len4,ptr PETSC_F90_2PTR_PARAM(ptrd));
357:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
358:   return(0);
359: }

361: PetscErrorCode  F90Array4dAccess(F90Array4d *ptr,MPI_Datatype type,void **array PETSC_F90_2PTR_PROTO(ptrd))
362: {
364:   if (type == MPIU_SCALAR) {
365:     f90array4daccessscalar_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
366:   } else if (type == MPIU_REAL) {
367:     f90array4daccessreal_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
368:   } else if (type == MPIU_INT) {
369:     f90array4daccessint_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
370:   } else if (type == MPIU_FORTRANADDR) {
371:     f90array4daccessfortranaddr_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
372:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
373:   return(0);
374: }

376: PetscErrorCode  F90Array4dDestroy(F90Array4d *ptr,MPI_Datatype type PETSC_F90_2PTR_PROTO(ptrd))
377: {
379:   if (type == MPIU_SCALAR) {
380:     f90array4ddestroyscalar_(ptr PETSC_F90_2PTR_PARAM(ptrd));
381:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
382:   return(0);
383: }

385: /*************************************************************************/
386: #if defined(PETSC_HAVE_FORTRAN_CAPS)
387: #define f90array1dgetaddrscalar_            F90ARRAY1DGETADDRSCALAR
388: #define f90array1dgetaddrreal_              F90ARRAY1DGETADDRREAL
389: #define f90array1dgetaddrint_               F90ARRAY1DGETADDRINT
390: #define f90array1dgetaddrfortranaddr_       F90ARRAY1DGETADDRFORTRANADDR
391: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
392: #define f90array1dgetaddrscalar_            f90array1dgetaddrscalar
393: #define f90array1dgetaddrreal_              f90array1dgetaddrreal
394: #define f90array1dgetaddrint_               f90array1dgetaddrint
395: #define f90array1dgetaddrfortranaddr_       f90array1dgetaddrfortranaddr
396: #endif

398: PETSC_EXTERN void PETSC_STDCALL f90array1dgetaddrscalar_(void *array, PetscFortranAddr *address)
399: {
400:   *address = (PetscFortranAddr)array;
401: }
402: PETSC_EXTERN void PETSC_STDCALL f90array1dgetaddrreal_(void *array, PetscFortranAddr *address)
403: {
404:   *address = (PetscFortranAddr)array;
405: }
406: PETSC_EXTERN void PETSC_STDCALL f90array1dgetaddrint_(void *array, PetscFortranAddr *address)
407: {
408:   *address = (PetscFortranAddr)array;
409: }
410: PETSC_EXTERN void PETSC_STDCALL f90array1dgetaddrfortranaddr_(void *array, PetscFortranAddr *address)
411: {
412:   *address = (PetscFortranAddr)array;
413: }

415: /*************************************************************************/
416: #if defined(PETSC_HAVE_FORTRAN_CAPS)
417: #define f90array2dgetaddrscalar_            F90ARRAY2DGETADDRSCALAR
418: #define f90array2dgetaddrreal_              F90ARRAY2DGETADDRREAL
419: #define f90array2dgetaddrint_               F90ARRAY2DGETADDRINT
420: #define f90array2dgetaddrfortranaddr_       F90ARRAY2DGETADDRFORTRANADDR
421: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
422: #define f90array2dgetaddrscalar_            f90array2dgetaddrscalar
423: #define f90array2dgetaddrreal_              f90array2dgetaddrreal
424: #define f90array2dgetaddrint_               f90array2dgetaddrint
425: #define f90array2dgetaddrfortranaddr_       f90array2dgetaddrfortranaddr
426: #endif

428: PETSC_EXTERN void PETSC_STDCALL f90array2dgetaddrscalar_(void *array, PetscFortranAddr *address)
429: {
430:   *address = (PetscFortranAddr)array;
431: }
432: PETSC_EXTERN void PETSC_STDCALL f90array2dgetaddrreal_(void *array, PetscFortranAddr *address)
433: {
434:   *address = (PetscFortranAddr)array;
435: }
436: PETSC_EXTERN void PETSC_STDCALL f90array2dgetaddrint_(void *array, PetscFortranAddr *address)
437: {
438:   *address = (PetscFortranAddr)array;
439: }
440: PETSC_EXTERN void PETSC_STDCALL f90array2dgetaddrfortranaddr_(void *array, PetscFortranAddr *address)
441: {
442:   *address = (PetscFortranAddr)array;
443: }

445: /*************************************************************************/
446: #if defined(PETSC_HAVE_FORTRAN_CAPS)
447: #define f90array3dgetaddrscalar_            F90ARRAY3DGETADDRSCALAR
448: #define f90array3dgetaddrreal_              F90ARRAY3DGETADDRREAL
449: #define f90array3dgetaddrint_               F90ARRAY3DGETADDRINT
450: #define f90array3dgetaddrfortranaddr_       F90ARRAY3DGETADDRFORTRANADDR
451: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
452: #define f90array3dgetaddrscalar_            f90array3dgetaddrscalar
453: #define f90array3dgetaddrreal_              f90array3dgetaddrreal
454: #define f90array3dgetaddrint_               f90array3dgetaddrint
455: #define f90array3dgetaddrfortranaddr_       f90array3dgetaddrfortranaddr
456: #endif

458: PETSC_EXTERN void PETSC_STDCALL f90array3dgetaddrscalar_(void *array, PetscFortranAddr *address)
459: {
460:   *address = (PetscFortranAddr)array;
461: }
462: PETSC_EXTERN void PETSC_STDCALL f90array3dgetaddrreal_(void *array, PetscFortranAddr *address)
463: {
464:   *address = (PetscFortranAddr)array;
465: }
466: PETSC_EXTERN void PETSC_STDCALL f90array3dgetaddrint_(void *array, PetscFortranAddr *address)
467: {
468:   *address = (PetscFortranAddr)array;
469: }
470: PETSC_EXTERN void PETSC_STDCALL f90array3dgetaddrfortranaddr_(void *array, PetscFortranAddr *address)
471: {
472:   *address = (PetscFortranAddr)array;
473: }

475: /*************************************************************************/
476: #if defined(PETSC_HAVE_FORTRAN_CAPS)
477: #define f90array4dgetaddrscalar_            F90ARRAY4DGETADDRSCALAR
478: #define f90array4dgetaddrreal_              F90ARRAY4DGETADDRREAL
479: #define f90array4dgetaddrint_               F90ARRAY4DGETADDRINT
480: #define f90array4dgetaddrfortranaddr_       F90ARRAY4DGETADDRFORTRANADDR
481: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
482: #define f90array4dgetaddrscalar_            f90array4dgetaddrscalar
483: #define f90array4dgetaddrreal_              f90array4dgetaddrreal
484: #define f90array4dgetaddrint_               f90array4dgetaddrint
485: #define f90array4dgetaddrfortranaddr_       f90array4dgetaddrfortranaddr
486: #endif

488: PETSC_EXTERN void PETSC_STDCALL f90array4dgetaddrscalar_(void *array, PetscFortranAddr *address)
489: {
490:   *address = (PetscFortranAddr)array;
491: }
492: PETSC_EXTERN void PETSC_STDCALL f90array4dgetaddrreal_(void *array, PetscFortranAddr *address)
493: {
494:   *address = (PetscFortranAddr)array;
495: }
496: PETSC_EXTERN void PETSC_STDCALL f90array4dgetaddrint_(void *array, PetscFortranAddr *address)
497: {
498:   *address = (PetscFortranAddr)array;
499: }
500: PETSC_EXTERN void PETSC_STDCALL f90array4dgetaddrfortranaddr_(void *array, PetscFortranAddr *address)
501: {
502:   *address = (PetscFortranAddr)array;
503: }