1 #include <petsc/private/fortranimpl.h> 2 3 /*MC 4 PetscFortranAddr - a variable type in Fortran that can hold a 5 regular C pointer. 6 7 Note: 8 Used, for example, as the file argument in `PetscFOpen()` 9 10 Level: beginner 11 12 .seealso: `PetscOffset`, `PetscInt` 13 M*/ 14 /*MC 15 PetscOffset - a variable type in Fortran used with `VecGetArray()` 16 and `ISGetIndices()` 17 18 Level: beginner 19 20 .seealso: `PetscFortranAddr`, `PetscInt` 21 M*/ 22 23 /* 24 This is code for translating PETSc memory addresses to integer offsets 25 for Fortran. 26 */ 27 char *PETSC_NULL_CHARACTER_Fortran = 0; 28 void *PETSC_NULL_INTEGER_Fortran = 0; 29 void *PETSC_NULL_SCALAR_Fortran = 0; 30 void *PETSC_NULL_DOUBLE_Fortran = 0; 31 void *PETSC_NULL_REAL_Fortran = 0; 32 void *PETSC_NULL_BOOL_Fortran = 0; 33 EXTERN_C_BEGIN 34 void (*PETSC_NULL_FUNCTION_Fortran)(void) = 0; 35 EXTERN_C_END 36 void *PETSC_NULL_MPI_COMM_Fortran = 0; 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 %zu Fortran %zu\n",tmp1,tmp3); 66 PETSCABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB); 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 PetscScalar *work; 120 PetscContainer container; 121 122 PetscCall(PetscMalloc1(N+align,&work)); 123 124 /* recompute shift for newly allocated space */ 125 tmp3 = (size_t) work; 126 if (tmp3 > tmp1) { /* C is bigger than Fortran */ 127 shift = (align*sizeof(PetscScalar) - (PetscInt)((tmp3 - tmp1) % (align*sizeof(PetscScalar)))) % (align*sizeof(PetscScalar)); 128 } else { 129 shift = (PetscInt)((tmp1 - tmp3) % (align*sizeof(PetscScalar))); 130 } 131 132 /* shift work by that number of bytes */ 133 work = (PetscScalar*)(((char*)work) + shift); 134 PetscCall(PetscArraycpy(work,addr,N)); 135 136 /* store in the first location in addr how much you shift it */ 137 ((PetscInt*)addr)[0] = shift; 138 139 PetscCall(PetscContainerCreate(PETSC_COMM_SELF,&container)); 140 PetscCall(PetscContainerSetPointer(container,addr)); 141 PetscCall(PetscObjectCompose(obj,"GetArrayPtr",(PetscObject)container)); 142 143 tmp3 = (size_t) work; 144 if (tmp3 > tmp1) { /* C is bigger than Fortran */ 145 tmp2 = (tmp3 - tmp1)/sizeof(PetscScalar); 146 itmp2 = (size_t) tmp2; 147 shift = (align*sizeof(PetscScalar) - (PetscInt)((tmp3 - tmp1) % (align*sizeof(PetscScalar)))) % (align*sizeof(PetscScalar)); 148 } else { 149 tmp2 = (tmp1 - tmp3)/sizeof(PetscScalar); 150 itmp2 = -((size_t) tmp2); 151 shift = (PetscInt)((tmp1 - tmp3) % (align*sizeof(PetscScalar))); 152 } 153 if (shift) { 154 (*PetscErrorPrintf)("PetscScalarAddressToFortran:C and Fortran arrays are\n"); 155 (*PetscErrorPrintf)("not commonly aligned.\n"); 156 (*PetscErrorPrintf)("Locations/sizeof(PetscScalar): C %g Fortran %g\n",(double)(((PetscReal)tmp3)/(PetscReal)sizeof(PetscScalar)),(double)(((PetscReal)tmp1)/(PetscReal)sizeof(PetscScalar))); 157 PETSCABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB); 158 } 159 PetscCall(PetscInfo(obj,"Efficiency warning, copying array in XXXGetArray() due\n\ 160 to alignment differences between C and Fortran\n")); 161 } 162 *res = itmp2; 163 return 0; 164 } 165 166 /* 167 obj - the PETSc object where the scalar pointer came from 168 base - the Fortran array address 169 addr - the Fortran offset from base 170 N - the amount of data 171 172 lx - the array space that is to be passed to XXXXRestoreArray() 173 */ 174 PetscErrorCode PetscScalarAddressFromFortran(PetscObject obj,PetscScalar *base,size_t addr,PetscInt N,PetscScalar **lx) 175 { 176 PetscInt shift; 177 PetscContainer container; 178 PetscScalar *tlx; 179 180 PetscCall(PetscObjectQuery(obj,"GetArrayPtr",(PetscObject*)&container)); 181 if (container) { 182 PetscCall(PetscContainerGetPointer(container,(void**)lx)); 183 tlx = base + addr; 184 185 shift = *(PetscInt*)*lx; 186 PetscCall(PetscArraycpy(*lx,tlx,N)); 187 tlx = (PetscScalar*)(((char*)tlx) - shift); 188 189 PetscCall(PetscFree(tlx)); 190 PetscCall(PetscContainerDestroy(&container)); 191 PetscCall(PetscObjectCompose(obj,"GetArrayPtr",0)); 192 } else { 193 *lx = base + addr; 194 } 195 return 0; 196 } 197 198 #if defined(PETSC_HAVE_FORTRAN_CAPS) 199 #define petscisinfornanscalar_ PETSCISINFORNANSCALAR 200 #define petscisinfornanreal_ PETSCISINFORNANREAL 201 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 202 #define petscisinfornanscalar_ petscisinfornanscalar 203 #define petscisinfornanreal_ petscisinfornanreal 204 #endif 205 206 PETSC_EXTERN PetscBool petscisinfornanscalar_(PetscScalar *v) 207 { 208 return (PetscBool) PetscIsInfOrNanScalar(*v); 209 } 210 211 PETSC_EXTERN PetscBool petscisinfornanreal_(PetscReal *v) 212 { 213 return (PetscBool) PetscIsInfOrNanReal(*v); 214 } 215