xref: /petsc/src/mat/interface/ftn-custom/zmatrixf.c (revision 22688b47e4cfe7550dca274eea02c9b2d86fd385)
1ce0a2cd1SBarry Smith #include "private/fortranimpl.h"
2f4e70085SSatish Balay #include "petscmat.h"
3f4e70085SSatish Balay 
4f4e70085SSatish Balay #if defined(PETSC_HAVE_FORTRAN_CAPS)
57d6bfa3bSBarry Smith #define matdestroymatrices_              MATDESTROYMATRICES
65dffd610SBarry Smith #define matgetfactor_                    MATGETFACTOR
735bd34faSBarry Smith #define matfactorgetsolverpackage_       MATFACTORGETSOLVERPACKAGE
8f4e70085SSatish Balay #define matgetrowij_                     MATGETROWIJ
9f4e70085SSatish Balay #define matrestorerowij_                 MATRESTOREROWIJ
10f4e70085SSatish Balay #define matgetrow_                       MATGETROW
11f4e70085SSatish Balay #define matrestorerow_                   MATRESTOREROW
12f4e70085SSatish Balay #define matview_                         MATVIEW
13f4e70085SSatish Balay #define matgetarray_                     MATGETARRAY
14f4e70085SSatish Balay #define matrestorearray_                 MATRESTOREARRAY
15f4e70085SSatish Balay #define matconvert_                      MATCONVERT
16f4e70085SSatish Balay #define matgetsubmatrices_               MATGETSUBMATRICES
17f4e70085SSatish Balay #define matzerorows_                     MATZEROROWS
18f4e70085SSatish Balay #define matzerorowsis_                   MATZEROROWSIS
19f4e70085SSatish Balay #define matzerorowslocal_                MATZEROROWSLOCAL
20f4e70085SSatish Balay #define matzerorowslocalis_              MATZEROROWSLOCALIS
211eea217eSSatish Balay #define matsetoptionsprefix_             MATSETOPTIONSPREFIX
227c54600cSBarry Smith #define matgetvecs_                      MATGETVECS
235ba43861SSatish Balay #define matnullspaceremove_              MATNULLSPACEREMOV
245ba43861SSatish Balay #define matgetinfo_                      MATGETINFO
25*22688b47SSatish Balay #define matlufactor_                     MATLUFACTOR
26*22688b47SSatish Balay #define matilufactor_                    MATILUFACTOR
27*22688b47SSatish Balay #define matlufactorsymbolic_             MATLUFACTORSYMBOLIC
28*22688b47SSatish Balay #define matlufactornumeric_              MATLUFACTORNUMERIC
29*22688b47SSatish Balay #define matcholeskyfactor_               MATCHOLESKYFACTOR
30*22688b47SSatish Balay #define matcholeskyfactorsymbolic_       MATCHOLESKYFACTORSYMBOLIC
31*22688b47SSatish Balay #define matcholeskyfactornumeric_        MATCHOLESKYFACTORNUMERIC
32*22688b47SSatish Balay #define matilufactorsymbolic_            MATILUFACTORSYMBOLIC
33*22688b47SSatish Balay #define maticcfactorsymbolic_            MATICCFACTORSYMBOLIC
34*22688b47SSatish Balay #define maticcfactor_                    MATICCFACTOR
35*22688b47SSatish Balay #define matfactorinfoinitialize_         MATFACTORINFOINITIALIZE
36f4e70085SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
377d6bfa3bSBarry Smith #define matdestroymatrices_              matdestroymatrices_
385dffd610SBarry Smith #define matgetfactor_                    matgetfactor
3935bd34faSBarry Smith #define matfactorgetsolverpackage_       matfactorgetsolverpackage
407c54600cSBarry Smith #define matgetvecs_                      matgetvecs
41f4e70085SSatish Balay #define matgetrowij_                     matgetrowij
42f4e70085SSatish Balay #define matrestorerowij_                 matrestorerowij
43f4e70085SSatish Balay #define matgetrow_                       matgetrow
44f4e70085SSatish Balay #define matrestorerow_                   matrestorerow
45f4e70085SSatish Balay #define matview_                         matview
46f4e70085SSatish Balay #define matgetarray_                     matgetarray
47f4e70085SSatish Balay #define matrestorearray_                 matrestorearray
48f4e70085SSatish Balay #define matconvert_                      matconvert
49f4e70085SSatish Balay #define matgetsubmatrices_               matgetsubmatrices
50f4e70085SSatish Balay #define matzerorows_                     matzerorows
51f4e70085SSatish Balay #define matzerorowsis_                   matzerorowsis
52f4e70085SSatish Balay #define matzerorowslocal_                matzerorowslocal
53f4e70085SSatish Balay #define matzerorowslocalis_              matzerorowslocalis
541eea217eSSatish Balay #define matsetoptionsprefix_             matsetoptionsprefix
55812c3f48SMatthew Knepley #define matnullspaceremove_              matnullspaceremove
565ba43861SSatish Balay #define matgetinfo_                      matgetinfo
57*22688b47SSatish Balay #define matlufactor_                     matlufactor
58*22688b47SSatish Balay #define matilufactor_                    matilufactor
59*22688b47SSatish Balay #define matlufactorsymbolic_             matlufactorsymbolic
60*22688b47SSatish Balay #define matlufactornumeric_              matlufactornumeric
61*22688b47SSatish Balay #define matcholeskyfactor_               matcholeskyfactor
62*22688b47SSatish Balay #define matcholeskyfactorsymbolic_       matcholeskyfactorsymbolic
63*22688b47SSatish Balay #define matcholeskyfactornumeric_        matcholeskyfactornumeric
64*22688b47SSatish Balay #define matilufactorsymbolic_            matilufactorsymbolic
65*22688b47SSatish Balay #define maticcfactorsymbolic_            maticcfactorsymbolic
66*22688b47SSatish Balay #define maticcfactor_                    maticcfactor
67*22688b47SSatish Balay #define matfactorinfoinitialize_         matfactorinfoinitialize
68f4e70085SSatish Balay #endif
69f4e70085SSatish Balay 
70f4e70085SSatish Balay EXTERN_C_BEGIN
71f4e70085SSatish Balay 
727c54600cSBarry Smith void PETSC_STDCALL   matgetvecs_(Mat *mat,Vec *right,Vec *left, int *ierr )
737c54600cSBarry Smith {
747c54600cSBarry Smith   CHKFORTRANNULLOBJECT(right);
757c54600cSBarry Smith   CHKFORTRANNULLOBJECT(left);
767c54600cSBarry Smith   *ierr = MatGetVecs(*mat,right,left);
777c54600cSBarry Smith }
787c54600cSBarry Smith 
798f7157efSSatish Balay void PETSC_STDCALL matgetrowij_(Mat *B,PetscInt *shift,PetscTruth *sym,PetscTruth *blockcompressed,PetscInt *n,PetscInt *ia,size_t *iia,
808f7157efSSatish Balay                                 PetscInt *ja,size_t *jja,PetscTruth *done,PetscErrorCode *ierr)
81f4e70085SSatish Balay {
82f4e70085SSatish Balay   PetscInt *IA,*JA;
838f7157efSSatish Balay   *ierr = MatGetRowIJ(*B,*shift,*sym,*blockcompressed,n,&IA,&JA,done);if (*ierr) return;
84f4e70085SSatish Balay   *iia  = PetscIntAddressToFortran(ia,IA);
85f4e70085SSatish Balay   *jja  = PetscIntAddressToFortran(ja,JA);
86f4e70085SSatish Balay }
87f4e70085SSatish Balay 
888f7157efSSatish Balay void PETSC_STDCALL matrestorerowij_(Mat *B,PetscInt *shift,PetscTruth *sym,PetscTruth *blockcompressed, PetscInt *n,PetscInt *ia,size_t *iia,
898f7157efSSatish Balay                                     PetscInt *ja,size_t *jja,PetscTruth *done,PetscErrorCode *ierr)
90f4e70085SSatish Balay {
91f4e70085SSatish Balay   PetscInt *IA = PetscIntAddressFromFortran(ia,*iia),*JA = PetscIntAddressFromFortran(ja,*jja);
928f7157efSSatish Balay   *ierr = MatRestoreRowIJ(*B,*shift,*sym,*blockcompressed,n,&IA,&JA,done);
93f4e70085SSatish Balay }
94f4e70085SSatish Balay 
95f4e70085SSatish Balay /*
96f4e70085SSatish Balay    This is a poor way of storing the column and value pointers
97f4e70085SSatish Balay   generated by MatGetRow() to be returned with MatRestoreRow()
98f4e70085SSatish Balay   but there is not natural,good place else to store them. Hence
99f4e70085SSatish Balay   Fortran programmers can only have one outstanding MatGetRows()
100f4e70085SSatish Balay   at a time.
101f4e70085SSatish Balay */
102f4e70085SSatish Balay static PetscErrorCode    matgetrowactive = 0;
103f4e70085SSatish Balay static const PetscInt    *my_ocols = 0;
104f4e70085SSatish Balay static const PetscScalar *my_ovals = 0;
105f4e70085SSatish Balay 
106f4e70085SSatish Balay void PETSC_STDCALL matgetrow_(Mat *mat,PetscInt *row,PetscInt *ncols,PetscInt *cols,PetscScalar *vals,PetscErrorCode *ierr)
107f4e70085SSatish Balay {
108f4e70085SSatish Balay   const PetscInt    **oocols = &my_ocols;
109f4e70085SSatish Balay   const PetscScalar **oovals = &my_ovals;
110f4e70085SSatish Balay 
111f4e70085SSatish Balay   if (matgetrowactive) {
112f4e70085SSatish Balay      PetscError(__LINE__,"MatGetRow_Fortran",__FILE__,__SDIR__,1,0,
113f4e70085SSatish Balay                "Cannot have two MatGetRow() active simultaneously\n\
114f4e70085SSatish Balay                call MatRestoreRow() before calling MatGetRow() a second time");
115f4e70085SSatish Balay      *ierr = 1;
116f4e70085SSatish Balay      return;
117f4e70085SSatish Balay   }
118f4e70085SSatish Balay 
119f4e70085SSatish Balay   CHKFORTRANNULLINTEGER(cols); if (!cols) oocols = PETSC_NULL;
120f4e70085SSatish Balay   CHKFORTRANNULLSCALAR(vals);  if (!vals) oovals = PETSC_NULL;
121f4e70085SSatish Balay 
122f4e70085SSatish Balay   *ierr = MatGetRow(*mat,*row,ncols,oocols,oovals);
123f4e70085SSatish Balay   if (*ierr) return;
124f4e70085SSatish Balay 
125f4e70085SSatish Balay   if (oocols) { *ierr = PetscMemcpy(cols,my_ocols,(*ncols)*sizeof(PetscInt)); if (*ierr) return;}
126f4e70085SSatish Balay   if (oovals) { *ierr = PetscMemcpy(vals,my_ovals,(*ncols)*sizeof(PetscScalar)); if (*ierr) return; }
127f4e70085SSatish Balay   matgetrowactive = 1;
128f4e70085SSatish Balay }
129f4e70085SSatish Balay 
130f4e70085SSatish Balay void PETSC_STDCALL matrestorerow_(Mat *mat,PetscInt *row,PetscInt *ncols,PetscInt *cols,PetscScalar *vals,PetscErrorCode *ierr)
131f4e70085SSatish Balay {
132f4e70085SSatish Balay   const PetscInt         **oocols = &my_ocols;
133f4e70085SSatish Balay   const PetscScalar **oovals = &my_ovals;
134f4e70085SSatish Balay   if (!matgetrowactive) {
135f4e70085SSatish Balay      PetscError(__LINE__,"MatRestoreRow_Fortran",__FILE__,__SDIR__,1,0,
136f4e70085SSatish Balay                "Must call MatGetRow() first");
137f4e70085SSatish Balay      *ierr = 1;
138f4e70085SSatish Balay      return;
139f4e70085SSatish Balay   }
140f4e70085SSatish Balay   CHKFORTRANNULLINTEGER(cols); if (!cols) oocols = PETSC_NULL;
141f4e70085SSatish Balay   CHKFORTRANNULLSCALAR(vals);  if (!vals) oovals = PETSC_NULL;
142f4e70085SSatish Balay 
143f4e70085SSatish Balay   *ierr = MatRestoreRow(*mat,*row,ncols,oocols,oovals);
144f4e70085SSatish Balay   matgetrowactive = 0;
145f4e70085SSatish Balay }
146f4e70085SSatish Balay 
147f4e70085SSatish Balay void PETSC_STDCALL matview_(Mat *mat,PetscViewer *vin,PetscErrorCode *ierr)
148f4e70085SSatish Balay {
149f4e70085SSatish Balay   PetscViewer v;
150f4e70085SSatish Balay   PetscPatchDefaultViewers_Fortran(vin,v);
151f4e70085SSatish Balay   *ierr = MatView(*mat,v);
152f4e70085SSatish Balay }
153f4e70085SSatish Balay 
154f4e70085SSatish Balay void PETSC_STDCALL matgetarray_(Mat *mat,PetscScalar *fa,size_t *ia,PetscErrorCode *ierr)
155f4e70085SSatish Balay {
156f4e70085SSatish Balay   PetscScalar *mm;
157f4e70085SSatish Balay   PetscInt    m,n;
158f4e70085SSatish Balay 
159f4e70085SSatish Balay   *ierr = MatGetArray(*mat,&mm); if (*ierr) return;
160f4e70085SSatish Balay   *ierr = MatGetSize(*mat,&m,&n);  if (*ierr) return;
161f91d1997SBarry Smith   *ierr = PetscScalarAddressToFortran((PetscObject)*mat,1,fa,mm,m*n,ia); if (*ierr) return;
162f4e70085SSatish Balay }
163f4e70085SSatish Balay 
164f4e70085SSatish Balay void PETSC_STDCALL matrestorearray_(Mat *mat,PetscScalar *fa,size_t *ia,PetscErrorCode *ierr)
165f4e70085SSatish Balay {
166f4e70085SSatish Balay   PetscScalar          *lx;
167f4e70085SSatish Balay   PetscInt                  m,n;
168f4e70085SSatish Balay 
169f4e70085SSatish Balay   *ierr = MatGetSize(*mat,&m,&n); if (*ierr) return;
170f4e70085SSatish Balay   *ierr = PetscScalarAddressFromFortran((PetscObject)*mat,fa,*ia,m*n,&lx);if (*ierr) return;
171f4e70085SSatish Balay   *ierr = MatRestoreArray(*mat,&lx);if (*ierr) return;
172f4e70085SSatish Balay }
173f4e70085SSatish Balay 
17435bd34faSBarry Smith void PETSC_STDCALL matfactorgetsolverpackage_(Mat *mat,CHAR name PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
17535bd34faSBarry Smith {
17635bd34faSBarry Smith   const char *tname;
17735bd34faSBarry Smith 
17835bd34faSBarry Smith   *ierr = MatFactorGetSolverPackage(*mat,&tname);if (*ierr) return;
17935bd34faSBarry Smith   if (name != PETSC_NULL_CHARACTER_Fortran) {
18035bd34faSBarry Smith     *ierr = PetscStrncpy(name,tname,len);if (*ierr) return;
18135bd34faSBarry Smith   }
18235bd34faSBarry Smith   FIXRETURNCHAR(PETSC_TRUE,name,len);
18335bd34faSBarry Smith }
18435bd34faSBarry Smith 
1855dffd610SBarry Smith void PETSC_STDCALL matgetfactor_(Mat *mat,CHAR outtype PETSC_MIXED_LEN(len),MatFactorType ftype,Mat *M,PetscErrorCode *ierr PETSC_END_LEN(len))
1865dffd610SBarry Smith {
1875dffd610SBarry Smith   char *t;
1885dffd610SBarry Smith   FIXCHAR(outtype,len,t);
189c911e420SBarry Smith   *ierr = MatGetFactor(*mat,t,ftype,M);
1905dffd610SBarry Smith   FREECHAR(outtype,t);
1915dffd610SBarry Smith }
1925dffd610SBarry Smith 
193f4e70085SSatish Balay void PETSC_STDCALL matconvert_(Mat *mat,CHAR outtype PETSC_MIXED_LEN(len),MatReuse *reuse,Mat *M,PetscErrorCode *ierr PETSC_END_LEN(len))
194f4e70085SSatish Balay {
195f4e70085SSatish Balay   char *t;
196f4e70085SSatish Balay   FIXCHAR(outtype,len,t);
197f4e70085SSatish Balay   *ierr = MatConvert(*mat,t,*reuse,M);
198f4e70085SSatish Balay   FREECHAR(outtype,t);
199f4e70085SSatish Balay }
200f4e70085SSatish Balay 
201f4e70085SSatish Balay /*
202f4e70085SSatish Balay     MatGetSubmatrices() is slightly different from C since the
203f4e70085SSatish Balay     Fortran provides the array to hold the submatrix objects,while in C that
204f4e70085SSatish Balay     array is allocated by the MatGetSubmatrices()
205f4e70085SSatish Balay */
206f4e70085SSatish Balay void PETSC_STDCALL matgetsubmatrices_(Mat *mat,PetscInt *n,IS *isrow,IS *iscol,MatReuse *scall,Mat *smat,PetscErrorCode *ierr)
207f4e70085SSatish Balay {
208f4e70085SSatish Balay   Mat *lsmat;
209f4e70085SSatish Balay   PetscInt i;
210f4e70085SSatish Balay 
211f4e70085SSatish Balay   if (*scall == MAT_INITIAL_MATRIX) {
212f4e70085SSatish Balay     *ierr = MatGetSubMatrices(*mat,*n,isrow,iscol,*scall,&lsmat);
213f4e70085SSatish Balay     for (i=0; i<*n; i++) {
214f4e70085SSatish Balay       smat[i] = lsmat[i];
215f4e70085SSatish Balay     }
216f4e70085SSatish Balay     *ierr = PetscFree(lsmat);
217f4e70085SSatish Balay   } else {
218f4e70085SSatish Balay     *ierr = MatGetSubMatrices(*mat,*n,isrow,iscol,*scall,&smat);
219f4e70085SSatish Balay   }
220f4e70085SSatish Balay }
221f4e70085SSatish Balay 
2227d6bfa3bSBarry Smith /*
2237d6bfa3bSBarry Smith     MatDestroyMatrices() is slightly different from C since the
2247d6bfa3bSBarry Smith     Fortran provides the array to hold the submatrix objects,while in C that
2257d6bfa3bSBarry Smith     array is allocated by the MatGetSubmatrices()
2267d6bfa3bSBarry Smith */
2277d6bfa3bSBarry Smith void PETSC_STDCALL matdestroymatrices_(Mat *mat,PetscInt *n,Mat *smat,PetscErrorCode *ierr)
2287d6bfa3bSBarry Smith {
2297d6bfa3bSBarry Smith   PetscInt i;
2307d6bfa3bSBarry Smith 
2317d6bfa3bSBarry Smith   for (i=0; i<*n; i++) {
2327d6bfa3bSBarry Smith     *ierr = MatDestroy(smat[i]);if (*ierr) return;
2337d6bfa3bSBarry Smith   }
2347d6bfa3bSBarry Smith }
2357d6bfa3bSBarry Smith 
236f4e70085SSatish Balay void PETSC_STDCALL matzerorows_(Mat *mat,PetscInt *numRows,PetscInt *rows,PetscScalar *diag,PetscErrorCode *ierr)
237f4e70085SSatish Balay {
238f4e70085SSatish Balay   *ierr = MatZeroRows(*mat,*numRows,rows,*diag);
239f4e70085SSatish Balay }
240f4e70085SSatish Balay 
241f4e70085SSatish Balay void PETSC_STDCALL matzerorowsis_(Mat *mat,IS *is,PetscScalar *diag,PetscErrorCode *ierr)
242f4e70085SSatish Balay {
243f4e70085SSatish Balay   *ierr = MatZeroRowsIS(*mat,*is,*diag);
244f4e70085SSatish Balay }
245f4e70085SSatish Balay 
246f4e70085SSatish Balay void PETSC_STDCALL matzerorowslocal_(Mat *mat,PetscInt *numRows,PetscInt *rows,PetscScalar *diag,PetscErrorCode *ierr)
247f4e70085SSatish Balay {
248f4e70085SSatish Balay   *ierr = MatZeroRowsLocal(*mat,*numRows,rows,*diag);
249f4e70085SSatish Balay }
250f4e70085SSatish Balay 
251f4e70085SSatish Balay void PETSC_STDCALL matzerorowslocalis_(Mat *mat,IS *is,PetscScalar *diag,PetscErrorCode *ierr)
252f4e70085SSatish Balay {
253f4e70085SSatish Balay   *ierr = MatZeroRowsLocalIS(*mat,*is,*diag);
254f4e70085SSatish Balay }
255f4e70085SSatish Balay 
2561eea217eSSatish Balay 
2571eea217eSSatish Balay void PETSC_STDCALL matsetoptionsprefix_(Mat *mat,CHAR prefix PETSC_MIXED_LEN(len),
2581eea217eSSatish Balay                                         PetscErrorCode *ierr PETSC_END_LEN(len))
2591eea217eSSatish Balay {
2601eea217eSSatish Balay   char *t;
2611eea217eSSatish Balay 
2621eea217eSSatish Balay   FIXCHAR(prefix,len,t);
2631eea217eSSatish Balay   *ierr = MatSetOptionsPrefix(*mat,t);
2641eea217eSSatish Balay   FREECHAR(prefix,t);
2651eea217eSSatish Balay }
2661eea217eSSatish Balay 
267812c3f48SMatthew Knepley void PETSC_STDCALL matnullspaceremove_(MatNullSpace *sp,Vec *vec,Vec *out,PetscErrorCode *ierr)
268812c3f48SMatthew Knepley {
269812c3f48SMatthew Knepley   CHKFORTRANNULLOBJECT(out);
270812c3f48SMatthew Knepley   *ierr = MatNullSpaceRemove(*sp,*vec,out);
271812c3f48SMatthew Knepley }
2721eea217eSSatish Balay 
2735ba43861SSatish Balay void PETSC_STDCALL   matgetinfo_(Mat *mat,MatInfoType *flag,MatInfo *info, int *__ierr )
2745ba43861SSatish Balay {
2755ba43861SSatish Balay   *__ierr = MatGetInfo(*mat,*flag,info);
2765ba43861SSatish Balay }
2775ba43861SSatish Balay 
278*22688b47SSatish Balay void PETSC_STDCALL   matlufactor_(Mat *mat,IS *row,IS *col, MatFactorInfo *info, int *__ierr )
279*22688b47SSatish Balay {
280*22688b47SSatish Balay   *__ierr = MatLUFactor(*mat,*row,*col,info);
281*22688b47SSatish Balay }
282*22688b47SSatish Balay 
283*22688b47SSatish Balay void PETSC_STDCALL   matilufactor_(Mat *mat,IS *row,IS *col, MatFactorInfo *info, int *__ierr )
284*22688b47SSatish Balay {
285*22688b47SSatish Balay   *__ierr = MatILUFactor(*mat,*row,*col,info);
286*22688b47SSatish Balay }
287*22688b47SSatish Balay 
288*22688b47SSatish Balay void PETSC_STDCALL   matlufactorsymbolic_(Mat *fact,Mat *mat,IS *row,IS *col, MatFactorInfo *info, int *__ierr )
289*22688b47SSatish Balay {
290*22688b47SSatish Balay   *__ierr = MatLUFactorSymbolic(*fact,*mat,*row,*col,info);
291*22688b47SSatish Balay }
292*22688b47SSatish Balay 
293*22688b47SSatish Balay void PETSC_STDCALL   matlufactornumeric_(Mat *fact,Mat *mat, MatFactorInfo *info, int *__ierr )
294*22688b47SSatish Balay {
295*22688b47SSatish Balay   *__ierr = MatLUFactorNumeric(*fact,*mat,info);
296*22688b47SSatish Balay }
297*22688b47SSatish Balay 
298*22688b47SSatish Balay void PETSC_STDCALL   matcholeskyfactor_(Mat *mat,IS *perm, MatFactorInfo *info, int *__ierr )
299*22688b47SSatish Balay {
300*22688b47SSatish Balay   *__ierr = MatCholeskyFactor(*mat,*perm,info);
301*22688b47SSatish Balay }
302*22688b47SSatish Balay 
303*22688b47SSatish Balay void PETSC_STDCALL   matcholeskyfactorsymbolic_(Mat *fact,Mat *mat,IS *perm, MatFactorInfo *info, int *__ierr )
304*22688b47SSatish Balay {
305*22688b47SSatish Balay   *__ierr = MatCholeskyFactorSymbolic(*fact,*mat,*perm,info);
306*22688b47SSatish Balay }
307*22688b47SSatish Balay 
308*22688b47SSatish Balay void PETSC_STDCALL   matcholeskyfactornumeric_(Mat *fact,Mat *mat, MatFactorInfo *info, int *__ierr )
309*22688b47SSatish Balay {
310*22688b47SSatish Balay   *__ierr = MatCholeskyFactorNumeric(*fact,*mat,info);
311*22688b47SSatish Balay }
312*22688b47SSatish Balay 
313*22688b47SSatish Balay void PETSC_STDCALL   matilufactorsymbolic_(Mat *fact,Mat *mat,IS *row,IS *col, MatFactorInfo *info, int *__ierr )
314*22688b47SSatish Balay {
315*22688b47SSatish Balay   *__ierr = MatILUFactorSymbolic(*fact,*mat,*row,*col,info);
316*22688b47SSatish Balay }
317*22688b47SSatish Balay 
318*22688b47SSatish Balay void PETSC_STDCALL   maticcfactorsymbolic_(Mat *fact,Mat *mat,IS *perm, MatFactorInfo *info, int *__ierr )
319*22688b47SSatish Balay {
320*22688b47SSatish Balay   *__ierr = MatICCFactorSymbolic(*fact,*mat,*perm,info);
321*22688b47SSatish Balay }
322*22688b47SSatish Balay 
323*22688b47SSatish Balay void PETSC_STDCALL   maticcfactor_(Mat *mat,IS *row, MatFactorInfo* info, int *__ierr )
324*22688b47SSatish Balay {
325*22688b47SSatish Balay   *__ierr = MatICCFactor(*mat,*row,info);
326*22688b47SSatish Balay }
327*22688b47SSatish Balay 
328*22688b47SSatish Balay void PETSC_STDCALL   matfactorinfoinitialize_(MatFactorInfo *info, int *__ierr )
329*22688b47SSatish Balay {
330*22688b47SSatish Balay   *__ierr = MatFactorInfoInitialize(info);
331*22688b47SSatish Balay }
332*22688b47SSatish Balay 
333f4e70085SSatish Balay EXTERN_C_END
334