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