1af0996ceSBarry Smith #include <petsc/private/fortranimpl.h> 29306f9a3SSatish Balay 39306f9a3SSatish Balay /*MC 49306f9a3SSatish Balay PetscFortranAddr - a variable type in Fortran that can hold a 59306f9a3SSatish Balay regular C pointer. 69306f9a3SSatish Balay 795452b02SPatrick Sanan Notes: 895452b02SPatrick Sanan Used, for example, as the file argument in PetscFOpen() 99306f9a3SSatish Balay 109306f9a3SSatish Balay Level: beginner 119306f9a3SSatish Balay 129306f9a3SSatish Balay .seealso: PetscOffset, PetscInt 139306f9a3SSatish Balay M*/ 149306f9a3SSatish Balay /*MC 159306f9a3SSatish Balay PetscOffset - a variable type in Fortran used with VecGetArray() 169306f9a3SSatish Balay and ISGetIndices() 179306f9a3SSatish Balay 189306f9a3SSatish Balay Level: beginner 199306f9a3SSatish Balay 209306f9a3SSatish Balay .seealso: PetscFortranAddr, PetscInt 219306f9a3SSatish Balay M*/ 229306f9a3SSatish Balay 239306f9a3SSatish Balay /* 249306f9a3SSatish Balay This is code for translating PETSc memory addresses to integer offsets 259306f9a3SSatish Balay for Fortran. 269306f9a3SSatish Balay */ 279306f9a3SSatish Balay char *PETSC_NULL_CHARACTER_Fortran = 0; 289306f9a3SSatish Balay void *PETSC_NULL_INTEGER_Fortran = 0; 299306f9a3SSatish Balay void *PETSC_NULL_SCALAR_Fortran = 0; 309306f9a3SSatish Balay void *PETSC_NULL_DOUBLE_Fortran = 0; 319306f9a3SSatish Balay void *PETSC_NULL_REAL_Fortran = 0; 325c550465SJed Brown void *PETSC_NULL_BOOL_Fortran = 0; 339306f9a3SSatish Balay EXTERN_C_BEGIN 349306f9a3SSatish Balay void (*PETSC_NULL_FUNCTION_Fortran)(void) = 0; 359306f9a3SSatish Balay EXTERN_C_END 3699e0435eSBarry Smith 378ea3bf28SBarry Smith size_t PetscIntAddressToFortran(const PetscInt *base,const PetscInt *addr) 389306f9a3SSatish Balay { 399306f9a3SSatish Balay size_t tmp1 = (size_t) base,tmp2 = 0; 409306f9a3SSatish Balay size_t tmp3 = (size_t) addr; 419306f9a3SSatish Balay size_t itmp2; 429306f9a3SSatish Balay 439306f9a3SSatish Balay #if !defined(PETSC_HAVE_CRAY90_POINTER) 449306f9a3SSatish Balay if (tmp3 > tmp1) { 459306f9a3SSatish Balay tmp2 = (tmp3 - tmp1)/sizeof(PetscInt); 469306f9a3SSatish Balay itmp2 = (size_t) tmp2; 479306f9a3SSatish Balay } else { 489306f9a3SSatish Balay tmp2 = (tmp1 - tmp3)/sizeof(PetscInt); 499306f9a3SSatish Balay itmp2 = -((size_t) tmp2); 509306f9a3SSatish Balay } 519306f9a3SSatish Balay #else 529306f9a3SSatish Balay if (tmp3 > tmp1) { 539306f9a3SSatish Balay tmp2 = (tmp3 - tmp1); 549306f9a3SSatish Balay itmp2 = (size_t) tmp2; 559306f9a3SSatish Balay } else { 569306f9a3SSatish Balay tmp2 = (tmp1 - tmp3); 579306f9a3SSatish Balay itmp2 = -((size_t) tmp2); 589306f9a3SSatish Balay } 599306f9a3SSatish Balay #endif 609306f9a3SSatish Balay 619306f9a3SSatish Balay if (base + itmp2 != addr) { 629306f9a3SSatish Balay (*PetscErrorPrintf)("PetscIntAddressToFortran:C and Fortran arrays are\n"); 639306f9a3SSatish Balay (*PetscErrorPrintf)("not commonly aligned or are too far apart to be indexed \n"); 649306f9a3SSatish Balay (*PetscErrorPrintf)("by an integer. Locations: C %uld Fortran %uld\n",tmp1,tmp3); 65*41e02c4dSJunchao Zhang PETSCABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB); 669306f9a3SSatish Balay } 679306f9a3SSatish Balay return itmp2; 689306f9a3SSatish Balay } 699306f9a3SSatish Balay 708ea3bf28SBarry Smith PetscInt *PetscIntAddressFromFortran(const PetscInt *base,size_t addr) 719306f9a3SSatish Balay { 728ea3bf28SBarry Smith return (PetscInt *)(base + addr); 739306f9a3SSatish Balay } 749306f9a3SSatish Balay 759306f9a3SSatish Balay /* 769306f9a3SSatish Balay obj - PETSc object on which request is made 779306f9a3SSatish Balay base - Fortran array address 789306f9a3SSatish Balay addr - C array address 799306f9a3SSatish Balay res - will contain offset from C to Fortran 809306f9a3SSatish Balay shift - number of bytes that prevent base and addr from being commonly aligned 819306f9a3SSatish Balay N - size of the array 829306f9a3SSatish Balay 83f91d1997SBarry Smith align indicates alignment relative to PetscScalar, 1 means aligned on PetscScalar, 2 means aligned on 2 PetscScalar etc 849306f9a3SSatish Balay */ 85f91d1997SBarry Smith PetscErrorCode PetscScalarAddressToFortran(PetscObject obj,PetscInt align,PetscScalar *base,PetscScalar *addr,PetscInt N,size_t *res) 869306f9a3SSatish Balay { 87e366c363SBarry Smith size_t tmp1 = (size_t) base,tmp2; 889306f9a3SSatish Balay size_t tmp3 = (size_t) addr; 899306f9a3SSatish Balay size_t itmp2; 909306f9a3SSatish Balay PetscInt shift; 919306f9a3SSatish Balay 929306f9a3SSatish Balay #if !defined(PETSC_HAVE_CRAY90_POINTER) 939306f9a3SSatish Balay if (tmp3 > tmp1) { /* C is bigger than Fortran */ 949306f9a3SSatish Balay tmp2 = (tmp3 - tmp1)/sizeof(PetscScalar); 959306f9a3SSatish Balay itmp2 = (size_t) tmp2; 96f91d1997SBarry Smith shift = (align*sizeof(PetscScalar) - (PetscInt)((tmp3 - tmp1) % (align*sizeof(PetscScalar)))) % (align*sizeof(PetscScalar)); 979306f9a3SSatish Balay } else { 989306f9a3SSatish Balay tmp2 = (tmp1 - tmp3)/sizeof(PetscScalar); 999306f9a3SSatish Balay itmp2 = -((size_t) tmp2); 100f91d1997SBarry Smith shift = (PetscInt)((tmp1 - tmp3) % (align*sizeof(PetscScalar))); 1019306f9a3SSatish Balay } 1029306f9a3SSatish Balay #else 1039306f9a3SSatish Balay if (tmp3 > tmp1) { /* C is bigger than Fortran */ 1049306f9a3SSatish Balay tmp2 = (tmp3 - tmp1); 1059306f9a3SSatish Balay itmp2 = (size_t) tmp2; 1069306f9a3SSatish Balay } else { 1079306f9a3SSatish Balay tmp2 = (tmp1 - tmp3); 1089306f9a3SSatish Balay itmp2 = -((size_t) tmp2); 1099306f9a3SSatish Balay } 1109306f9a3SSatish Balay shift = 0; 1119306f9a3SSatish Balay #endif 1129306f9a3SSatish Balay 1139306f9a3SSatish Balay if (shift) { 1149306f9a3SSatish Balay /* 1159306f9a3SSatish Balay Fortran and C not PetscScalar aligned,recover by copying values into 1169306f9a3SSatish Balay memory that is aligned with the Fortran 1179306f9a3SSatish Balay */ 1189306f9a3SSatish Balay PetscErrorCode ierr; 1199306f9a3SSatish Balay PetscScalar *work; 120776b82aeSLisandro Dalcin PetscContainer container; 1219306f9a3SSatish Balay 122854ce69bSBarry Smith ierr = PetscMalloc1(N+align,&work);CHKERRQ(ierr); 123f91d1997SBarry Smith 124f91d1997SBarry Smith /* recompute shift for newly allocated space */ 125f91d1997SBarry Smith tmp3 = (size_t) work; 126f91d1997SBarry Smith if (tmp3 > tmp1) { /* C is bigger than Fortran */ 127f91d1997SBarry Smith shift = (align*sizeof(PetscScalar) - (PetscInt)((tmp3 - tmp1) % (align*sizeof(PetscScalar)))) % (align*sizeof(PetscScalar)); 128f91d1997SBarry Smith } else { 129f91d1997SBarry Smith shift = (PetscInt)((tmp1 - tmp3) % (align*sizeof(PetscScalar))); 130f91d1997SBarry Smith } 1319306f9a3SSatish Balay 1329306f9a3SSatish Balay /* shift work by that number of bytes */ 1339306f9a3SSatish Balay work = (PetscScalar*)(((char*)work) + shift); 134580bdb30SBarry Smith ierr = PetscArraycpy(work,addr,N);CHKERRQ(ierr); 1359306f9a3SSatish Balay 1369306f9a3SSatish Balay /* store in the first location in addr how much you shift it */ 1379306f9a3SSatish Balay ((PetscInt*)addr)[0] = shift; 1389306f9a3SSatish Balay 139776b82aeSLisandro Dalcin ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 140776b82aeSLisandro Dalcin ierr = PetscContainerSetPointer(container,addr);CHKERRQ(ierr); 1419306f9a3SSatish Balay ierr = PetscObjectCompose(obj,"GetArrayPtr",(PetscObject)container);CHKERRQ(ierr); 1429306f9a3SSatish Balay 1439306f9a3SSatish Balay tmp3 = (size_t) work; 1449306f9a3SSatish Balay if (tmp3 > tmp1) { /* C is bigger than Fortran */ 1459306f9a3SSatish Balay tmp2 = (tmp3 - tmp1)/sizeof(PetscScalar); 1469306f9a3SSatish Balay itmp2 = (size_t) tmp2; 147f91d1997SBarry Smith shift = (align*sizeof(PetscScalar) - (PetscInt)((tmp3 - tmp1) % (align*sizeof(PetscScalar)))) % (align*sizeof(PetscScalar)); 1489306f9a3SSatish Balay } else { 1499306f9a3SSatish Balay tmp2 = (tmp1 - tmp3)/sizeof(PetscScalar); 1509306f9a3SSatish Balay itmp2 = -((size_t) tmp2); 151f91d1997SBarry Smith shift = (PetscInt)((tmp1 - tmp3) % (align*sizeof(PetscScalar))); 1529306f9a3SSatish Balay } 1539306f9a3SSatish Balay if (shift) { 1549306f9a3SSatish Balay (*PetscErrorPrintf)("PetscScalarAddressToFortran:C and Fortran arrays are\n"); 1559306f9a3SSatish Balay (*PetscErrorPrintf)("not commonly aligned.\n"); 156bcaeba4dSBarry Smith (*PetscErrorPrintf)("Locations/sizeof(PetscScalar): C %f Fortran %f\n",((PetscReal)tmp3)/(PetscReal)sizeof(PetscScalar),((PetscReal)tmp1)/(PetscReal)sizeof(PetscScalar)); 157*41e02c4dSJunchao Zhang PETSCABORT(PETSC_COMM_WORLD,PETSC_ERR_PLIB); 1589306f9a3SSatish Balay } 1591e2582c4SBarry Smith ierr = PetscInfo(obj,"Efficiency warning, copying array in XXXGetArray() due\n\ 160ae15b995SBarry Smith to alignment differences between C and Fortran\n");CHKERRQ(ierr); 1619306f9a3SSatish Balay } 1629306f9a3SSatish Balay *res = itmp2; 1639306f9a3SSatish Balay return 0; 1649306f9a3SSatish Balay } 1659306f9a3SSatish Balay 1669306f9a3SSatish Balay /* 1679306f9a3SSatish Balay obj - the PETSc object where the scalar pointer came from 1689306f9a3SSatish Balay base - the Fortran array address 1699306f9a3SSatish Balay addr - the Fortran offset from base 1709306f9a3SSatish Balay N - the amount of data 1719306f9a3SSatish Balay 1729306f9a3SSatish Balay lx - the array space that is to be passed to XXXXRestoreArray() 1739306f9a3SSatish Balay */ 1749306f9a3SSatish Balay PetscErrorCode PetscScalarAddressFromFortran(PetscObject obj,PetscScalar *base,size_t addr,PetscInt N,PetscScalar **lx) 1759306f9a3SSatish Balay { 1769306f9a3SSatish Balay PetscErrorCode ierr; 1779306f9a3SSatish Balay PetscInt shift; 178776b82aeSLisandro Dalcin PetscContainer container; 1799306f9a3SSatish Balay PetscScalar *tlx; 1809306f9a3SSatish Balay 1819306f9a3SSatish Balay ierr = PetscObjectQuery(obj,"GetArrayPtr",(PetscObject*)&container);CHKERRQ(ierr); 1829306f9a3SSatish Balay if (container) { 183776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void**)lx);CHKERRQ(ierr); 1849306f9a3SSatish Balay tlx = base + addr; 1859306f9a3SSatish Balay 1869306f9a3SSatish Balay shift = *(PetscInt*)*lx; 187580bdb30SBarry Smith ierr = PetscArraycpy(*lx,tlx,N);CHKERRQ(ierr); 1889306f9a3SSatish Balay tlx = (PetscScalar*)(((char*)tlx) - shift); 189a297a907SKarl Rupp 1909306f9a3SSatish Balay ierr = PetscFree(tlx);CHKERRQ(ierr); 1916bf464f9SBarry Smith ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 1929306f9a3SSatish Balay ierr = PetscObjectCompose(obj,"GetArrayPtr",0);CHKERRQ(ierr); 1939306f9a3SSatish Balay } else { 1949306f9a3SSatish Balay *lx = base + addr; 1959306f9a3SSatish Balay } 1969306f9a3SSatish Balay return 0; 1979306f9a3SSatish Balay } 1989306f9a3SSatish Balay 19983886165SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS) 200f66fdb6dSSatish Balay #define petscisinfornanscalar_ PETSCISINFORNANSCALAR 201f66fdb6dSSatish Balay #define petscisinfornanreal_ PETSCISINFORNANREAL 20283886165SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 203f66fdb6dSSatish Balay #define petscisinfornanscalar_ petscisinfornanscalar 204f66fdb6dSSatish Balay #define petscisinfornanreal_ petscisinfornanreal 20583886165SBarry Smith #endif 20683886165SBarry Smith 2078cc058d9SJed Brown PETSC_EXTERN PetscBool PETSC_STDCALL petscisinfornanscalar_(PetscScalar *v) 20883886165SBarry Smith { 209ace3abfcSBarry Smith return (PetscBool) PetscIsInfOrNanScalar(*v); 210f66fdb6dSSatish Balay } 211f66fdb6dSSatish Balay 2128cc058d9SJed Brown PETSC_EXTERN PetscBool PETSC_STDCALL petscisinfornanreal_(PetscReal *v) 213f66fdb6dSSatish Balay { 214ace3abfcSBarry Smith return (PetscBool) PetscIsInfOrNanReal(*v); 21583886165SBarry Smith } 21683886165SBarry Smith 21783886165SBarry Smith 2189306f9a3SSatish Balay 219