Actual source code: f90_cwrap.c
petsc-3.9.0 2018-04-07
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: The MPI documentation in multiple places says that one can never us
18: Fortran MPI_Datatypes in C (or vis-versa) but this is problematic since users could never
19: call C routines from Fortran that have MPI_Datatype arguments. Jed states that the Fortran
20: MPI_Datatypes will always be available in C if the MPI was built to support Fortran. This function
21: relys on this.
22: @*/
23: PetscErrorCode PetscMPIFortranDatatypeToC(MPI_Fint unit,MPI_Datatype *dtype)
24: {
25: MPI_Datatype ftype;
28: ftype = MPI_Type_f2c(unit);
29: if (ftype == MPI_INTEGER) *dtype = MPI_INT;
30: else if (ftype == MPI_INTEGER8) *dtype = MPIU_INT64;
31: else if (ftype == MPI_DOUBLE_PRECISION) *dtype = MPI_DOUBLE;
32: else if (ftype == MPI_COMPLEX16) *dtype = MPI_C_DOUBLE_COMPLEX;
33: else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unknown Fortran MPI_Datatype");
34: return(0);
35: }
37: /*************************************************************************/
39: #if defined(PETSC_HAVE_FORTRAN_CAPS)
40: #define f90array1dcreatescalar_ F90ARRAY1DCREATESCALAR
41: #define f90array1daccessscalar_ F90ARRAY1DACCESSSCALAR
42: #define f90array1ddestroyscalar_ F90ARRAY1DDESTROYSCALAR
43: #define f90array1dcreatereal_ F90ARRAY1DCREATEREAL
44: #define f90array1daccessreal_ F90ARRAY1DACCESSREAL
45: #define f90array1ddestroyreal_ F90ARRAY1DDESTROYREAL
46: #define f90array1dcreateint_ F90ARRAY1DCREATEINT
47: #define f90array1daccessint_ F90ARRAY1DACCESSINT
48: #define f90array1ddestroyint_ F90ARRAY1DDESTROYINT
49: #define f90array1dcreatefortranaddr_ F90ARRAY1DCREATEFORTRANADDR
50: #define f90array1daccessfortranaddr_ F90ARRAY1DACCESSFORTRANADDR
51: #define f90array1ddestroyfortranaddr_ F90ARRAY1DDESTROYFORTRANADDR
52: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
53: #define f90array1dcreatescalar_ f90array1dcreatescalar
54: #define f90array1daccessscalar_ f90array1daccessscalar
55: #define f90array1ddestroyscalar_ f90array1ddestroyscalar
56: #define f90array1dcreatereal_ f90array1dcreatereal
57: #define f90array1daccessreal_ f90array1daccessreal
58: #define f90array1ddestroyreal_ f90array1ddestroyreal
59: #define f90array1dcreateint_ f90array1dcreateint
60: #define f90array1daccessint_ f90array1daccessint
61: #define f90array1ddestroyint_ f90array1ddestroyint
62: #define f90array1dcreatefortranaddr_ f90array1dcreatefortranaddr
63: #define f90array1daccessfortranaddr_ f90array1daccessfortranaddr
64: #define f90array1ddestroyfortranaddr_ f90array1ddestroyfortranaddr
65: #endif
67: PETSC_EXTERN void PETSC_STDCALL f90array1dcreatescalar_(void *,PetscInt *,PetscInt *,F90Array1d * PETSC_F90_2PTR_PROTO_NOVAR);
68: PETSC_EXTERN void PETSC_STDCALL f90array1daccessscalar_(F90Array1d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
69: PETSC_EXTERN void PETSC_STDCALL f90array1ddestroyscalar_(F90Array1d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
70: PETSC_EXTERN void PETSC_STDCALL f90array1dcreatereal_(void *,PetscInt *,PetscInt *,F90Array1d * PETSC_F90_2PTR_PROTO_NOVAR);
71: PETSC_EXTERN void PETSC_STDCALL f90array1daccessreal_(F90Array1d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
72: PETSC_EXTERN void PETSC_STDCALL f90array1ddestroyreal_(F90Array1d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
73: PETSC_EXTERN void PETSC_STDCALL f90array1dcreateint_(void *,PetscInt *,PetscInt *,F90Array1d * PETSC_F90_2PTR_PROTO_NOVAR);
74: PETSC_EXTERN void PETSC_STDCALL f90array1daccessint_(F90Array1d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
75: PETSC_EXTERN void PETSC_STDCALL f90array1ddestroyint_(F90Array1d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
76: PETSC_EXTERN void PETSC_STDCALL f90array1dcreatefortranaddr_(void *,PetscInt *,PetscInt *,F90Array1d * PETSC_F90_2PTR_PROTO_NOVAR);
77: PETSC_EXTERN void PETSC_STDCALL f90array1daccessfortranaddr_(F90Array1d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
78: PETSC_EXTERN void PETSC_STDCALL f90array1ddestroyfortranaddr_(F90Array1d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
80: PetscErrorCode F90Array1dCreate(void *array,MPI_Datatype type,PetscInt start,PetscInt len,F90Array1d *ptr PETSC_F90_2PTR_PROTO(ptrd))
81: {
83: if (type == MPIU_SCALAR) {
84: if (!len) array = PETSC_NULL_SCALAR_Fortran;
85: f90array1dcreatescalar_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
86: } else if (type == MPIU_REAL) {
87: if (!len) array = PETSC_NULL_REAL_Fortran;
88: f90array1dcreatereal_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
89: } else if (type == MPIU_INT) {
90: if (!len) array = PETSC_NULL_INTEGER_Fortran;
91: f90array1dcreateint_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
92: } else if (type == MPIU_FORTRANADDR) {
93: f90array1dcreatefortranaddr_(array,&start,&len,ptr PETSC_F90_2PTR_PARAM(ptrd));
94: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
95: return(0);
96: }
98: PetscErrorCode F90Array1dAccess(F90Array1d *ptr,MPI_Datatype type,void **array PETSC_F90_2PTR_PROTO(ptrd))
99: {
101: if (type == MPIU_SCALAR) {
102: f90array1daccessscalar_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
103: if (*array == PETSC_NULL_SCALAR_Fortran) *array = 0;
104: } else if (type == MPIU_REAL) {
105: f90array1daccessreal_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
106: if (*array == PETSC_NULL_REAL_Fortran) *array = 0;
107: } else if (type == MPIU_INT) {
108: f90array1daccessint_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
109: if (*array == PETSC_NULL_INTEGER_Fortran) *array = 0;
110: } else if (type == MPIU_FORTRANADDR) {
111: f90array1daccessfortranaddr_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
112: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
113: return(0);
114: }
116: PetscErrorCode F90Array1dDestroy(F90Array1d *ptr,MPI_Datatype type PETSC_F90_2PTR_PROTO(ptrd))
117: {
119: if (type == MPIU_SCALAR) {
120: f90array1ddestroyscalar_(ptr PETSC_F90_2PTR_PARAM(ptrd));
121: } else if (type == MPIU_REAL) {
122: f90array1ddestroyreal_(ptr PETSC_F90_2PTR_PARAM(ptrd));
123: } else if (type == MPIU_INT) {
124: f90array1ddestroyint_(ptr PETSC_F90_2PTR_PARAM(ptrd));
125: } else if (type == MPIU_FORTRANADDR) {
126: f90array1ddestroyfortranaddr_(ptr PETSC_F90_2PTR_PARAM(ptrd));
127: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
128: return(0);
129: }
131: /*************************************************************************/
133: #if defined(PETSC_HAVE_FORTRAN_CAPS)
134: #define f90array2dcreatescalar_ F90ARRAY2DCREATESCALAR
135: #define f90array2daccessscalar_ F90ARRAY2DACCESSSCALAR
136: #define f90array2ddestroyscalar_ F90ARRAY2DDESTROYSCALAR
137: #define f90array2dcreatereal_ F90ARRAY2DCREATEREAL
138: #define f90array2daccessreal_ F90ARRAY2DACCESSREAL
139: #define f90array2ddestroyreal_ F90ARRAY2DDESTROYREAL
140: #define f90array2dcreateint_ F90ARRAY2DCREATEINT
141: #define f90array2daccessint_ F90ARRAY2DACCESSINT
142: #define f90array2ddestroyint_ F90ARRAY2DDESTROYINT
143: #define f90array2dcreatefortranaddr_ F90ARRAY2DCREATEFORTRANADDR
144: #define f90array2daccessfortranaddr_ F90ARRAY2DACCESSFORTRANADDR
145: #define f90array2ddestroyfortranaddr_ F90ARRAY2DDESTROYFORTRANADDR
146: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
147: #define f90array2dcreatescalar_ f90array2dcreatescalar
148: #define f90array2daccessscalar_ f90array2daccessscalar
149: #define f90array2ddestroyscalar_ f90array2ddestroyscalar
150: #define f90array2dcreatereal_ f90array2dcreatereal
151: #define f90array2daccessreal_ f90array2daccessreal
152: #define f90array2ddestroyreal_ f90array2ddestroyreal
153: #define f90array2dcreateint_ f90array2dcreateint
154: #define f90array2daccessint_ f90array2daccessint
155: #define f90array2ddestroyint_ f90array2ddestroyint
156: #define f90array2dcreatefortranaddr_ f90array2dcreatefortranaddr
157: #define f90array2daccessfortranaddr_ f90array2daccessfortranaddr
158: #define f90array2ddestroyfortranaddr_ f90array2ddestroyfortranaddr
159: #endif
161: PETSC_EXTERN void PETSC_STDCALL f90array2dcreatescalar_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array2d * PETSC_F90_2PTR_PROTO_NOVAR);
162: PETSC_EXTERN void PETSC_STDCALL f90array2daccessscalar_(F90Array2d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
163: PETSC_EXTERN void PETSC_STDCALL f90array2ddestroyscalar_(F90Array2d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
164: PETSC_EXTERN void PETSC_STDCALL f90array2dcreatereal_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array2d * PETSC_F90_2PTR_PROTO_NOVAR);
165: PETSC_EXTERN void PETSC_STDCALL f90array2daccessreal_(F90Array2d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
166: PETSC_EXTERN void PETSC_STDCALL f90array2ddestroyreal_(F90Array2d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
167: PETSC_EXTERN void PETSC_STDCALL f90array2dcreateint_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array2d * PETSC_F90_2PTR_PROTO_NOVAR);
168: PETSC_EXTERN void PETSC_STDCALL f90array2daccessint_(F90Array2d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
169: PETSC_EXTERN void PETSC_STDCALL f90array2ddestroyint_(F90Array2d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
170: PETSC_EXTERN void PETSC_STDCALL f90array2dcreatefortranaddr_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array2d * PETSC_F90_2PTR_PROTO_NOVAR);
171: PETSC_EXTERN void PETSC_STDCALL f90array2daccessfortranaddr_(F90Array2d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
172: PETSC_EXTERN void PETSC_STDCALL f90array2ddestroyfortranaddr_(F90Array2d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
174: PetscErrorCode F90Array2dCreate(void *array,MPI_Datatype type,PetscInt start1,PetscInt len1,PetscInt start2,PetscInt len2,F90Array2d *ptr PETSC_F90_2PTR_PROTO(ptrd))
175: {
177: if (type == MPIU_SCALAR) {
178: f90array2dcreatescalar_(array,&start1,&len1,&start2,&len2,ptr PETSC_F90_2PTR_PARAM(ptrd));
179: } else if (type == MPIU_REAL) {
180: f90array2dcreatereal_(array,&start1,&len1,&start2,&len2,ptr PETSC_F90_2PTR_PARAM(ptrd));
181: } else if (type == MPIU_INT) {
182: f90array2dcreateint_(array,&start1,&len1,&start2,&len2,ptr PETSC_F90_2PTR_PARAM(ptrd));
183: } else if (type == MPIU_FORTRANADDR) {
184: f90array2dcreatefortranaddr_(array,&start1,&len1,&start2,&len2,ptr PETSC_F90_2PTR_PARAM(ptrd));
185: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
186: return(0);
187: }
189: PetscErrorCode F90Array2dAccess(F90Array2d *ptr,MPI_Datatype type,void **array PETSC_F90_2PTR_PROTO(ptrd))
190: {
192: if (type == MPIU_SCALAR) {
193: f90array2daccessscalar_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
194: } else if (type == MPIU_REAL) {
195: f90array2daccessreal_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
196: } else if (type == MPIU_INT) {
197: f90array2daccessint_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
198: } else if (type == MPIU_FORTRANADDR) {
199: f90array2daccessfortranaddr_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
200: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
201: return(0);
202: }
204: PetscErrorCode F90Array2dDestroy(F90Array2d *ptr,MPI_Datatype type PETSC_F90_2PTR_PROTO(ptrd))
205: {
207: if (type == MPIU_SCALAR) {
208: f90array2ddestroyscalar_(ptr PETSC_F90_2PTR_PARAM(ptrd));
209: } else if (type == MPIU_REAL) {
210: f90array2ddestroyreal_(ptr PETSC_F90_2PTR_PARAM(ptrd));
211: } else if (type == MPIU_INT) {
212: f90array2ddestroyint_(ptr PETSC_F90_2PTR_PARAM(ptrd));
213: } else if (type == MPIU_FORTRANADDR) {
214: f90array2ddestroyfortranaddr_(ptr PETSC_F90_2PTR_PARAM(ptrd));
215: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
216: return(0);
217: }
219: /*************************************************************************/
221: #if defined(PETSC_HAVE_FORTRAN_CAPS)
222: #define f90array3dcreatescalar_ F90ARRAY3DCREATESCALAR
223: #define f90array3daccessscalar_ F90ARRAY3DACCESSSCALAR
224: #define f90array3ddestroyscalar_ F90ARRAY3DDESTROYSCALAR
225: #define f90array3dcreatereal_ F90ARRAY3DCREATEREAL
226: #define f90array3daccessreal_ F90ARRAY3DACCESSREAL
227: #define f90array3ddestroyreal_ F90ARRAY3DDESTROYREAL
228: #define f90array3dcreateint_ F90ARRAY3DCREATEINT
229: #define f90array3daccessint_ F90ARRAY3DACCESSINT
230: #define f90array3ddestroyint_ F90ARRAY3DDESTROYINT
231: #define f90array3dcreatefortranaddr_ F90ARRAY3DCREATEFORTRANADDR
232: #define f90array3daccessfortranaddr_ F90ARRAY3DACCESSFORTRANADDR
233: #define f90array3ddestroyfortranaddr_ F90ARRAY3DDESTROYFORTRANADDR
234: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
235: #define f90array3dcreatescalar_ f90array3dcreatescalar
236: #define f90array3daccessscalar_ f90array3daccessscalar
237: #define f90array3ddestroyscalar_ f90array3ddestroyscalar
238: #define f90array3dcreatereal_ f90array3dcreatereal
239: #define f90array3daccessreal_ f90array3daccessreal
240: #define f90array3ddestroyreal_ f90array3ddestroyreal
241: #define f90array3dcreateint_ f90array3dcreateint
242: #define f90array3daccessint_ f90array3daccessint
243: #define f90array3ddestroyint_ f90array3ddestroyint
244: #define f90array3dcreatefortranaddr_ f90array3dcreatefortranaddr
245: #define f90array3daccessfortranaddr_ f90array3daccessfortranaddr
246: #define f90array3ddestroyfortranaddr_ f90array3ddestroyfortranaddr
247: #endif
249: PETSC_EXTERN void PETSC_STDCALL f90array3dcreatescalar_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array3d * PETSC_F90_2PTR_PROTO_NOVAR);
250: PETSC_EXTERN void PETSC_STDCALL f90array3daccessscalar_(F90Array3d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
251: PETSC_EXTERN void PETSC_STDCALL f90array3ddestroyscalar_(F90Array3d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
252: PETSC_EXTERN void PETSC_STDCALL f90array3dcreatereal_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array3d * PETSC_F90_2PTR_PROTO_NOVAR);
253: PETSC_EXTERN void PETSC_STDCALL f90array3daccessreal_(F90Array3d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
254: PETSC_EXTERN void PETSC_STDCALL f90array3ddestroyreal_(F90Array3d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
255: PETSC_EXTERN void PETSC_STDCALL f90array3dcreateint_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array3d * PETSC_F90_2PTR_PROTO_NOVAR);
256: PETSC_EXTERN void PETSC_STDCALL f90array3daccessint_(F90Array3d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
257: PETSC_EXTERN void PETSC_STDCALL f90array3ddestroyint_(F90Array3d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
258: PETSC_EXTERN void PETSC_STDCALL f90array3dcreatefortranaddr_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,F90Array3d * PETSC_F90_2PTR_PROTO_NOVAR);
259: PETSC_EXTERN void PETSC_STDCALL f90array3daccessfortranaddr_(F90Array3d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
260: PETSC_EXTERN void PETSC_STDCALL f90array3ddestroyfortranaddr_(F90Array3d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
262: 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))
263: {
265: if (type == MPIU_SCALAR) {
266: f90array3dcreatescalar_(array,&start1,&len1,&start2,&len2,&start3,&len3,ptr PETSC_F90_2PTR_PARAM(ptrd));
267: } else if (type == MPIU_REAL) {
268: f90array3dcreatereal_(array,&start1,&len1,&start2,&len2,&start3,&len3,ptr PETSC_F90_2PTR_PARAM(ptrd));
269: } else if (type == MPIU_INT) {
270: f90array3dcreateint_(array,&start1,&len1,&start2,&len2,&start3,&len3,ptr PETSC_F90_2PTR_PARAM(ptrd));
271: } else if (type == MPIU_FORTRANADDR) {
272: f90array3dcreatefortranaddr_(array,&start1,&len1,&start2,&len2,&start3,&len3,ptr PETSC_F90_2PTR_PARAM(ptrd));
273: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
274: return(0);
275: }
277: PetscErrorCode F90Array3dAccess(F90Array3d *ptr,MPI_Datatype type,void **array PETSC_F90_2PTR_PROTO(ptrd))
278: {
280: if (type == MPIU_SCALAR) {
281: f90array3daccessscalar_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
282: } else if (type == MPIU_REAL) {
283: f90array3daccessreal_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
284: } else if (type == MPIU_INT) {
285: f90array3daccessint_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
286: } else if (type == MPIU_FORTRANADDR) {
287: f90array3daccessfortranaddr_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
288: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
289: return(0);
290: }
292: PetscErrorCode F90Array3dDestroy(F90Array3d *ptr,MPI_Datatype type PETSC_F90_2PTR_PROTO(ptrd))
293: {
295: if (type == MPIU_SCALAR) {
296: f90array3ddestroyscalar_(ptr PETSC_F90_2PTR_PARAM(ptrd));
297: } else if (type == MPIU_REAL) {
298: f90array3ddestroyreal_(ptr PETSC_F90_2PTR_PARAM(ptrd));
299: } else if (type == MPIU_INT) {
300: f90array3ddestroyint_(ptr PETSC_F90_2PTR_PARAM(ptrd));
301: } else if (type == MPIU_FORTRANADDR) {
302: f90array3ddestroyfortranaddr_(ptr PETSC_F90_2PTR_PARAM(ptrd));
303: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
304: return(0);
305: }
307: /*************************************************************************/
308: #if defined(PETSC_HAVE_FORTRAN_CAPS)
309: #define f90array4dcreatescalar_ F90ARRAY4DCREATESCALAR
310: #define f90array4daccessscalar_ F90ARRAY4DACCESSSCALAR
311: #define f90array4ddestroyscalar_ F90ARRAY4DDESTROYSCALAR
312: #define f90array4dcreatereal_ F90ARRAY4DCREATEREAL
313: #define f90array4daccessreal_ F90ARRAY4DACCESSREAL
314: #define f90array4ddestroyreal_ F90ARRAY4DDESTROYREAL
315: #define f90array4dcreateint_ F90ARRAY4DCREATEINT
316: #define f90array4daccessint_ F90ARRAY4DACCESSINT
317: #define f90array4ddestroyint_ F90ARRAY4DDESTROYINT
318: #define f90array4dcreatefortranaddr_ F90ARRAY4DCREATEFORTRANADDR
319: #define f90array4daccessfortranaddr_ F90ARRAY4DACCESSFORTRANADDR
320: #define f90array4ddestroyfortranaddr_ F90ARRAY4DDESTROYFORTRANADDR
321: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
322: #define f90array4dcreatescalar_ f90array4dcreatescalar
323: #define f90array4daccessscalar_ f90array4daccessscalar
324: #define f90array4ddestroyscalar_ f90array4ddestroyscalar
325: #define f90array4dcreatereal_ f90array4dcreatereal
326: #define f90array4daccessreal_ f90array4daccessreal
327: #define f90array4ddestroyreal_ f90array4ddestroyreal
328: #define f90array4dcreateint_ f90array4dcreateint
329: #define f90array4daccessint_ f90array4daccessint
330: #define f90array4ddestroyint_ f90array4ddestroyint
331: #define f90array4dcreatefortranaddr_ f90array4dcreatefortranaddr
332: #define f90array4daccessfortranaddr_ f90array4daccessfortranaddr
333: #define f90array4ddestroyfortranaddr_ f90array4ddestroyfortranaddr
334: #endif
336: PETSC_EXTERN void PETSC_STDCALL f90array4dcreatescalar_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt*,PetscInt*,F90Array4d * PETSC_F90_2PTR_PROTO_NOVAR);
337: PETSC_EXTERN void PETSC_STDCALL f90array4daccessscalar_(F90Array4d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
338: PETSC_EXTERN void PETSC_STDCALL f90array4ddestroyscalar_(F90Array4d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
339: PETSC_EXTERN void PETSC_STDCALL f90array4dcreatereal_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt*,PetscInt*,F90Array4d * PETSC_F90_2PTR_PROTO_NOVAR);
340: PETSC_EXTERN void PETSC_STDCALL f90array4daccessreal_(F90Array4d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
341: PETSC_EXTERN void PETSC_STDCALL f90array4ddestroyreal_(F90Array4d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
342: PETSC_EXTERN void PETSC_STDCALL f90array4dcreateint_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt*,PetscInt*,F90Array4d * PETSC_F90_2PTR_PROTO_NOVAR);
343: PETSC_EXTERN void PETSC_STDCALL f90array4daccessint_(F90Array4d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
344: PETSC_EXTERN void PETSC_STDCALL f90array4ddestroyint_(F90Array4d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
345: PETSC_EXTERN void PETSC_STDCALL f90array4dcreatefortranaddr_(void *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt *,PetscInt*,PetscInt*,F90Array4d * PETSC_F90_2PTR_PROTO_NOVAR);
346: PETSC_EXTERN void PETSC_STDCALL f90array4daccessfortranaddr_(F90Array4d*,void** PETSC_F90_2PTR_PROTO_NOVAR);
347: PETSC_EXTERN void PETSC_STDCALL f90array4ddestroyfortranaddr_(F90Array4d *ptr PETSC_F90_2PTR_PROTO_NOVAR);
349: 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))
350: {
352: if (type == MPIU_SCALAR) {
353: f90array4dcreatescalar_(array,&start1,&len1,&start2,&len2,&start3,&len3,&start4,&len4,ptr PETSC_F90_2PTR_PARAM(ptrd));
354: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
355: return(0);
356: }
358: PetscErrorCode F90Array4dAccess(F90Array4d *ptr,MPI_Datatype type,void **array PETSC_F90_2PTR_PROTO(ptrd))
359: {
361: if (type == MPIU_SCALAR) {
362: f90array4daccessscalar_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
363: } else if (type == MPIU_REAL) {
364: f90array4daccessreal_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
365: } else if (type == MPIU_INT) {
366: f90array4daccessint_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
367: } else if (type == MPIU_FORTRANADDR) {
368: f90array4daccessfortranaddr_(ptr,array PETSC_F90_2PTR_PARAM(ptrd));
369: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
370: return(0);
371: }
373: PetscErrorCode F90Array4dDestroy(F90Array4d *ptr,MPI_Datatype type PETSC_F90_2PTR_PROTO(ptrd))
374: {
376: if (type == MPIU_SCALAR) {
377: f90array4ddestroyscalar_(ptr PETSC_F90_2PTR_PARAM(ptrd));
378: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported MPI_Datatype");
379: return(0);
380: }
382: /*************************************************************************/
383: #if defined(PETSC_HAVE_FORTRAN_CAPS)
384: #define f90array1dgetaddrscalar_ F90ARRAY1DGETADDRSCALAR
385: #define f90array1dgetaddrreal_ F90ARRAY1DGETADDRREAL
386: #define f90array1dgetaddrint_ F90ARRAY1DGETADDRINT
387: #define f90array1dgetaddrfortranaddr_ F90ARRAY1DGETADDRFORTRANADDR
388: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
389: #define f90array1dgetaddrscalar_ f90array1dgetaddrscalar
390: #define f90array1dgetaddrreal_ f90array1dgetaddrreal
391: #define f90array1dgetaddrint_ f90array1dgetaddrint
392: #define f90array1dgetaddrfortranaddr_ f90array1dgetaddrfortranaddr
393: #endif
395: PETSC_EXTERN void PETSC_STDCALL f90array1dgetaddrscalar_(void *array, PetscFortranAddr *address)
396: {
397: *address = (PetscFortranAddr)array;
398: }
399: PETSC_EXTERN void PETSC_STDCALL f90array1dgetaddrreal_(void *array, PetscFortranAddr *address)
400: {
401: *address = (PetscFortranAddr)array;
402: }
403: PETSC_EXTERN void PETSC_STDCALL f90array1dgetaddrint_(void *array, PetscFortranAddr *address)
404: {
405: *address = (PetscFortranAddr)array;
406: }
407: PETSC_EXTERN void PETSC_STDCALL f90array1dgetaddrfortranaddr_(void *array, PetscFortranAddr *address)
408: {
409: *address = (PetscFortranAddr)array;
410: }
412: /*************************************************************************/
413: #if defined(PETSC_HAVE_FORTRAN_CAPS)
414: #define f90array2dgetaddrscalar_ F90ARRAY2DGETADDRSCALAR
415: #define f90array2dgetaddrreal_ F90ARRAY2DGETADDRREAL
416: #define f90array2dgetaddrint_ F90ARRAY2DGETADDRINT
417: #define f90array2dgetaddrfortranaddr_ F90ARRAY2DGETADDRFORTRANADDR
418: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
419: #define f90array2dgetaddrscalar_ f90array2dgetaddrscalar
420: #define f90array2dgetaddrreal_ f90array2dgetaddrreal
421: #define f90array2dgetaddrint_ f90array2dgetaddrint
422: #define f90array2dgetaddrfortranaddr_ f90array2dgetaddrfortranaddr
423: #endif
425: PETSC_EXTERN void PETSC_STDCALL f90array2dgetaddrscalar_(void *array, PetscFortranAddr *address)
426: {
427: *address = (PetscFortranAddr)array;
428: }
429: PETSC_EXTERN void PETSC_STDCALL f90array2dgetaddrreal_(void *array, PetscFortranAddr *address)
430: {
431: *address = (PetscFortranAddr)array;
432: }
433: PETSC_EXTERN void PETSC_STDCALL f90array2dgetaddrint_(void *array, PetscFortranAddr *address)
434: {
435: *address = (PetscFortranAddr)array;
436: }
437: PETSC_EXTERN void PETSC_STDCALL f90array2dgetaddrfortranaddr_(void *array, PetscFortranAddr *address)
438: {
439: *address = (PetscFortranAddr)array;
440: }
442: /*************************************************************************/
443: #if defined(PETSC_HAVE_FORTRAN_CAPS)
444: #define f90array3dgetaddrscalar_ F90ARRAY3DGETADDRSCALAR
445: #define f90array3dgetaddrreal_ F90ARRAY3DGETADDRREAL
446: #define f90array3dgetaddrint_ F90ARRAY3DGETADDRINT
447: #define f90array3dgetaddrfortranaddr_ F90ARRAY3DGETADDRFORTRANADDR
448: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
449: #define f90array3dgetaddrscalar_ f90array3dgetaddrscalar
450: #define f90array3dgetaddrreal_ f90array3dgetaddrreal
451: #define f90array3dgetaddrint_ f90array3dgetaddrint
452: #define f90array3dgetaddrfortranaddr_ f90array3dgetaddrfortranaddr
453: #endif
455: PETSC_EXTERN void PETSC_STDCALL f90array3dgetaddrscalar_(void *array, PetscFortranAddr *address)
456: {
457: *address = (PetscFortranAddr)array;
458: }
459: PETSC_EXTERN void PETSC_STDCALL f90array3dgetaddrreal_(void *array, PetscFortranAddr *address)
460: {
461: *address = (PetscFortranAddr)array;
462: }
463: PETSC_EXTERN void PETSC_STDCALL f90array3dgetaddrint_(void *array, PetscFortranAddr *address)
464: {
465: *address = (PetscFortranAddr)array;
466: }
467: PETSC_EXTERN void PETSC_STDCALL f90array3dgetaddrfortranaddr_(void *array, PetscFortranAddr *address)
468: {
469: *address = (PetscFortranAddr)array;
470: }
472: /*************************************************************************/
473: #if defined(PETSC_HAVE_FORTRAN_CAPS)
474: #define f90array4dgetaddrscalar_ F90ARRAY4DGETADDRSCALAR
475: #define f90array4dgetaddrreal_ F90ARRAY4DGETADDRREAL
476: #define f90array4dgetaddrint_ F90ARRAY4DGETADDRINT
477: #define f90array4dgetaddrfortranaddr_ F90ARRAY4DGETADDRFORTRANADDR
478: #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
479: #define f90array4dgetaddrscalar_ f90array4dgetaddrscalar
480: #define f90array4dgetaddrreal_ f90array4dgetaddrreal
481: #define f90array4dgetaddrint_ f90array4dgetaddrint
482: #define f90array4dgetaddrfortranaddr_ f90array4dgetaddrfortranaddr
483: #endif
485: PETSC_EXTERN void PETSC_STDCALL f90array4dgetaddrscalar_(void *array, PetscFortranAddr *address)
486: {
487: *address = (PetscFortranAddr)array;
488: }
489: PETSC_EXTERN void PETSC_STDCALL f90array4dgetaddrreal_(void *array, PetscFortranAddr *address)
490: {
491: *address = (PetscFortranAddr)array;
492: }
493: PETSC_EXTERN void PETSC_STDCALL f90array4dgetaddrint_(void *array, PetscFortranAddr *address)
494: {
495: *address = (PetscFortranAddr)array;
496: }
497: PETSC_EXTERN void PETSC_STDCALL f90array4dgetaddrfortranaddr_(void *array, PetscFortranAddr *address)
498: {
499: *address = (PetscFortranAddr)array;
500: }