xref: /petsc/src/sys/ftn-custom/zutils.c (revision 5b6bfdb9644f185dbf5e5a09b808ec241507e1e7)
1 #include <petsc/private/fortranimpl.h>
2 
3 void *PETSCNULLPOINTERADDRESS = NULL;
4 
5 /*MC
6    PetscFortranAddr - a variable type in Fortran that can hold a
7      regular C pointer.
8 
9    Notes: Used, for example, as the file argument in PetscFOpen()
10 
11    Level: beginner
12 
13 .seealso:  PetscOffset, PetscInt
14 M*/
15 /*MC
16    PetscOffset - a variable type in Fortran used with VecGetArray()
17      and ISGetIndices()
18 
19    Level: beginner
20 
21 .seealso:  PetscFortranAddr, PetscInt
22 M*/
23 
24 /*
25     This is code for translating PETSc memory addresses to integer offsets
26     for Fortran.
27 */
28 char *PETSC_NULL_CHARACTER_Fortran = 0;
29 void *PETSC_NULL_INTEGER_Fortran   = 0;
30 void *PETSC_NULL_SCALAR_Fortran    = 0;
31 void *PETSC_NULL_DOUBLE_Fortran    = 0;
32 void *PETSC_NULL_REAL_Fortran      = 0;
33 void *PETSC_NULL_BOOL_Fortran      = 0;
34 EXTERN_C_BEGIN
35 void (*PETSC_NULL_FUNCTION_Fortran)(void) = 0;
36 EXTERN_C_END
37 
38 size_t PetscIntAddressToFortran(const PetscInt *base,const PetscInt *addr)
39 {
40   size_t tmp1 = (size_t) base,tmp2 = 0;
41   size_t tmp3 = (size_t) addr;
42   size_t itmp2;
43 
44 #if !defined(PETSC_HAVE_CRAY90_POINTER)
45   if (tmp3 > tmp1) {
46     tmp2  = (tmp3 - tmp1)/sizeof(PetscInt);
47     itmp2 = (size_t) tmp2;
48   } else {
49     tmp2  = (tmp1 - tmp3)/sizeof(PetscInt);
50     itmp2 = -((size_t) tmp2);
51   }
52 #else
53   if (tmp3 > tmp1) {
54     tmp2  = (tmp3 - tmp1);
55     itmp2 = (size_t) tmp2;
56   } else {
57     tmp2  = (tmp1 - tmp3);
58     itmp2 = -((size_t) tmp2);
59   }
60 #endif
61 
62   if (base + itmp2 != addr) {
63     (*PetscErrorPrintf)("PetscIntAddressToFortran:C and Fortran arrays are\n");
64     (*PetscErrorPrintf)("not commonly aligned or are too far apart to be indexed \n");
65     (*PetscErrorPrintf)("by an integer. Locations: C %uld Fortran %uld\n",tmp1,tmp3);
66     MPI_Abort(PETSC_COMM_WORLD,1);
67   }
68   return itmp2;
69 }
70 
71 PetscInt *PetscIntAddressFromFortran(const PetscInt *base,size_t addr)
72 {
73   return (PetscInt *)(base + addr);
74 }
75 
76 /*
77        obj - PETSc object on which request is made
78        base - Fortran array address
79        addr - C array address
80        res  - will contain offset from C to Fortran
81        shift - number of bytes that prevent base and addr from being commonly aligned
82        N - size of the array
83 
84        align indicates alignment relative to PetscScalar, 1 means aligned on PetscScalar, 2 means aligned on 2 PetscScalar etc
85 */
86 PetscErrorCode PetscScalarAddressToFortran(PetscObject obj,PetscInt align,PetscScalar *base,PetscScalar *addr,PetscInt N,size_t *res)
87 {
88   size_t   tmp1 = (size_t) base,tmp2;
89   size_t   tmp3 = (size_t) addr;
90   size_t   itmp2;
91   PetscInt shift;
92 
93 #if !defined(PETSC_HAVE_CRAY90_POINTER)
94   if (tmp3 > tmp1) {  /* C is bigger than Fortran */
95     tmp2  = (tmp3 - tmp1)/sizeof(PetscScalar);
96     itmp2 = (size_t) tmp2;
97     shift = (align*sizeof(PetscScalar) - (PetscInt)((tmp3 - tmp1) % (align*sizeof(PetscScalar)))) % (align*sizeof(PetscScalar));
98   } else {
99     tmp2  = (tmp1 - tmp3)/sizeof(PetscScalar);
100     itmp2 = -((size_t) tmp2);
101     shift = (PetscInt)((tmp1 - tmp3) % (align*sizeof(PetscScalar)));
102   }
103 #else
104   if (tmp3 > tmp1) {  /* C is bigger than Fortran */
105     tmp2  = (tmp3 - tmp1);
106     itmp2 = (size_t) tmp2;
107   } else {
108     tmp2  = (tmp1 - tmp3);
109     itmp2 = -((size_t) tmp2);
110   }
111   shift = 0;
112 #endif
113 
114   if (shift) {
115     /*
116         Fortran and C not PetscScalar aligned,recover by copying values into
117         memory that is aligned with the Fortran
118     */
119     PetscErrorCode ierr;
120     PetscScalar    *work;
121     PetscContainer container;
122 
123     ierr = PetscMalloc1(N+align,&work);CHKERRQ(ierr);
124 
125     /* recompute shift for newly allocated space */
126     tmp3 = (size_t) work;
127     if (tmp3 > tmp1) {  /* C is bigger than Fortran */
128       shift = (align*sizeof(PetscScalar) - (PetscInt)((tmp3 - tmp1) % (align*sizeof(PetscScalar)))) % (align*sizeof(PetscScalar));
129     } else {
130       shift = (PetscInt)((tmp1 - tmp3) % (align*sizeof(PetscScalar)));
131     }
132 
133     /* shift work by that number of bytes */
134     work = (PetscScalar*)(((char*)work) + shift);
135     ierr = PetscMemcpy(work,addr,N*sizeof(PetscScalar));CHKERRQ(ierr);
136 
137     /* store in the first location in addr how much you shift it */
138     ((PetscInt*)addr)[0] = shift;
139 
140     ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
141     ierr = PetscContainerSetPointer(container,addr);CHKERRQ(ierr);
142     ierr = PetscObjectCompose(obj,"GetArrayPtr",(PetscObject)container);CHKERRQ(ierr);
143 
144     tmp3 = (size_t) work;
145     if (tmp3 > tmp1) {  /* C is bigger than Fortran */
146       tmp2  = (tmp3 - tmp1)/sizeof(PetscScalar);
147       itmp2 = (size_t) tmp2;
148       shift = (align*sizeof(PetscScalar) - (PetscInt)((tmp3 - tmp1) % (align*sizeof(PetscScalar)))) % (align*sizeof(PetscScalar));
149     } else {
150       tmp2  = (tmp1 - tmp3)/sizeof(PetscScalar);
151       itmp2 = -((size_t) tmp2);
152       shift = (PetscInt)((tmp1 - tmp3) % (align*sizeof(PetscScalar)));
153     }
154     if (shift) {
155       (*PetscErrorPrintf)("PetscScalarAddressToFortran:C and Fortran arrays are\n");
156       (*PetscErrorPrintf)("not commonly aligned.\n");
157       (*PetscErrorPrintf)("Locations/sizeof(PetscScalar): C %f Fortran %f\n",((PetscReal)tmp3)/(PetscReal)sizeof(PetscScalar),((PetscReal)tmp1)/(PetscReal)sizeof(PetscScalar));
158       MPI_Abort(PETSC_COMM_WORLD,1);
159     }
160     ierr = PetscInfo(obj,"Efficiency warning, copying array in XXXGetArray() due\n\
161     to alignment differences between C and Fortran\n");CHKERRQ(ierr);
162   }
163   *res = itmp2;
164   return 0;
165 }
166 
167 /*
168     obj - the PETSc object where the scalar pointer came from
169     base - the Fortran array address
170     addr - the Fortran offset from base
171     N    - the amount of data
172 
173     lx   - the array space that is to be passed to XXXXRestoreArray()
174 */
175 PetscErrorCode PetscScalarAddressFromFortran(PetscObject obj,PetscScalar *base,size_t addr,PetscInt N,PetscScalar **lx)
176 {
177   PetscErrorCode ierr;
178   PetscInt       shift;
179   PetscContainer container;
180   PetscScalar    *tlx;
181 
182   ierr = PetscObjectQuery(obj,"GetArrayPtr",(PetscObject*)&container);CHKERRQ(ierr);
183   if (container) {
184     ierr = PetscContainerGetPointer(container,(void**)lx);CHKERRQ(ierr);
185     tlx  = base + addr;
186 
187     shift = *(PetscInt*)*lx;
188     ierr  = PetscMemcpy(*lx,tlx,N*sizeof(PetscScalar));CHKERRQ(ierr);
189     tlx   = (PetscScalar*)(((char*)tlx) - shift);
190 
191     ierr = PetscFree(tlx);CHKERRQ(ierr);
192     ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
193     ierr = PetscObjectCompose(obj,"GetArrayPtr",0);CHKERRQ(ierr);
194   } else {
195     *lx = base + addr;
196   }
197   return 0;
198 }
199 
200 #if defined(PETSC_HAVE_FORTRAN_CAPS)
201 #define petscisinfornanscalar_          PETSCISINFORNANSCALAR
202 #define petscisinfornanreal_            PETSCISINFORNANREAL
203 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
204 #define petscisinfornanscalar_          petscisinfornanscalar
205 #define petscisinfornanreal_            petscisinfornanreal
206 #endif
207 
208 PETSC_EXTERN PetscBool PETSC_STDCALL petscisinfornanscalar_(PetscScalar *v)
209 {
210   return (PetscBool) PetscIsInfOrNanScalar(*v);
211 }
212 
213 PETSC_EXTERN PetscBool PETSC_STDCALL petscisinfornanreal_(PetscReal *v)
214 {
215   return (PetscBool) PetscIsInfOrNanReal(*v);
216 }
217 
218 
219 
220