xref: /petsc/src/mat/interface/ftn-custom/zmatrixf.c (revision a58c3bc3391eee32bc3fd94ac7edeea38fe57aae)
1 #include "zpetsc.h"
2 #include "petscmat.h"
3 
4 #if defined(PETSC_HAVE_FORTRAN_CAPS)
5 #define matgetrowij_                     MATGETROWIJ
6 #define matrestorerowij_                 MATRESTOREROWIJ
7 #define matgetrow_                       MATGETROW
8 #define matrestorerow_                   MATRESTOREROW
9 #define matview_                         MATVIEW
10 #define matgetarray_                     MATGETARRAY
11 #define matrestorearray_                 MATRESTOREARRAY
12 #define mattranspose_                    MATTRANSPOSE
13 #define matconvert_                      MATCONVERT
14 #define matgetsubmatrices_               MATGETSUBMATRICES
15 #define matzerorows_                     MATZEROROWS
16 #define matzerorowsis_                   MATZEROROWSIS
17 #define matzerorowslocal_                MATZEROROWSLOCAL
18 #define matzerorowslocalis_              MATZEROROWSLOCALIS
19 #define matsetoptionsprefix_             MATSETOPTIONSPREFIX
20 #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
21 #define matgetrowij_                     matgetrowij
22 #define matrestorerowij_                 matrestorerowij
23 #define matgetrow_                       matgetrow
24 #define matrestorerow_                   matrestorerow
25 #define matview_                         matview
26 #define matgetarray_                     matgetarray
27 #define matrestorearray_                 matrestorearray
28 #define mattranspose_                    mattranspose
29 #define matconvert_                      matconvert
30 #define matgetsubmatrices_               matgetsubmatrices
31 #define matzerorows_                     matzerorows
32 #define matzerorowsis_                   matzerorowsis
33 #define matzerorowslocal_                matzerorowslocal
34 #define matzerorowslocalis_              matzerorowslocalis
35 #define matsetoptionsprefix_             matsetoptionsprefix
36 #endif
37 
38 EXTERN_C_BEGIN
39 
40 void PETSC_STDCALL matgetrowij_(Mat *B,PetscInt *shift,PetscTruth *sym,PetscTruth *blockcompressed,PetscInt *n,PetscInt *ia,size_t *iia,
41                                 PetscInt *ja,size_t *jja,PetscTruth *done,PetscErrorCode *ierr)
42 {
43   PetscInt *IA,*JA;
44   *ierr = MatGetRowIJ(*B,*shift,*sym,*blockcompressed,n,&IA,&JA,done);if (*ierr) return;
45   *iia  = PetscIntAddressToFortran(ia,IA);
46   *jja  = PetscIntAddressToFortran(ja,JA);
47 }
48 
49 void PETSC_STDCALL matrestorerowij_(Mat *B,PetscInt *shift,PetscTruth *sym,PetscTruth *blockcompressed, PetscInt *n,PetscInt *ia,size_t *iia,
50                                     PetscInt *ja,size_t *jja,PetscTruth *done,PetscErrorCode *ierr)
51 {
52   PetscInt *IA = PetscIntAddressFromFortran(ia,*iia),*JA = PetscIntAddressFromFortran(ja,*jja);
53   *ierr = MatRestoreRowIJ(*B,*shift,*sym,*blockcompressed,n,&IA,&JA,done);
54 }
55 
56 /*
57    This is a poor way of storing the column and value pointers
58   generated by MatGetRow() to be returned with MatRestoreRow()
59   but there is not natural,good place else to store them. Hence
60   Fortran programmers can only have one outstanding MatGetRows()
61   at a time.
62 */
63 static PetscErrorCode    matgetrowactive = 0;
64 static const PetscInt    *my_ocols = 0;
65 static const PetscScalar *my_ovals = 0;
66 
67 void PETSC_STDCALL matgetrow_(Mat *mat,PetscInt *row,PetscInt *ncols,PetscInt *cols,PetscScalar *vals,PetscErrorCode *ierr)
68 {
69   const PetscInt    **oocols = &my_ocols;
70   const PetscScalar **oovals = &my_ovals;
71 
72   if (matgetrowactive) {
73      PetscError(__LINE__,"MatGetRow_Fortran",__FILE__,__SDIR__,1,0,
74                "Cannot have two MatGetRow() active simultaneously\n\
75                call MatRestoreRow() before calling MatGetRow() a second time");
76      *ierr = 1;
77      return;
78   }
79 
80   CHKFORTRANNULLINTEGER(cols); if (!cols) oocols = PETSC_NULL;
81   CHKFORTRANNULLSCALAR(vals);  if (!vals) oovals = PETSC_NULL;
82 
83   *ierr = MatGetRow(*mat,*row,ncols,oocols,oovals);
84   if (*ierr) return;
85 
86   if (oocols) { *ierr = PetscMemcpy(cols,my_ocols,(*ncols)*sizeof(PetscInt)); if (*ierr) return;}
87   if (oovals) { *ierr = PetscMemcpy(vals,my_ovals,(*ncols)*sizeof(PetscScalar)); if (*ierr) return; }
88   matgetrowactive = 1;
89 }
90 
91 void PETSC_STDCALL matrestorerow_(Mat *mat,PetscInt *row,PetscInt *ncols,PetscInt *cols,PetscScalar *vals,PetscErrorCode *ierr)
92 {
93   const PetscInt         **oocols = &my_ocols;
94   const PetscScalar **oovals = &my_ovals;
95   if (!matgetrowactive) {
96      PetscError(__LINE__,"MatRestoreRow_Fortran",__FILE__,__SDIR__,1,0,
97                "Must call MatGetRow() first");
98      *ierr = 1;
99      return;
100   }
101   CHKFORTRANNULLINTEGER(cols); if (!cols) oocols = PETSC_NULL;
102   CHKFORTRANNULLSCALAR(vals);  if (!vals) oovals = PETSC_NULL;
103 
104   *ierr = MatRestoreRow(*mat,*row,ncols,oocols,oovals);
105   matgetrowactive = 0;
106 }
107 
108 void PETSC_STDCALL matview_(Mat *mat,PetscViewer *vin,PetscErrorCode *ierr)
109 {
110   PetscViewer v;
111   PetscPatchDefaultViewers_Fortran(vin,v);
112   *ierr = MatView(*mat,v);
113 }
114 
115 void PETSC_STDCALL matgetarray_(Mat *mat,PetscScalar *fa,size_t *ia,PetscErrorCode *ierr)
116 {
117   PetscScalar *mm;
118   PetscInt    m,n;
119 
120   *ierr = MatGetArray(*mat,&mm); if (*ierr) return;
121   *ierr = MatGetSize(*mat,&m,&n);  if (*ierr) return;
122   *ierr = PetscScalarAddressToFortran((PetscObject)*mat,1,fa,mm,m*n,ia); if (*ierr) return;
123 }
124 
125 void PETSC_STDCALL matrestorearray_(Mat *mat,PetscScalar *fa,size_t *ia,PetscErrorCode *ierr)
126 {
127   PetscScalar          *lx;
128   PetscInt                  m,n;
129 
130   *ierr = MatGetSize(*mat,&m,&n); if (*ierr) return;
131   *ierr = PetscScalarAddressFromFortran((PetscObject)*mat,fa,*ia,m*n,&lx);if (*ierr) return;
132   *ierr = MatRestoreArray(*mat,&lx);if (*ierr) return;
133 }
134 
135 void PETSC_STDCALL mattranspose_(Mat *mat,Mat *B,PetscErrorCode *ierr)
136 {
137   CHKFORTRANNULLOBJECT(B);
138   *ierr = MatTranspose(*mat,B);
139 }
140 
141 void PETSC_STDCALL matconvert_(Mat *mat,CHAR outtype PETSC_MIXED_LEN(len),MatReuse *reuse,Mat *M,PetscErrorCode *ierr PETSC_END_LEN(len))
142 {
143   char *t;
144   FIXCHAR(outtype,len,t);
145   *ierr = MatConvert(*mat,t,*reuse,M);
146   FREECHAR(outtype,t);
147 }
148 
149 /*
150     MatGetSubmatrices() is slightly different from C since the
151     Fortran provides the array to hold the submatrix objects,while in C that
152     array is allocated by the MatGetSubmatrices()
153 */
154 void PETSC_STDCALL matgetsubmatrices_(Mat *mat,PetscInt *n,IS *isrow,IS *iscol,MatReuse *scall,Mat *smat,PetscErrorCode *ierr)
155 {
156   Mat *lsmat;
157   PetscInt i;
158 
159   if (*scall == MAT_INITIAL_MATRIX) {
160     *ierr = MatGetSubMatrices(*mat,*n,isrow,iscol,*scall,&lsmat);
161     for (i=0; i<*n; i++) {
162       smat[i] = lsmat[i];
163     }
164     *ierr = PetscFree(lsmat);
165   } else {
166     *ierr = MatGetSubMatrices(*mat,*n,isrow,iscol,*scall,&smat);
167   }
168 }
169 
170 void PETSC_STDCALL matzerorows_(Mat *mat,PetscInt *numRows,PetscInt *rows,PetscScalar *diag,PetscErrorCode *ierr)
171 {
172   *ierr = MatZeroRows(*mat,*numRows,rows,*diag);
173 }
174 
175 void PETSC_STDCALL matzerorowsis_(Mat *mat,IS *is,PetscScalar *diag,PetscErrorCode *ierr)
176 {
177   *ierr = MatZeroRowsIS(*mat,*is,*diag);
178 }
179 
180 void PETSC_STDCALL matzerorowslocal_(Mat *mat,PetscInt *numRows,PetscInt *rows,PetscScalar *diag,PetscErrorCode *ierr)
181 {
182   *ierr = MatZeroRowsLocal(*mat,*numRows,rows,*diag);
183 }
184 
185 void PETSC_STDCALL matzerorowslocalis_(Mat *mat,IS *is,PetscScalar *diag,PetscErrorCode *ierr)
186 {
187   *ierr = MatZeroRowsLocalIS(*mat,*is,*diag);
188 }
189 
190 
191 void PETSC_STDCALL matsetoptionsprefix_(Mat *mat,CHAR prefix PETSC_MIXED_LEN(len),
192                                         PetscErrorCode *ierr PETSC_END_LEN(len))
193 {
194   char *t;
195 
196   FIXCHAR(prefix,len,t);
197   *ierr = MatSetOptionsPrefix(*mat,t);
198   FREECHAR(prefix,t);
199 }
200 
201 
202 EXTERN_C_END
203