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