16dd63270SBarry Smith #include <petscmat.h> 26dd63270SBarry Smith #include <petsc/private/ftnimpl.h> 36dd63270SBarry Smith 46dd63270SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS) 5*f14a4c02SBarry Smith #define matgetrow_ MATGETROW 6*f14a4c02SBarry Smith #define matrestorerow_ MATRESTOREROW 76dd63270SBarry Smith #define matmpiaijgetseqaij_ MATMPIAIJGETSEQAIJ 86dd63270SBarry Smith #define matmpiaijrestoreseqaij_ MATMPIAIJRESTORESEQAIJ 96dd63270SBarry Smith #define matdensegetarray1d_ MATDENSEGETARRAY1D 106dd63270SBarry Smith #define matdenserestorearray1d_ MATDENSERESTOREARRAY1D 116dd63270SBarry Smith #define matdensegetarrayread1d_ MATDENSEGETARRAYREAD1D 126dd63270SBarry Smith #define matdenserestorearrayread1d_ MATDENSERESTOREARRAYREAD1D 136dd63270SBarry Smith #define matdensegetarraywrite1d_ MATDENSEGETARRAYWRITE1D 146dd63270SBarry Smith #define matdenserestorearraywrite1d_ MATDENSERESTOREARRAYWRITE1D 156dd63270SBarry Smith #define matdensegetarray2d_ MATDENSEGETARRAY2D 166dd63270SBarry Smith #define matdenserestorearray2d_ MATDENSERESTOREARRAY2D 176dd63270SBarry Smith #define matdensegetarrayread2d_ MATDENSEGETARRAYREAD2D 186dd63270SBarry Smith #define matdenserestorearrayread2d_ MATDENSERESTOREARRAYREAD2D 196dd63270SBarry Smith #define matdensegetarraywrite2d_ MATDENSEGETARRAYWRITE2D 206dd63270SBarry Smith #define matdenserestorearraywrite2d_ MATDENSERESTOREARRAYWRITE2D 216dd63270SBarry Smith #define matdensegetcolumn_ MATDENSEGETCOLUMN 226dd63270SBarry Smith #define matdenserestorecolumn_ MATDENSERESTORECOLUMN 236dd63270SBarry Smith #define matseqaijgetarray_ MATSEQAIJGETARRAY 246dd63270SBarry Smith #define matseqaijrestorearray_ MATSEQAIJRESTOREARRAY 256dd63270SBarry Smith #define matgetghosts_ MATGETGHOSTS 266dd63270SBarry Smith #define matgetrowij_ MATGETROWIJ 276dd63270SBarry Smith #define matrestorerowij_ MATRESTOREROWIJ 286dd63270SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE) 29*f14a4c02SBarry Smith #define matgetrow_ matgetrow 30*f14a4c02SBarry Smith #define matrestorerow_ matrestorerow 316dd63270SBarry Smith #define matmpiaijgetseqaij_ matmpiaijgetseqaij 326dd63270SBarry Smith #define matmpiaijrestoreseqaij_ matmpiaijrestoreseqaij 336dd63270SBarry Smith #define matdensegetarray1d_ matdensegetarray1d 346dd63270SBarry Smith #define matdenserestorearray1d_ matdenserestorearray1d 356dd63270SBarry Smith #define matdensegetarrayread1d_ matdensegetarrayread1d 366dd63270SBarry Smith #define matdenserestorearrayread1d_ matdenserestorearrayread1d 376dd63270SBarry Smith #define matdensegetarraywrite1d_ matdensegetarraywrite1d 386dd63270SBarry Smith #define matdenserestorearraywrite1d_ matdenserestorearraywrite1d 396dd63270SBarry Smith #define matdensegetarray2d_ matdensegetarray2d 406dd63270SBarry Smith #define matdenserestorearray2d_ matdenserestorearray2d 416dd63270SBarry Smith #define matdensegetarrayread2d_ matdensegetarrayread2d 426dd63270SBarry Smith #define matdenserestorearrayread2d_ matdenserestorearrayread2d 436dd63270SBarry Smith #define matdensegetarraywrite2d_ matdensegetarraywrite2d 446dd63270SBarry Smith #define matdenserestorearraywrite2d_ matdenserestorearraywrite2d 456dd63270SBarry Smith #define matdensegetcolumn_ matdensegetcolumn 466dd63270SBarry Smith #define matdenserestorecolumn_ matdenserestorecolumn 476dd63270SBarry Smith #define matseqaijgetarray_ matseqaijgetarray 486dd63270SBarry Smith #define matseqaijrestorearray_ matseqaijrestorearray 496dd63270SBarry Smith #define matgetghosts_ matgetghosts 506dd63270SBarry Smith #define matgetrowij_ matgetrowij 516dd63270SBarry Smith #define matrestorerowij_ matrestorerowij 526dd63270SBarry Smith #endif 536dd63270SBarry Smith 54*f14a4c02SBarry Smith PETSC_EXTERN void matgetrow_(Mat *B, PetscInt *row, PetscInt *N, F90Array1d *ia, F90Array1d *a, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(iad) PETSC_F90_2PTR_PROTO(jad)) 55*f14a4c02SBarry Smith { 56*f14a4c02SBarry Smith PetscInt n; 57*f14a4c02SBarry Smith const PetscInt *II = NULL; 58*f14a4c02SBarry Smith const PetscScalar *A = NULL; 59*f14a4c02SBarry Smith 60*f14a4c02SBarry Smith if (FORTRANNULLINTEGERPOINTER(ia) && FORTRANNULLSCALARPOINTER(a)) { 61*f14a4c02SBarry Smith *ierr = MatGetRow(*B, *row, &n, NULL, NULL); 62*f14a4c02SBarry Smith } else if (FORTRANNULLINTEGERPOINTER(ia)) { 63*f14a4c02SBarry Smith *ierr = MatGetRow(*B, *row, &n, NULL, &A); 64*f14a4c02SBarry Smith } else if (FORTRANNULLSCALARPOINTER(a)) { 65*f14a4c02SBarry Smith *ierr = MatGetRow(*B, *row, &n, &II, NULL); 66*f14a4c02SBarry Smith } else { 67*f14a4c02SBarry Smith *ierr = MatGetRow(*B, *row, &n, &II, &A); 68*f14a4c02SBarry Smith } 69*f14a4c02SBarry Smith if (*ierr) return; 70*f14a4c02SBarry Smith if (II) { *ierr = F90Array1dCreate((void *)II, MPIU_INT, 1, n, ia PETSC_F90_2PTR_PARAM(iad)); } 71*f14a4c02SBarry Smith if (A) { *ierr = F90Array1dCreate((void *)A, MPIU_SCALAR, 1, n, a PETSC_F90_2PTR_PARAM(jad)); } 72*f14a4c02SBarry Smith if (!FORTRANNULLINTEGER(N)) *N = n; 73*f14a4c02SBarry Smith } 74*f14a4c02SBarry Smith PETSC_EXTERN void matrestorerow_(Mat *B, PetscInt *row, PetscInt *N, F90Array1d *ia, F90Array1d *a, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(iad) PETSC_F90_2PTR_PROTO(jad)) 75*f14a4c02SBarry Smith { 76*f14a4c02SBarry Smith const PetscInt *IA = NULL; 77*f14a4c02SBarry Smith const PetscScalar *A = NULL; 78*f14a4c02SBarry Smith PetscInt n; 79*f14a4c02SBarry Smith 80*f14a4c02SBarry Smith if (FORTRANNULLINTEGERPOINTER(ia) && FORTRANNULLSCALARPOINTER(a)) return; 81*f14a4c02SBarry Smith if (!FORTRANNULLINTEGERPOINTER(ia)) { 82*f14a4c02SBarry Smith *ierr = F90Array1dAccess(ia, MPIU_INT, (void **)&IA PETSC_F90_2PTR_PARAM(iad)); 83*f14a4c02SBarry Smith if (*ierr) return; 84*f14a4c02SBarry Smith *ierr = F90Array1dDestroy(ia, MPIU_INT PETSC_F90_2PTR_PARAM(iad)); 85*f14a4c02SBarry Smith if (*ierr) return; 86*f14a4c02SBarry Smith } 87*f14a4c02SBarry Smith if (!FORTRANNULLSCALARPOINTER(a)) { 88*f14a4c02SBarry Smith *ierr = F90Array1dAccess(a, MPIU_SCALAR, (void **)&A PETSC_F90_2PTR_PARAM(jad)); 89*f14a4c02SBarry Smith if (*ierr) return; 90*f14a4c02SBarry Smith *ierr = F90Array1dDestroy(a, MPIU_INT PETSC_F90_2PTR_PARAM(jad)); 91*f14a4c02SBarry Smith if (*ierr) return; 92*f14a4c02SBarry Smith } 93*f14a4c02SBarry Smith if (FORTRANNULLINTEGERPOINTER(ia)) { 94*f14a4c02SBarry Smith *ierr = MatRestoreRow(*B, *row, &n, NULL, &A); 95*f14a4c02SBarry Smith } else if (FORTRANNULLSCALARPOINTER(a)) { 96*f14a4c02SBarry Smith *ierr = MatRestoreRow(*B, *row, &n, &IA, NULL); 97*f14a4c02SBarry Smith } else { 98*f14a4c02SBarry Smith *ierr = MatRestoreRow(*B, *row, &n, &IA, &A); 99*f14a4c02SBarry Smith } 100*f14a4c02SBarry Smith } 1016dd63270SBarry Smith PETSC_EXTERN void matgetghosts_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 1026dd63270SBarry Smith { 1036dd63270SBarry Smith const PetscInt *ghosts; 1046dd63270SBarry Smith PetscInt N; 1056dd63270SBarry Smith 1066dd63270SBarry Smith *ierr = MatGetGhosts(*mat, &N, &ghosts); 1076dd63270SBarry Smith if (*ierr) return; 1086dd63270SBarry Smith *ierr = F90Array1dCreate((PetscInt *)ghosts, MPIU_INT, 1, N, ptr PETSC_F90_2PTR_PARAM(ptrd)); 1096dd63270SBarry Smith } 1106dd63270SBarry Smith PETSC_EXTERN void matdensegetarray2d_(Mat *mat, F90Array2d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 1116dd63270SBarry Smith { 1126dd63270SBarry Smith PetscScalar *fa; 1136dd63270SBarry Smith PetscInt m, N, lda; 1146dd63270SBarry Smith *ierr = MatDenseGetArray(*mat, &fa); 1156dd63270SBarry Smith if (*ierr) return; 1166dd63270SBarry Smith *ierr = MatGetLocalSize(*mat, &m, NULL); 1176dd63270SBarry Smith if (*ierr) return; 1186dd63270SBarry Smith *ierr = MatGetSize(*mat, NULL, &N); 1196dd63270SBarry Smith if (*ierr) return; 1206dd63270SBarry Smith *ierr = MatDenseGetLDA(*mat, &lda); 1216dd63270SBarry Smith if (*ierr) return; 1226dd63270SBarry Smith if (m != lda) { // TODO: add F90Array2dLDACreate() 1236dd63270SBarry Smith *ierr = PetscError(((PetscObject)*mat)->comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_ARG_BADPTR, PETSC_ERROR_INITIAL, "Array lda %" PetscInt_FMT " must match number of local rows %" PetscInt_FMT, lda, m); 1246dd63270SBarry Smith return; 1256dd63270SBarry Smith } 1266dd63270SBarry Smith *ierr = F90Array2dCreate(fa, MPIU_SCALAR, 1, m, 1, N, ptr PETSC_F90_2PTR_PARAM(ptrd)); 1276dd63270SBarry Smith } 1286dd63270SBarry Smith PETSC_EXTERN void matdenserestorearray2d_(Mat *mat, F90Array2d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 1296dd63270SBarry Smith { 1306dd63270SBarry Smith PetscScalar *fa; 1316dd63270SBarry Smith *ierr = F90Array2dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd)); 1326dd63270SBarry Smith if (*ierr) return; 1336dd63270SBarry Smith *ierr = F90Array2dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd)); 1346dd63270SBarry Smith if (*ierr) return; 1356dd63270SBarry Smith *ierr = MatDenseRestoreArray(*mat, &fa); 1366dd63270SBarry Smith } 1376dd63270SBarry Smith PETSC_EXTERN void matdensegetarrayread2d_(Mat *mat, F90Array2d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 1386dd63270SBarry Smith { 1396dd63270SBarry Smith const PetscScalar *fa; 1406dd63270SBarry Smith PetscInt m, N, lda; 1416dd63270SBarry Smith *ierr = MatDenseGetArrayRead(*mat, &fa); 1426dd63270SBarry Smith if (*ierr) return; 1436dd63270SBarry Smith *ierr = MatGetLocalSize(*mat, &m, NULL); 1446dd63270SBarry Smith if (*ierr) return; 1456dd63270SBarry Smith *ierr = MatGetSize(*mat, NULL, &N); 1466dd63270SBarry Smith if (*ierr) return; 1476dd63270SBarry Smith *ierr = MatDenseGetLDA(*mat, &lda); 1486dd63270SBarry Smith if (*ierr) return; 1496dd63270SBarry Smith if (m != lda) { // TODO: add F90Array2dLDACreate() 1506dd63270SBarry Smith *ierr = PetscError(((PetscObject)*mat)->comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_ARG_BADPTR, PETSC_ERROR_INITIAL, "Array lda %" PetscInt_FMT " must match number of local rows %" PetscInt_FMT, lda, m); 1516dd63270SBarry Smith return; 1526dd63270SBarry Smith } 1536dd63270SBarry Smith *ierr = F90Array2dCreate((void **)fa, MPIU_SCALAR, 1, m, 1, N, ptr PETSC_F90_2PTR_PARAM(ptrd)); 1546dd63270SBarry Smith } 1556dd63270SBarry Smith PETSC_EXTERN void matdenserestorearrayread2d_(Mat *mat, F90Array2d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 1566dd63270SBarry Smith { 1576dd63270SBarry Smith const PetscScalar *fa; 1586dd63270SBarry Smith *ierr = F90Array2dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd)); 1596dd63270SBarry Smith if (*ierr) return; 1606dd63270SBarry Smith *ierr = F90Array2dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd)); 1616dd63270SBarry Smith if (*ierr) return; 1626dd63270SBarry Smith *ierr = MatDenseRestoreArrayRead(*mat, &fa); 1636dd63270SBarry Smith } 1646dd63270SBarry Smith PETSC_EXTERN void matdensegetarraywrite2d_(Mat *mat, F90Array2d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 1656dd63270SBarry Smith { 1666dd63270SBarry Smith PetscScalar *fa; 1676dd63270SBarry Smith PetscInt m, N, lda; 1686dd63270SBarry Smith *ierr = MatDenseGetArrayWrite(*mat, &fa); 1696dd63270SBarry Smith if (*ierr) return; 1706dd63270SBarry Smith *ierr = MatGetLocalSize(*mat, &m, NULL); 1716dd63270SBarry Smith if (*ierr) return; 1726dd63270SBarry Smith *ierr = MatGetSize(*mat, NULL, &N); 1736dd63270SBarry Smith if (*ierr) return; 1746dd63270SBarry Smith *ierr = MatDenseGetLDA(*mat, &lda); 1756dd63270SBarry Smith if (*ierr) return; 1766dd63270SBarry Smith if (m != lda) { // TODO: add F90Array2dLDACreate() 1776dd63270SBarry Smith *ierr = PetscError(((PetscObject)*mat)->comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_ARG_BADPTR, PETSC_ERROR_INITIAL, "Array lda %" PetscInt_FMT " must match number of local rows %" PetscInt_FMT, lda, m); 1786dd63270SBarry Smith return; 1796dd63270SBarry Smith } 1806dd63270SBarry Smith *ierr = F90Array2dCreate(fa, MPIU_SCALAR, 1, m, 1, N, ptr PETSC_F90_2PTR_PARAM(ptrd)); 1816dd63270SBarry Smith } 1826dd63270SBarry Smith PETSC_EXTERN void matdenserestorearraywrite2d_(Mat *mat, F90Array2d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 1836dd63270SBarry Smith { 1846dd63270SBarry Smith PetscScalar *fa; 1856dd63270SBarry Smith *ierr = F90Array2dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd)); 1866dd63270SBarry Smith if (*ierr) return; 1876dd63270SBarry Smith *ierr = F90Array2dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd)); 1886dd63270SBarry Smith if (*ierr) return; 1896dd63270SBarry Smith *ierr = MatDenseRestoreArrayWrite(*mat, &fa); 1906dd63270SBarry Smith } 1916dd63270SBarry Smith PETSC_EXTERN void matdensegetarray1d_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 1926dd63270SBarry Smith { 1936dd63270SBarry Smith PetscScalar *fa; 1946dd63270SBarry Smith PetscInt m, N, lda; 1956dd63270SBarry Smith *ierr = MatDenseGetArray(*mat, &fa); 1966dd63270SBarry Smith if (*ierr) return; 1976dd63270SBarry Smith *ierr = MatGetLocalSize(*mat, &m, NULL); 1986dd63270SBarry Smith if (*ierr) return; 1996dd63270SBarry Smith *ierr = MatGetSize(*mat, NULL, &N); 2006dd63270SBarry Smith if (*ierr) return; 2016dd63270SBarry Smith *ierr = MatDenseGetLDA(*mat, &lda); 2026dd63270SBarry Smith if (*ierr) return; 2036dd63270SBarry Smith if (m != lda) { // TODO: add F90Array1dLDACreate() 2046dd63270SBarry Smith *ierr = PetscError(((PetscObject)*mat)->comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_ARG_BADPTR, PETSC_ERROR_INITIAL, "Array lda %" PetscInt_FMT " must match number of local rows %" PetscInt_FMT, lda, m); 2056dd63270SBarry Smith return; 2066dd63270SBarry Smith } 2076dd63270SBarry Smith *ierr = F90Array1dCreate(fa, MPIU_SCALAR, 1, m * N, ptr PETSC_F90_2PTR_PARAM(ptrd)); 2086dd63270SBarry Smith } 2096dd63270SBarry Smith PETSC_EXTERN void matdenserestorearray1d_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 2106dd63270SBarry Smith { 2116dd63270SBarry Smith PetscScalar *fa; 2126dd63270SBarry Smith *ierr = F90Array1dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd)); 2136dd63270SBarry Smith if (*ierr) return; 2146dd63270SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd)); 2156dd63270SBarry Smith if (*ierr) return; 2166dd63270SBarry Smith *ierr = MatDenseRestoreArray(*mat, &fa); 2176dd63270SBarry Smith } 2186dd63270SBarry Smith PETSC_EXTERN void matdensegetarrayread1d_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 2196dd63270SBarry Smith { 2206dd63270SBarry Smith const PetscScalar *fa; 2216dd63270SBarry Smith PetscInt m, N, lda; 2226dd63270SBarry Smith *ierr = MatDenseGetArrayRead(*mat, &fa); 2236dd63270SBarry Smith if (*ierr) return; 2246dd63270SBarry Smith *ierr = MatGetLocalSize(*mat, &m, NULL); 2256dd63270SBarry Smith if (*ierr) return; 2266dd63270SBarry Smith *ierr = MatGetSize(*mat, NULL, &N); 2276dd63270SBarry Smith if (*ierr) return; 2286dd63270SBarry Smith *ierr = MatDenseGetLDA(*mat, &lda); 2296dd63270SBarry Smith if (*ierr) return; 2306dd63270SBarry Smith if (m != lda) { // TODO: add F90Array1dLDACreate() 2316dd63270SBarry Smith *ierr = PetscError(((PetscObject)*mat)->comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_ARG_BADPTR, PETSC_ERROR_INITIAL, "Array lda %" PetscInt_FMT " must match number of local rows %" PetscInt_FMT, lda, m); 2326dd63270SBarry Smith return; 2336dd63270SBarry Smith } 2346dd63270SBarry Smith *ierr = F90Array1dCreate((void **)fa, MPIU_SCALAR, 1, m * N, ptr PETSC_F90_2PTR_PARAM(ptrd)); 2356dd63270SBarry Smith } 2366dd63270SBarry Smith PETSC_EXTERN void matdenserestorearrayread1d_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 2376dd63270SBarry Smith { 2386dd63270SBarry Smith const PetscScalar *fa; 2396dd63270SBarry Smith *ierr = F90Array1dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd)); 2406dd63270SBarry Smith if (*ierr) return; 2416dd63270SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd)); 2426dd63270SBarry Smith if (*ierr) return; 2436dd63270SBarry Smith *ierr = MatDenseRestoreArrayRead(*mat, &fa); 2446dd63270SBarry Smith } 2456dd63270SBarry Smith PETSC_EXTERN void matdensegetarraywrite1d_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 2466dd63270SBarry Smith { 2476dd63270SBarry Smith PetscScalar *fa; 2486dd63270SBarry Smith PetscInt m, N, lda; 2496dd63270SBarry Smith *ierr = MatDenseGetArrayWrite(*mat, &fa); 2506dd63270SBarry Smith if (*ierr) return; 2516dd63270SBarry Smith *ierr = MatGetLocalSize(*mat, &m, NULL); 2526dd63270SBarry Smith if (*ierr) return; 2536dd63270SBarry Smith *ierr = MatGetSize(*mat, NULL, &N); 2546dd63270SBarry Smith if (*ierr) return; 2556dd63270SBarry Smith *ierr = MatDenseGetLDA(*mat, &lda); 2566dd63270SBarry Smith if (*ierr) return; 2576dd63270SBarry Smith if (m != lda) { // TODO: add F90Array1dLDACreate() 2586dd63270SBarry Smith *ierr = PetscError(((PetscObject)*mat)->comm, __LINE__, PETSC_FUNCTION_NAME, __FILE__, PETSC_ERR_ARG_BADPTR, PETSC_ERROR_INITIAL, "Array lda %" PetscInt_FMT " must match number of local rows %" PetscInt_FMT, lda, m); 2596dd63270SBarry Smith return; 2606dd63270SBarry Smith } 2616dd63270SBarry Smith *ierr = F90Array1dCreate(fa, MPIU_SCALAR, 1, m * N, ptr PETSC_F90_2PTR_PARAM(ptrd)); 2626dd63270SBarry Smith } 2636dd63270SBarry Smith PETSC_EXTERN void matdenserestorearraywrite1d_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 2646dd63270SBarry Smith { 2656dd63270SBarry Smith PetscScalar *fa; 2666dd63270SBarry Smith *ierr = F90Array1dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd)); 2676dd63270SBarry Smith if (*ierr) return; 2686dd63270SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd)); 2696dd63270SBarry Smith if (*ierr) return; 2706dd63270SBarry Smith *ierr = MatDenseRestoreArrayWrite(*mat, &fa); 2716dd63270SBarry Smith } 2726dd63270SBarry Smith PETSC_EXTERN void matdensegetcolumn_(Mat *mat, PetscInt *col, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 2736dd63270SBarry Smith { 2746dd63270SBarry Smith PetscScalar *fa; 2756dd63270SBarry Smith PetscInt m; 2766dd63270SBarry Smith *ierr = MatDenseGetColumn(*mat, *col, &fa); 2776dd63270SBarry Smith if (*ierr) return; 2786dd63270SBarry Smith *ierr = MatGetLocalSize(*mat, &m, NULL); 2796dd63270SBarry Smith if (*ierr) return; 2806dd63270SBarry Smith *ierr = F90Array1dCreate(fa, MPIU_SCALAR, 1, m, ptr PETSC_F90_2PTR_PARAM(ptrd)); 2816dd63270SBarry Smith } 2826dd63270SBarry Smith PETSC_EXTERN void matdenserestorecolumn_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 2836dd63270SBarry Smith { 2846dd63270SBarry Smith PetscScalar *fa; 2856dd63270SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd)); 2866dd63270SBarry Smith if (*ierr) return; 2876dd63270SBarry Smith *ierr = MatDenseRestoreColumn(*mat, &fa); 2886dd63270SBarry Smith } 2896dd63270SBarry Smith PETSC_EXTERN void matseqaijgetarray_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 2906dd63270SBarry Smith { 2916dd63270SBarry Smith PetscScalar *fa; 2926dd63270SBarry Smith PetscInt m, n; 2936dd63270SBarry Smith *ierr = MatSeqAIJGetArray(*mat, &fa); 2946dd63270SBarry Smith if (*ierr) return; 2956dd63270SBarry Smith *ierr = MatGetLocalSize(*mat, &m, &n); 2966dd63270SBarry Smith if (*ierr) return; 2976dd63270SBarry Smith *ierr = F90Array1dCreate(fa, MPIU_SCALAR, 1, m * n, ptr PETSC_F90_2PTR_PARAM(ptrd)); 2986dd63270SBarry Smith } 2996dd63270SBarry Smith PETSC_EXTERN void matseqaijrestorearray_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 3006dd63270SBarry Smith { 3016dd63270SBarry Smith PetscScalar *fa; 3026dd63270SBarry Smith *ierr = F90Array1dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd)); 3036dd63270SBarry Smith if (*ierr) return; 3046dd63270SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd)); 3056dd63270SBarry Smith if (*ierr) return; 3066dd63270SBarry Smith *ierr = MatSeqAIJRestoreArray(*mat, &fa); 3076dd63270SBarry Smith } 3086dd63270SBarry Smith PETSC_EXTERN void matgetrowij_(Mat *B, PetscInt *shift, PetscBool *sym, PetscBool *blockcompressed, PetscInt *n, F90Array1d *ia, F90Array1d *ja, PetscBool *done, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(iad) PETSC_F90_2PTR_PROTO(jad)) 3096dd63270SBarry Smith { 3106dd63270SBarry Smith const PetscInt *IA, *JA; 3116dd63270SBarry Smith *ierr = MatGetRowIJ(*B, *shift, *sym, *blockcompressed, n, &IA, &JA, done); 3126dd63270SBarry Smith if (*ierr) return; 3136dd63270SBarry Smith if (!*done) return; 3146dd63270SBarry Smith *ierr = F90Array1dCreate((PetscInt *)IA, MPIU_INT, 1, *n + 1, ia PETSC_F90_2PTR_PARAM(iad)); 3156dd63270SBarry Smith *ierr = F90Array1dCreate((PetscInt *)JA, MPIU_INT, 1, IA[*n], ja PETSC_F90_2PTR_PARAM(jad)); 3166dd63270SBarry Smith } 3176dd63270SBarry Smith PETSC_EXTERN void matrestorerowij_(Mat *B, PetscInt *shift, PetscBool *sym, PetscBool *blockcompressed, PetscInt *n, F90Array1d *ia, F90Array1d *ja, PetscBool *done, PetscErrorCode *ierr PETSC_F90_2PTR_PROTO(iad) PETSC_F90_2PTR_PROTO(jad)) 3186dd63270SBarry Smith { 3196dd63270SBarry Smith const PetscInt *IA, *JA; 3206dd63270SBarry Smith *ierr = F90Array1dAccess(ia, MPIU_INT, (void **)&IA PETSC_F90_2PTR_PARAM(iad)); 3216dd63270SBarry Smith if (*ierr) return; 3226dd63270SBarry Smith *ierr = F90Array1dDestroy(ia, MPIU_INT PETSC_F90_2PTR_PARAM(iad)); 3236dd63270SBarry Smith if (*ierr) return; 3246dd63270SBarry Smith *ierr = F90Array1dAccess(ja, MPIU_INT, (void **)&JA PETSC_F90_2PTR_PARAM(jad)); 3256dd63270SBarry Smith if (*ierr) return; 3266dd63270SBarry Smith *ierr = F90Array1dDestroy(ja, MPIU_INT PETSC_F90_2PTR_PARAM(jad)); 3276dd63270SBarry Smith if (*ierr) return; 3286dd63270SBarry Smith *ierr = MatRestoreRowIJ(*B, *shift, *sym, *blockcompressed, n, &IA, &JA, done); 3296dd63270SBarry Smith } 3306dd63270SBarry Smith PETSC_EXTERN void matmpiaijgetseqaij_(Mat *mat, Mat *A, Mat *B, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 3316dd63270SBarry Smith { 3326dd63270SBarry Smith const PetscInt *fa; 3336dd63270SBarry Smith PetscInt n; 3346dd63270SBarry Smith *ierr = MatMPIAIJGetSeqAIJ(*mat, A, B, &fa); 3356dd63270SBarry Smith if (*ierr) return; 3366dd63270SBarry Smith *ierr = MatGetLocalSize(*B, NULL, &n); 3376dd63270SBarry Smith if (*ierr) return; 3386dd63270SBarry Smith *ierr = F90Array1dCreate((void *)fa, MPIU_INT, 1, n, ptr PETSC_F90_2PTR_PARAM(ptrd)); 3396dd63270SBarry Smith } 3406dd63270SBarry Smith PETSC_EXTERN void matmpiaijrestoreseqaij_(Mat *mat, Mat *A, Mat *B, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd)) 3416dd63270SBarry Smith { 3426dd63270SBarry Smith *ierr = F90Array1dDestroy(ptr, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd)); 3436dd63270SBarry Smith if (*ierr) return; 3446dd63270SBarry Smith } 345