xref: /petsc/src/mat/interface/ftn-custom/zmatrixf90.c (revision 6dd63270497ad23dcf16ae500a87ff2b2a0b7474)
1*6dd63270SBarry Smith #include <petscmat.h>
2*6dd63270SBarry Smith #include <petsc/private/ftnimpl.h>
3*6dd63270SBarry Smith 
4*6dd63270SBarry Smith #if defined(PETSC_HAVE_FORTRAN_CAPS)
5*6dd63270SBarry Smith   #define matmpiaijgetseqaij_          MATMPIAIJGETSEQAIJ
6*6dd63270SBarry Smith   #define matmpiaijrestoreseqaij_      MATMPIAIJRESTORESEQAIJ
7*6dd63270SBarry Smith   #define matdensegetarray1d_          MATDENSEGETARRAY1D
8*6dd63270SBarry Smith   #define matdenserestorearray1d_      MATDENSERESTOREARRAY1D
9*6dd63270SBarry Smith   #define matdensegetarrayread1d_      MATDENSEGETARRAYREAD1D
10*6dd63270SBarry Smith   #define matdenserestorearrayread1d_  MATDENSERESTOREARRAYREAD1D
11*6dd63270SBarry Smith   #define matdensegetarraywrite1d_     MATDENSEGETARRAYWRITE1D
12*6dd63270SBarry Smith   #define matdenserestorearraywrite1d_ MATDENSERESTOREARRAYWRITE1D
13*6dd63270SBarry Smith   #define matdensegetarray2d_          MATDENSEGETARRAY2D
14*6dd63270SBarry Smith   #define matdenserestorearray2d_      MATDENSERESTOREARRAY2D
15*6dd63270SBarry Smith   #define matdensegetarrayread2d_      MATDENSEGETARRAYREAD2D
16*6dd63270SBarry Smith   #define matdenserestorearrayread2d_  MATDENSERESTOREARRAYREAD2D
17*6dd63270SBarry Smith   #define matdensegetarraywrite2d_     MATDENSEGETARRAYWRITE2D
18*6dd63270SBarry Smith   #define matdenserestorearraywrite2d_ MATDENSERESTOREARRAYWRITE2D
19*6dd63270SBarry Smith   #define matdensegetcolumn_           MATDENSEGETCOLUMN
20*6dd63270SBarry Smith   #define matdenserestorecolumn_       MATDENSERESTORECOLUMN
21*6dd63270SBarry Smith   #define matseqaijgetarray_           MATSEQAIJGETARRAY
22*6dd63270SBarry Smith   #define matseqaijrestorearray_       MATSEQAIJRESTOREARRAY
23*6dd63270SBarry Smith   #define matgetghosts_                MATGETGHOSTS
24*6dd63270SBarry Smith   #define matgetrowij_                 MATGETROWIJ
25*6dd63270SBarry Smith   #define matrestorerowij_             MATRESTOREROWIJ
26*6dd63270SBarry Smith #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
27*6dd63270SBarry Smith   #define matmpiaijgetseqaij_          matmpiaijgetseqaij
28*6dd63270SBarry Smith   #define matmpiaijrestoreseqaij_      matmpiaijrestoreseqaij
29*6dd63270SBarry Smith   #define matdensegetarray1d_          matdensegetarray1d
30*6dd63270SBarry Smith   #define matdenserestorearray1d_      matdenserestorearray1d
31*6dd63270SBarry Smith   #define matdensegetarrayread1d_      matdensegetarrayread1d
32*6dd63270SBarry Smith   #define matdenserestorearrayread1d_  matdenserestorearrayread1d
33*6dd63270SBarry Smith   #define matdensegetarraywrite1d_     matdensegetarraywrite1d
34*6dd63270SBarry Smith   #define matdenserestorearraywrite1d_ matdenserestorearraywrite1d
35*6dd63270SBarry Smith   #define matdensegetarray2d_          matdensegetarray2d
36*6dd63270SBarry Smith   #define matdenserestorearray2d_      matdenserestorearray2d
37*6dd63270SBarry Smith   #define matdensegetarrayread2d_      matdensegetarrayread2d
38*6dd63270SBarry Smith   #define matdenserestorearrayread2d_  matdenserestorearrayread2d
39*6dd63270SBarry Smith   #define matdensegetarraywrite2d_     matdensegetarraywrite2d
40*6dd63270SBarry Smith   #define matdenserestorearraywrite2d_ matdenserestorearraywrite2d
41*6dd63270SBarry Smith   #define matdensegetcolumn_           matdensegetcolumn
42*6dd63270SBarry Smith   #define matdenserestorecolumn_       matdenserestorecolumn
43*6dd63270SBarry Smith   #define matseqaijgetarray_           matseqaijgetarray
44*6dd63270SBarry Smith   #define matseqaijrestorearray_       matseqaijrestorearray
45*6dd63270SBarry Smith   #define matgetghosts_                matgetghosts
46*6dd63270SBarry Smith   #define matgetrowij_                 matgetrowij
47*6dd63270SBarry Smith   #define matrestorerowij_             matrestorerowij
48*6dd63270SBarry Smith #endif
49*6dd63270SBarry Smith 
50*6dd63270SBarry Smith PETSC_EXTERN void matgetghosts_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
51*6dd63270SBarry Smith {
52*6dd63270SBarry Smith   const PetscInt *ghosts;
53*6dd63270SBarry Smith   PetscInt        N;
54*6dd63270SBarry Smith 
55*6dd63270SBarry Smith   *ierr = MatGetGhosts(*mat, &N, &ghosts);
56*6dd63270SBarry Smith   if (*ierr) return;
57*6dd63270SBarry Smith   *ierr = F90Array1dCreate((PetscInt *)ghosts, MPIU_INT, 1, N, ptr PETSC_F90_2PTR_PARAM(ptrd));
58*6dd63270SBarry Smith }
59*6dd63270SBarry Smith PETSC_EXTERN void matdensegetarray2d_(Mat *mat, F90Array2d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
60*6dd63270SBarry Smith {
61*6dd63270SBarry Smith   PetscScalar *fa;
62*6dd63270SBarry Smith   PetscInt     m, N, lda;
63*6dd63270SBarry Smith   *ierr = MatDenseGetArray(*mat, &fa);
64*6dd63270SBarry Smith   if (*ierr) return;
65*6dd63270SBarry Smith   *ierr = MatGetLocalSize(*mat, &m, NULL);
66*6dd63270SBarry Smith   if (*ierr) return;
67*6dd63270SBarry Smith   *ierr = MatGetSize(*mat, NULL, &N);
68*6dd63270SBarry Smith   if (*ierr) return;
69*6dd63270SBarry Smith   *ierr = MatDenseGetLDA(*mat, &lda);
70*6dd63270SBarry Smith   if (*ierr) return;
71*6dd63270SBarry Smith   if (m != lda) { // TODO: add F90Array2dLDACreate()
72*6dd63270SBarry 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);
73*6dd63270SBarry Smith     return;
74*6dd63270SBarry Smith   }
75*6dd63270SBarry Smith   *ierr = F90Array2dCreate(fa, MPIU_SCALAR, 1, m, 1, N, ptr PETSC_F90_2PTR_PARAM(ptrd));
76*6dd63270SBarry Smith }
77*6dd63270SBarry Smith PETSC_EXTERN void matdenserestorearray2d_(Mat *mat, F90Array2d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
78*6dd63270SBarry Smith {
79*6dd63270SBarry Smith   PetscScalar *fa;
80*6dd63270SBarry Smith   *ierr = F90Array2dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd));
81*6dd63270SBarry Smith   if (*ierr) return;
82*6dd63270SBarry Smith   *ierr = F90Array2dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
83*6dd63270SBarry Smith   if (*ierr) return;
84*6dd63270SBarry Smith   *ierr = MatDenseRestoreArray(*mat, &fa);
85*6dd63270SBarry Smith }
86*6dd63270SBarry Smith PETSC_EXTERN void matdensegetarrayread2d_(Mat *mat, F90Array2d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
87*6dd63270SBarry Smith {
88*6dd63270SBarry Smith   const PetscScalar *fa;
89*6dd63270SBarry Smith   PetscInt           m, N, lda;
90*6dd63270SBarry Smith   *ierr = MatDenseGetArrayRead(*mat, &fa);
91*6dd63270SBarry Smith   if (*ierr) return;
92*6dd63270SBarry Smith   *ierr = MatGetLocalSize(*mat, &m, NULL);
93*6dd63270SBarry Smith   if (*ierr) return;
94*6dd63270SBarry Smith   *ierr = MatGetSize(*mat, NULL, &N);
95*6dd63270SBarry Smith   if (*ierr) return;
96*6dd63270SBarry Smith   *ierr = MatDenseGetLDA(*mat, &lda);
97*6dd63270SBarry Smith   if (*ierr) return;
98*6dd63270SBarry Smith   if (m != lda) { // TODO: add F90Array2dLDACreate()
99*6dd63270SBarry 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);
100*6dd63270SBarry Smith     return;
101*6dd63270SBarry Smith   }
102*6dd63270SBarry Smith   *ierr = F90Array2dCreate((void **)fa, MPIU_SCALAR, 1, m, 1, N, ptr PETSC_F90_2PTR_PARAM(ptrd));
103*6dd63270SBarry Smith }
104*6dd63270SBarry Smith PETSC_EXTERN void matdenserestorearrayread2d_(Mat *mat, F90Array2d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
105*6dd63270SBarry Smith {
106*6dd63270SBarry Smith   const PetscScalar *fa;
107*6dd63270SBarry Smith   *ierr = F90Array2dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd));
108*6dd63270SBarry Smith   if (*ierr) return;
109*6dd63270SBarry Smith   *ierr = F90Array2dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
110*6dd63270SBarry Smith   if (*ierr) return;
111*6dd63270SBarry Smith   *ierr = MatDenseRestoreArrayRead(*mat, &fa);
112*6dd63270SBarry Smith }
113*6dd63270SBarry Smith PETSC_EXTERN void matdensegetarraywrite2d_(Mat *mat, F90Array2d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
114*6dd63270SBarry Smith {
115*6dd63270SBarry Smith   PetscScalar *fa;
116*6dd63270SBarry Smith   PetscInt     m, N, lda;
117*6dd63270SBarry Smith   *ierr = MatDenseGetArrayWrite(*mat, &fa);
118*6dd63270SBarry Smith   if (*ierr) return;
119*6dd63270SBarry Smith   *ierr = MatGetLocalSize(*mat, &m, NULL);
120*6dd63270SBarry Smith   if (*ierr) return;
121*6dd63270SBarry Smith   *ierr = MatGetSize(*mat, NULL, &N);
122*6dd63270SBarry Smith   if (*ierr) return;
123*6dd63270SBarry Smith   *ierr = MatDenseGetLDA(*mat, &lda);
124*6dd63270SBarry Smith   if (*ierr) return;
125*6dd63270SBarry Smith   if (m != lda) { // TODO: add F90Array2dLDACreate()
126*6dd63270SBarry 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);
127*6dd63270SBarry Smith     return;
128*6dd63270SBarry Smith   }
129*6dd63270SBarry Smith   *ierr = F90Array2dCreate(fa, MPIU_SCALAR, 1, m, 1, N, ptr PETSC_F90_2PTR_PARAM(ptrd));
130*6dd63270SBarry Smith }
131*6dd63270SBarry Smith PETSC_EXTERN void matdenserestorearraywrite2d_(Mat *mat, F90Array2d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
132*6dd63270SBarry Smith {
133*6dd63270SBarry Smith   PetscScalar *fa;
134*6dd63270SBarry Smith   *ierr = F90Array2dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd));
135*6dd63270SBarry Smith   if (*ierr) return;
136*6dd63270SBarry Smith   *ierr = F90Array2dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
137*6dd63270SBarry Smith   if (*ierr) return;
138*6dd63270SBarry Smith   *ierr = MatDenseRestoreArrayWrite(*mat, &fa);
139*6dd63270SBarry Smith }
140*6dd63270SBarry Smith PETSC_EXTERN void matdensegetarray1d_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
141*6dd63270SBarry Smith {
142*6dd63270SBarry Smith   PetscScalar *fa;
143*6dd63270SBarry Smith   PetscInt     m, N, lda;
144*6dd63270SBarry Smith   *ierr = MatDenseGetArray(*mat, &fa);
145*6dd63270SBarry Smith   if (*ierr) return;
146*6dd63270SBarry Smith   *ierr = MatGetLocalSize(*mat, &m, NULL);
147*6dd63270SBarry Smith   if (*ierr) return;
148*6dd63270SBarry Smith   *ierr = MatGetSize(*mat, NULL, &N);
149*6dd63270SBarry Smith   if (*ierr) return;
150*6dd63270SBarry Smith   *ierr = MatDenseGetLDA(*mat, &lda);
151*6dd63270SBarry Smith   if (*ierr) return;
152*6dd63270SBarry Smith   if (m != lda) { // TODO: add F90Array1dLDACreate()
153*6dd63270SBarry 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);
154*6dd63270SBarry Smith     return;
155*6dd63270SBarry Smith   }
156*6dd63270SBarry Smith   *ierr = F90Array1dCreate(fa, MPIU_SCALAR, 1, m * N, ptr PETSC_F90_2PTR_PARAM(ptrd));
157*6dd63270SBarry Smith }
158*6dd63270SBarry Smith PETSC_EXTERN void matdenserestorearray1d_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
159*6dd63270SBarry Smith {
160*6dd63270SBarry Smith   PetscScalar *fa;
161*6dd63270SBarry Smith   *ierr = F90Array1dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd));
162*6dd63270SBarry Smith   if (*ierr) return;
163*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
164*6dd63270SBarry Smith   if (*ierr) return;
165*6dd63270SBarry Smith   *ierr = MatDenseRestoreArray(*mat, &fa);
166*6dd63270SBarry Smith }
167*6dd63270SBarry Smith PETSC_EXTERN void matdensegetarrayread1d_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
168*6dd63270SBarry Smith {
169*6dd63270SBarry Smith   const PetscScalar *fa;
170*6dd63270SBarry Smith   PetscInt           m, N, lda;
171*6dd63270SBarry Smith   *ierr = MatDenseGetArrayRead(*mat, &fa);
172*6dd63270SBarry Smith   if (*ierr) return;
173*6dd63270SBarry Smith   *ierr = MatGetLocalSize(*mat, &m, NULL);
174*6dd63270SBarry Smith   if (*ierr) return;
175*6dd63270SBarry Smith   *ierr = MatGetSize(*mat, NULL, &N);
176*6dd63270SBarry Smith   if (*ierr) return;
177*6dd63270SBarry Smith   *ierr = MatDenseGetLDA(*mat, &lda);
178*6dd63270SBarry Smith   if (*ierr) return;
179*6dd63270SBarry Smith   if (m != lda) { // TODO: add F90Array1dLDACreate()
180*6dd63270SBarry 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);
181*6dd63270SBarry Smith     return;
182*6dd63270SBarry Smith   }
183*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void **)fa, MPIU_SCALAR, 1, m * N, ptr PETSC_F90_2PTR_PARAM(ptrd));
184*6dd63270SBarry Smith }
185*6dd63270SBarry Smith PETSC_EXTERN void matdenserestorearrayread1d_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
186*6dd63270SBarry Smith {
187*6dd63270SBarry Smith   const PetscScalar *fa;
188*6dd63270SBarry Smith   *ierr = F90Array1dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd));
189*6dd63270SBarry Smith   if (*ierr) return;
190*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
191*6dd63270SBarry Smith   if (*ierr) return;
192*6dd63270SBarry Smith   *ierr = MatDenseRestoreArrayRead(*mat, &fa);
193*6dd63270SBarry Smith }
194*6dd63270SBarry Smith PETSC_EXTERN void matdensegetarraywrite1d_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
195*6dd63270SBarry Smith {
196*6dd63270SBarry Smith   PetscScalar *fa;
197*6dd63270SBarry Smith   PetscInt     m, N, lda;
198*6dd63270SBarry Smith   *ierr = MatDenseGetArrayWrite(*mat, &fa);
199*6dd63270SBarry Smith   if (*ierr) return;
200*6dd63270SBarry Smith   *ierr = MatGetLocalSize(*mat, &m, NULL);
201*6dd63270SBarry Smith   if (*ierr) return;
202*6dd63270SBarry Smith   *ierr = MatGetSize(*mat, NULL, &N);
203*6dd63270SBarry Smith   if (*ierr) return;
204*6dd63270SBarry Smith   *ierr = MatDenseGetLDA(*mat, &lda);
205*6dd63270SBarry Smith   if (*ierr) return;
206*6dd63270SBarry Smith   if (m != lda) { // TODO: add F90Array1dLDACreate()
207*6dd63270SBarry 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);
208*6dd63270SBarry Smith     return;
209*6dd63270SBarry Smith   }
210*6dd63270SBarry Smith   *ierr = F90Array1dCreate(fa, MPIU_SCALAR, 1, m * N, ptr PETSC_F90_2PTR_PARAM(ptrd));
211*6dd63270SBarry Smith }
212*6dd63270SBarry Smith PETSC_EXTERN void matdenserestorearraywrite1d_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
213*6dd63270SBarry Smith {
214*6dd63270SBarry Smith   PetscScalar *fa;
215*6dd63270SBarry Smith   *ierr = F90Array1dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd));
216*6dd63270SBarry Smith   if (*ierr) return;
217*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
218*6dd63270SBarry Smith   if (*ierr) return;
219*6dd63270SBarry Smith   *ierr = MatDenseRestoreArrayWrite(*mat, &fa);
220*6dd63270SBarry Smith }
221*6dd63270SBarry Smith PETSC_EXTERN void matdensegetcolumn_(Mat *mat, PetscInt *col, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
222*6dd63270SBarry Smith {
223*6dd63270SBarry Smith   PetscScalar *fa;
224*6dd63270SBarry Smith   PetscInt     m;
225*6dd63270SBarry Smith   *ierr = MatDenseGetColumn(*mat, *col, &fa);
226*6dd63270SBarry Smith   if (*ierr) return;
227*6dd63270SBarry Smith   *ierr = MatGetLocalSize(*mat, &m, NULL);
228*6dd63270SBarry Smith   if (*ierr) return;
229*6dd63270SBarry Smith   *ierr = F90Array1dCreate(fa, MPIU_SCALAR, 1, m, ptr PETSC_F90_2PTR_PARAM(ptrd));
230*6dd63270SBarry Smith }
231*6dd63270SBarry Smith PETSC_EXTERN void matdenserestorecolumn_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
232*6dd63270SBarry Smith {
233*6dd63270SBarry Smith   PetscScalar *fa;
234*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
235*6dd63270SBarry Smith   if (*ierr) return;
236*6dd63270SBarry Smith   *ierr = MatDenseRestoreColumn(*mat, &fa);
237*6dd63270SBarry Smith }
238*6dd63270SBarry Smith PETSC_EXTERN void matseqaijgetarray_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
239*6dd63270SBarry Smith {
240*6dd63270SBarry Smith   PetscScalar *fa;
241*6dd63270SBarry Smith   PetscInt     m, n;
242*6dd63270SBarry Smith   *ierr = MatSeqAIJGetArray(*mat, &fa);
243*6dd63270SBarry Smith   if (*ierr) return;
244*6dd63270SBarry Smith   *ierr = MatGetLocalSize(*mat, &m, &n);
245*6dd63270SBarry Smith   if (*ierr) return;
246*6dd63270SBarry Smith   *ierr = F90Array1dCreate(fa, MPIU_SCALAR, 1, m * n, ptr PETSC_F90_2PTR_PARAM(ptrd));
247*6dd63270SBarry Smith }
248*6dd63270SBarry Smith PETSC_EXTERN void matseqaijrestorearray_(Mat *mat, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
249*6dd63270SBarry Smith {
250*6dd63270SBarry Smith   PetscScalar *fa;
251*6dd63270SBarry Smith   *ierr = F90Array1dAccess(ptr, MPIU_SCALAR, (void **)&fa PETSC_F90_2PTR_PARAM(ptrd));
252*6dd63270SBarry Smith   if (*ierr) return;
253*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_SCALAR PETSC_F90_2PTR_PARAM(ptrd));
254*6dd63270SBarry Smith   if (*ierr) return;
255*6dd63270SBarry Smith   *ierr = MatSeqAIJRestoreArray(*mat, &fa);
256*6dd63270SBarry Smith }
257*6dd63270SBarry 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))
258*6dd63270SBarry Smith {
259*6dd63270SBarry Smith   const PetscInt *IA, *JA;
260*6dd63270SBarry Smith   *ierr = MatGetRowIJ(*B, *shift, *sym, *blockcompressed, n, &IA, &JA, done);
261*6dd63270SBarry Smith   if (*ierr) return;
262*6dd63270SBarry Smith   if (!*done) return;
263*6dd63270SBarry Smith   *ierr = F90Array1dCreate((PetscInt *)IA, MPIU_INT, 1, *n + 1, ia PETSC_F90_2PTR_PARAM(iad));
264*6dd63270SBarry Smith   *ierr = F90Array1dCreate((PetscInt *)JA, MPIU_INT, 1, IA[*n], ja PETSC_F90_2PTR_PARAM(jad));
265*6dd63270SBarry Smith }
266*6dd63270SBarry 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))
267*6dd63270SBarry Smith {
268*6dd63270SBarry Smith   const PetscInt *IA, *JA;
269*6dd63270SBarry Smith   *ierr = F90Array1dAccess(ia, MPIU_INT, (void **)&IA PETSC_F90_2PTR_PARAM(iad));
270*6dd63270SBarry Smith   if (*ierr) return;
271*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ia, MPIU_INT PETSC_F90_2PTR_PARAM(iad));
272*6dd63270SBarry Smith   if (*ierr) return;
273*6dd63270SBarry Smith   *ierr = F90Array1dAccess(ja, MPIU_INT, (void **)&JA PETSC_F90_2PTR_PARAM(jad));
274*6dd63270SBarry Smith   if (*ierr) return;
275*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ja, MPIU_INT PETSC_F90_2PTR_PARAM(jad));
276*6dd63270SBarry Smith   if (*ierr) return;
277*6dd63270SBarry Smith   *ierr = MatRestoreRowIJ(*B, *shift, *sym, *blockcompressed, n, &IA, &JA, done);
278*6dd63270SBarry Smith }
279*6dd63270SBarry Smith PETSC_EXTERN void matmpiaijgetseqaij_(Mat *mat, Mat *A, Mat *B, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
280*6dd63270SBarry Smith {
281*6dd63270SBarry Smith   const PetscInt *fa;
282*6dd63270SBarry Smith   PetscInt        n;
283*6dd63270SBarry Smith   *ierr = MatMPIAIJGetSeqAIJ(*mat, A, B, &fa);
284*6dd63270SBarry Smith   if (*ierr) return;
285*6dd63270SBarry Smith   *ierr = MatGetLocalSize(*B, NULL, &n);
286*6dd63270SBarry Smith   if (*ierr) return;
287*6dd63270SBarry Smith   *ierr = F90Array1dCreate((void *)fa, MPIU_INT, 1, n, ptr PETSC_F90_2PTR_PARAM(ptrd));
288*6dd63270SBarry Smith }
289*6dd63270SBarry Smith PETSC_EXTERN void matmpiaijrestoreseqaij_(Mat *mat, Mat *A, Mat *B, F90Array1d *ptr, int *ierr PETSC_F90_2PTR_PROTO(ptrd))
290*6dd63270SBarry Smith {
291*6dd63270SBarry Smith   *ierr = F90Array1dDestroy(ptr, MPIU_INT PETSC_F90_2PTR_PARAM(ptrd));
292*6dd63270SBarry Smith   if (*ierr) return;
293*6dd63270SBarry Smith }
294