Actual source code: zutils.c

  1: #include <petsc/private/fortranimpl.h>

  3: /*MC
  4:    PetscFortranAddr - a variable type in Fortran that can hold a
  5:      regular C pointer.

  7:    Note:
  8:     Used, for example, as the file argument in `PetscFOpen()`

 10:    Level: beginner

 12: .seealso:  `PetscOffset`, `PetscInt`
 13: M*/
 14: /*MC
 15:    PetscOffset - a variable type in Fortran used with `VecGetArray()`
 16:      and `ISGetIndices()`

 18:    Level: beginner

 20: .seealso:  `PetscFortranAddr`, `PetscInt`
 21: M*/

 23: /*
 24:     This is code for translating PETSc memory addresses to integer offsets
 25:     for Fortran.
 26: */
 27: char *PETSC_NULL_CHARACTER_Fortran = NULL;
 28: void *PETSC_NULL_INTEGER_Fortran   = NULL;
 29: void *PETSC_NULL_SCALAR_Fortran    = NULL;
 30: void *PETSC_NULL_DOUBLE_Fortran    = NULL;
 31: void *PETSC_NULL_REAL_Fortran      = NULL;
 32: void *PETSC_NULL_BOOL_Fortran      = NULL;
 33: EXTERN_C_BEGIN
 34: void (*PETSC_NULL_FUNCTION_Fortran)(void) = NULL;
 35: EXTERN_C_END
 36: void *PETSC_NULL_MPI_COMM_Fortran = NULL;

 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;

 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

 62:   if (base + itmp2 != addr) {
 63:     PetscCallAbort(PETSC_COMM_SELF, (*PetscErrorPrintf)("PetscIntAddressToFortran:C and Fortran arrays are\n"));
 64:     PetscCallAbort(PETSC_COMM_SELF, (*PetscErrorPrintf)("not commonly aligned or are too far apart to be indexed \n"));
 65:     PetscCallAbort(PETSC_COMM_SELF, (*PetscErrorPrintf)("by an integer. Locations: C %zu Fortran %zu\n", tmp1, tmp3));
 66:     PETSCABORT(PETSC_COMM_WORLD, PETSC_ERR_PLIB);
 67:   }
 68:   return itmp2;
 69: }

 71: PetscInt *PetscIntAddressFromFortran(const PetscInt *base, size_t addr)
 72: {
 73:   return (PetscInt *)(base + addr);
 74: }

 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

 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;

 93:   PetscFunctionBegin;
 94: #if !defined(PETSC_HAVE_CRAY90_POINTER)
 95:   if (tmp3 > tmp1) { /* C is bigger than Fortran */
 96:     tmp2  = (tmp3 - tmp1) / sizeof(PetscScalar);
 97:     itmp2 = (size_t)tmp2;
 98:     shift = (align * sizeof(PetscScalar) - (PetscInt)((tmp3 - tmp1) % (align * sizeof(PetscScalar)))) % (align * sizeof(PetscScalar));
 99:   } else {
100:     tmp2  = (tmp1 - tmp3) / sizeof(PetscScalar);
101:     itmp2 = -((size_t)tmp2);
102:     shift = (PetscInt)((tmp1 - tmp3) % (align * sizeof(PetscScalar)));
103:   }
104: #else
105:   if (tmp3 > tmp1) { /* C is bigger than Fortran */
106:     tmp2  = (tmp3 - tmp1);
107:     itmp2 = (size_t)tmp2;
108:   } else {
109:     tmp2  = (tmp1 - tmp3);
110:     itmp2 = -((size_t)tmp2);
111:   }
112:   shift = 0;
113: #endif

115:   if (shift) {
116:     /*
117:         Fortran and C not PetscScalar aligned,recover by copying values into
118:         memory that is aligned with the Fortran
119:     */
120:     PetscScalar   *work;
121:     PetscContainer container;

123:     PetscCall(PetscMalloc1(N + align, &work));

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:     }

133:     /* shift work by that number of bytes */
134:     work = (PetscScalar *)(((char *)work) + shift);
135:     PetscCall(PetscArraycpy(work, addr, N));

137:     /* store in the first location in addr how much you shift it */
138:     ((PetscInt *)addr)[0] = shift;

140:     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container));
141:     PetscCall(PetscContainerSetPointer(container, addr));
142:     PetscCall(PetscObjectCompose(obj, "GetArrayPtr", (PetscObject)container));

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:       PetscCall((*PetscErrorPrintf)("PetscScalarAddressToFortran:C and Fortran arrays are\n"));
156:       PetscCall((*PetscErrorPrintf)("not commonly aligned.\n"));
157:       PetscCall((*PetscErrorPrintf)("Locations/sizeof(PetscScalar): C %g Fortran %g\n", (double)(((PetscReal)tmp3) / (PetscReal)sizeof(PetscScalar)), (double)(((PetscReal)tmp1) / (PetscReal)sizeof(PetscScalar))));
158:       PETSCABORT(PETSC_COMM_WORLD, PETSC_ERR_PLIB);
159:     }
160:     PetscCall(PetscInfo(obj, "Efficiency warning, copying array in XXXGetArray() due\n\
161:     to alignment differences between C and Fortran\n"));
162:   }
163:   *res = itmp2;
164:   PetscFunctionReturn(PETSC_SUCCESS);
165: }

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

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:   PetscInt       shift;
178:   PetscContainer container;
179:   PetscScalar   *tlx;

181:   PetscFunctionBegin;
182:   PetscCall(PetscObjectQuery(obj, "GetArrayPtr", (PetscObject *)&container));
183:   if (container) {
184:     PetscCall(PetscContainerGetPointer(container, (void **)lx));
185:     tlx = base + addr;

187:     shift = *(PetscInt *)*lx;
188:     PetscCall(PetscArraycpy(*lx, tlx, N));
189:     tlx = (PetscScalar *)(((char *)tlx) - shift);

191:     PetscCall(PetscFree(tlx));
192:     PetscCall(PetscContainerDestroy(&container));
193:     PetscCall(PetscObjectCompose(obj, "GetArrayPtr", NULL));
194:   } else {
195:     *lx = base + addr;
196:   }
197:   PetscFunctionReturn(PETSC_SUCCESS);
198: }

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

208: PETSC_EXTERN PetscBool petscisinfornanscalar_(PetscScalar *v)
209: {
210:   return (PetscBool)PetscIsInfOrNanScalar(*v);
211: }

213: PETSC_EXTERN PetscBool petscisinfornanreal_(PetscReal *v)
214: {
215:   return (PetscBool)PetscIsInfOrNanReal(*v);
216: }