xref: /petsc/src/mat/interface/ftn-custom/zmatrixf.c (revision 5ba4386168b6862bc91681d5c2fac74d965dee3f)
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
23*5ba43861SSatish Balay #define matnullspaceremove_              MATNULLSPACEREMOV
24*5ba43861SSatish Balay #define matgetinfo_                      MATGETINFO
25f4e70085SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
267d6bfa3bSBarry Smith #define matdestroymatrices_              matdestroymatrices_
275dffd610SBarry Smith #define matgetfactor_                    matgetfactor
2835bd34faSBarry Smith #define matfactorgetsolverpackage_       matfactorgetsolverpackage
297c54600cSBarry Smith #define matgetvecs_                      matgetvecs
30f4e70085SSatish Balay #define matgetrowij_                     matgetrowij
31f4e70085SSatish Balay #define matrestorerowij_                 matrestorerowij
32f4e70085SSatish Balay #define matgetrow_                       matgetrow
33f4e70085SSatish Balay #define matrestorerow_                   matrestorerow
34f4e70085SSatish Balay #define matview_                         matview
35f4e70085SSatish Balay #define matgetarray_                     matgetarray
36f4e70085SSatish Balay #define matrestorearray_                 matrestorearray
37f4e70085SSatish Balay #define matconvert_                      matconvert
38f4e70085SSatish Balay #define matgetsubmatrices_               matgetsubmatrices
39f4e70085SSatish Balay #define matzerorows_                     matzerorows
40f4e70085SSatish Balay #define matzerorowsis_                   matzerorowsis
41f4e70085SSatish Balay #define matzerorowslocal_                matzerorowslocal
42f4e70085SSatish Balay #define matzerorowslocalis_              matzerorowslocalis
431eea217eSSatish Balay #define matsetoptionsprefix_             matsetoptionsprefix
44812c3f48SMatthew Knepley #define matnullspaceremove_              matnullspaceremove
45*5ba43861SSatish Balay #define matgetinfo_                      matgetinfo
46f4e70085SSatish Balay #endif
47f4e70085SSatish Balay 
48f4e70085SSatish Balay EXTERN_C_BEGIN
49f4e70085SSatish Balay 
507c54600cSBarry Smith void PETSC_STDCALL   matgetvecs_(Mat *mat,Vec *right,Vec *left, int *ierr )
517c54600cSBarry Smith {
527c54600cSBarry Smith   CHKFORTRANNULLOBJECT(right);
537c54600cSBarry Smith   CHKFORTRANNULLOBJECT(left);
547c54600cSBarry Smith   *ierr = MatGetVecs(*mat,right,left);
557c54600cSBarry Smith }
567c54600cSBarry Smith 
578f7157efSSatish Balay void PETSC_STDCALL matgetrowij_(Mat *B,PetscInt *shift,PetscTruth *sym,PetscTruth *blockcompressed,PetscInt *n,PetscInt *ia,size_t *iia,
588f7157efSSatish Balay                                 PetscInt *ja,size_t *jja,PetscTruth *done,PetscErrorCode *ierr)
59f4e70085SSatish Balay {
60f4e70085SSatish Balay   PetscInt *IA,*JA;
618f7157efSSatish Balay   *ierr = MatGetRowIJ(*B,*shift,*sym,*blockcompressed,n,&IA,&JA,done);if (*ierr) return;
62f4e70085SSatish Balay   *iia  = PetscIntAddressToFortran(ia,IA);
63f4e70085SSatish Balay   *jja  = PetscIntAddressToFortran(ja,JA);
64f4e70085SSatish Balay }
65f4e70085SSatish Balay 
668f7157efSSatish Balay void PETSC_STDCALL matrestorerowij_(Mat *B,PetscInt *shift,PetscTruth *sym,PetscTruth *blockcompressed, PetscInt *n,PetscInt *ia,size_t *iia,
678f7157efSSatish Balay                                     PetscInt *ja,size_t *jja,PetscTruth *done,PetscErrorCode *ierr)
68f4e70085SSatish Balay {
69f4e70085SSatish Balay   PetscInt *IA = PetscIntAddressFromFortran(ia,*iia),*JA = PetscIntAddressFromFortran(ja,*jja);
708f7157efSSatish Balay   *ierr = MatRestoreRowIJ(*B,*shift,*sym,*blockcompressed,n,&IA,&JA,done);
71f4e70085SSatish Balay }
72f4e70085SSatish Balay 
73f4e70085SSatish Balay /*
74f4e70085SSatish Balay    This is a poor way of storing the column and value pointers
75f4e70085SSatish Balay   generated by MatGetRow() to be returned with MatRestoreRow()
76f4e70085SSatish Balay   but there is not natural,good place else to store them. Hence
77f4e70085SSatish Balay   Fortran programmers can only have one outstanding MatGetRows()
78f4e70085SSatish Balay   at a time.
79f4e70085SSatish Balay */
80f4e70085SSatish Balay static PetscErrorCode    matgetrowactive = 0;
81f4e70085SSatish Balay static const PetscInt    *my_ocols = 0;
82f4e70085SSatish Balay static const PetscScalar *my_ovals = 0;
83f4e70085SSatish Balay 
84f4e70085SSatish Balay void PETSC_STDCALL matgetrow_(Mat *mat,PetscInt *row,PetscInt *ncols,PetscInt *cols,PetscScalar *vals,PetscErrorCode *ierr)
85f4e70085SSatish Balay {
86f4e70085SSatish Balay   const PetscInt    **oocols = &my_ocols;
87f4e70085SSatish Balay   const PetscScalar **oovals = &my_ovals;
88f4e70085SSatish Balay 
89f4e70085SSatish Balay   if (matgetrowactive) {
90f4e70085SSatish Balay      PetscError(__LINE__,"MatGetRow_Fortran",__FILE__,__SDIR__,1,0,
91f4e70085SSatish Balay                "Cannot have two MatGetRow() active simultaneously\n\
92f4e70085SSatish Balay                call MatRestoreRow() before calling MatGetRow() a second time");
93f4e70085SSatish Balay      *ierr = 1;
94f4e70085SSatish Balay      return;
95f4e70085SSatish Balay   }
96f4e70085SSatish Balay 
97f4e70085SSatish Balay   CHKFORTRANNULLINTEGER(cols); if (!cols) oocols = PETSC_NULL;
98f4e70085SSatish Balay   CHKFORTRANNULLSCALAR(vals);  if (!vals) oovals = PETSC_NULL;
99f4e70085SSatish Balay 
100f4e70085SSatish Balay   *ierr = MatGetRow(*mat,*row,ncols,oocols,oovals);
101f4e70085SSatish Balay   if (*ierr) return;
102f4e70085SSatish Balay 
103f4e70085SSatish Balay   if (oocols) { *ierr = PetscMemcpy(cols,my_ocols,(*ncols)*sizeof(PetscInt)); if (*ierr) return;}
104f4e70085SSatish Balay   if (oovals) { *ierr = PetscMemcpy(vals,my_ovals,(*ncols)*sizeof(PetscScalar)); if (*ierr) return; }
105f4e70085SSatish Balay   matgetrowactive = 1;
106f4e70085SSatish Balay }
107f4e70085SSatish Balay 
108f4e70085SSatish Balay void PETSC_STDCALL matrestorerow_(Mat *mat,PetscInt *row,PetscInt *ncols,PetscInt *cols,PetscScalar *vals,PetscErrorCode *ierr)
109f4e70085SSatish Balay {
110f4e70085SSatish Balay   const PetscInt         **oocols = &my_ocols;
111f4e70085SSatish Balay   const PetscScalar **oovals = &my_ovals;
112f4e70085SSatish Balay   if (!matgetrowactive) {
113f4e70085SSatish Balay      PetscError(__LINE__,"MatRestoreRow_Fortran",__FILE__,__SDIR__,1,0,
114f4e70085SSatish Balay                "Must call MatGetRow() first");
115f4e70085SSatish Balay      *ierr = 1;
116f4e70085SSatish Balay      return;
117f4e70085SSatish Balay   }
118f4e70085SSatish Balay   CHKFORTRANNULLINTEGER(cols); if (!cols) oocols = PETSC_NULL;
119f4e70085SSatish Balay   CHKFORTRANNULLSCALAR(vals);  if (!vals) oovals = PETSC_NULL;
120f4e70085SSatish Balay 
121f4e70085SSatish Balay   *ierr = MatRestoreRow(*mat,*row,ncols,oocols,oovals);
122f4e70085SSatish Balay   matgetrowactive = 0;
123f4e70085SSatish Balay }
124f4e70085SSatish Balay 
125f4e70085SSatish Balay void PETSC_STDCALL matview_(Mat *mat,PetscViewer *vin,PetscErrorCode *ierr)
126f4e70085SSatish Balay {
127f4e70085SSatish Balay   PetscViewer v;
128f4e70085SSatish Balay   PetscPatchDefaultViewers_Fortran(vin,v);
129f4e70085SSatish Balay   *ierr = MatView(*mat,v);
130f4e70085SSatish Balay }
131f4e70085SSatish Balay 
132f4e70085SSatish Balay void PETSC_STDCALL matgetarray_(Mat *mat,PetscScalar *fa,size_t *ia,PetscErrorCode *ierr)
133f4e70085SSatish Balay {
134f4e70085SSatish Balay   PetscScalar *mm;
135f4e70085SSatish Balay   PetscInt    m,n;
136f4e70085SSatish Balay 
137f4e70085SSatish Balay   *ierr = MatGetArray(*mat,&mm); if (*ierr) return;
138f4e70085SSatish Balay   *ierr = MatGetSize(*mat,&m,&n);  if (*ierr) return;
139f91d1997SBarry Smith   *ierr = PetscScalarAddressToFortran((PetscObject)*mat,1,fa,mm,m*n,ia); if (*ierr) return;
140f4e70085SSatish Balay }
141f4e70085SSatish Balay 
142f4e70085SSatish Balay void PETSC_STDCALL matrestorearray_(Mat *mat,PetscScalar *fa,size_t *ia,PetscErrorCode *ierr)
143f4e70085SSatish Balay {
144f4e70085SSatish Balay   PetscScalar          *lx;
145f4e70085SSatish Balay   PetscInt                  m,n;
146f4e70085SSatish Balay 
147f4e70085SSatish Balay   *ierr = MatGetSize(*mat,&m,&n); if (*ierr) return;
148f4e70085SSatish Balay   *ierr = PetscScalarAddressFromFortran((PetscObject)*mat,fa,*ia,m*n,&lx);if (*ierr) return;
149f4e70085SSatish Balay   *ierr = MatRestoreArray(*mat,&lx);if (*ierr) return;
150f4e70085SSatish Balay }
151f4e70085SSatish Balay 
15235bd34faSBarry Smith void PETSC_STDCALL matfactorgetsolverpackage_(Mat *mat,CHAR name PETSC_MIXED_LEN(len),PetscErrorCode *ierr PETSC_END_LEN(len))
15335bd34faSBarry Smith {
15435bd34faSBarry Smith   const char *tname;
15535bd34faSBarry Smith 
15635bd34faSBarry Smith   *ierr = MatFactorGetSolverPackage(*mat,&tname);if (*ierr) return;
15735bd34faSBarry Smith   if (name != PETSC_NULL_CHARACTER_Fortran) {
15835bd34faSBarry Smith     *ierr = PetscStrncpy(name,tname,len);if (*ierr) return;
15935bd34faSBarry Smith   }
16035bd34faSBarry Smith   FIXRETURNCHAR(PETSC_TRUE,name,len);
16135bd34faSBarry Smith }
16235bd34faSBarry Smith 
1635dffd610SBarry Smith void PETSC_STDCALL matgetfactor_(Mat *mat,CHAR outtype PETSC_MIXED_LEN(len),MatFactorType ftype,Mat *M,PetscErrorCode *ierr PETSC_END_LEN(len))
1645dffd610SBarry Smith {
1655dffd610SBarry Smith   char *t;
1665dffd610SBarry Smith   FIXCHAR(outtype,len,t);
167c911e420SBarry Smith   *ierr = MatGetFactor(*mat,t,ftype,M);
1685dffd610SBarry Smith   FREECHAR(outtype,t);
1695dffd610SBarry Smith }
1705dffd610SBarry Smith 
171f4e70085SSatish Balay void PETSC_STDCALL matconvert_(Mat *mat,CHAR outtype PETSC_MIXED_LEN(len),MatReuse *reuse,Mat *M,PetscErrorCode *ierr PETSC_END_LEN(len))
172f4e70085SSatish Balay {
173f4e70085SSatish Balay   char *t;
174f4e70085SSatish Balay   FIXCHAR(outtype,len,t);
175f4e70085SSatish Balay   *ierr = MatConvert(*mat,t,*reuse,M);
176f4e70085SSatish Balay   FREECHAR(outtype,t);
177f4e70085SSatish Balay }
178f4e70085SSatish Balay 
179f4e70085SSatish Balay /*
180f4e70085SSatish Balay     MatGetSubmatrices() is slightly different from C since the
181f4e70085SSatish Balay     Fortran provides the array to hold the submatrix objects,while in C that
182f4e70085SSatish Balay     array is allocated by the MatGetSubmatrices()
183f4e70085SSatish Balay */
184f4e70085SSatish Balay void PETSC_STDCALL matgetsubmatrices_(Mat *mat,PetscInt *n,IS *isrow,IS *iscol,MatReuse *scall,Mat *smat,PetscErrorCode *ierr)
185f4e70085SSatish Balay {
186f4e70085SSatish Balay   Mat *lsmat;
187f4e70085SSatish Balay   PetscInt i;
188f4e70085SSatish Balay 
189f4e70085SSatish Balay   if (*scall == MAT_INITIAL_MATRIX) {
190f4e70085SSatish Balay     *ierr = MatGetSubMatrices(*mat,*n,isrow,iscol,*scall,&lsmat);
191f4e70085SSatish Balay     for (i=0; i<*n; i++) {
192f4e70085SSatish Balay       smat[i] = lsmat[i];
193f4e70085SSatish Balay     }
194f4e70085SSatish Balay     *ierr = PetscFree(lsmat);
195f4e70085SSatish Balay   } else {
196f4e70085SSatish Balay     *ierr = MatGetSubMatrices(*mat,*n,isrow,iscol,*scall,&smat);
197f4e70085SSatish Balay   }
198f4e70085SSatish Balay }
199f4e70085SSatish Balay 
2007d6bfa3bSBarry Smith /*
2017d6bfa3bSBarry Smith     MatDestroyMatrices() is slightly different from C since the
2027d6bfa3bSBarry Smith     Fortran provides the array to hold the submatrix objects,while in C that
2037d6bfa3bSBarry Smith     array is allocated by the MatGetSubmatrices()
2047d6bfa3bSBarry Smith */
2057d6bfa3bSBarry Smith void PETSC_STDCALL matdestroymatrices_(Mat *mat,PetscInt *n,Mat *smat,PetscErrorCode *ierr)
2067d6bfa3bSBarry Smith {
2077d6bfa3bSBarry Smith   PetscInt i;
2087d6bfa3bSBarry Smith 
2097d6bfa3bSBarry Smith   for (i=0; i<*n; i++) {
2107d6bfa3bSBarry Smith     *ierr = MatDestroy(smat[i]);if (*ierr) return;
2117d6bfa3bSBarry Smith   }
2127d6bfa3bSBarry Smith }
2137d6bfa3bSBarry Smith 
214f4e70085SSatish Balay void PETSC_STDCALL matzerorows_(Mat *mat,PetscInt *numRows,PetscInt *rows,PetscScalar *diag,PetscErrorCode *ierr)
215f4e70085SSatish Balay {
216f4e70085SSatish Balay   *ierr = MatZeroRows(*mat,*numRows,rows,*diag);
217f4e70085SSatish Balay }
218f4e70085SSatish Balay 
219f4e70085SSatish Balay void PETSC_STDCALL matzerorowsis_(Mat *mat,IS *is,PetscScalar *diag,PetscErrorCode *ierr)
220f4e70085SSatish Balay {
221f4e70085SSatish Balay   *ierr = MatZeroRowsIS(*mat,*is,*diag);
222f4e70085SSatish Balay }
223f4e70085SSatish Balay 
224f4e70085SSatish Balay void PETSC_STDCALL matzerorowslocal_(Mat *mat,PetscInt *numRows,PetscInt *rows,PetscScalar *diag,PetscErrorCode *ierr)
225f4e70085SSatish Balay {
226f4e70085SSatish Balay   *ierr = MatZeroRowsLocal(*mat,*numRows,rows,*diag);
227f4e70085SSatish Balay }
228f4e70085SSatish Balay 
229f4e70085SSatish Balay void PETSC_STDCALL matzerorowslocalis_(Mat *mat,IS *is,PetscScalar *diag,PetscErrorCode *ierr)
230f4e70085SSatish Balay {
231f4e70085SSatish Balay   *ierr = MatZeroRowsLocalIS(*mat,*is,*diag);
232f4e70085SSatish Balay }
233f4e70085SSatish Balay 
2341eea217eSSatish Balay 
2351eea217eSSatish Balay void PETSC_STDCALL matsetoptionsprefix_(Mat *mat,CHAR prefix PETSC_MIXED_LEN(len),
2361eea217eSSatish Balay                                         PetscErrorCode *ierr PETSC_END_LEN(len))
2371eea217eSSatish Balay {
2381eea217eSSatish Balay   char *t;
2391eea217eSSatish Balay 
2401eea217eSSatish Balay   FIXCHAR(prefix,len,t);
2411eea217eSSatish Balay   *ierr = MatSetOptionsPrefix(*mat,t);
2421eea217eSSatish Balay   FREECHAR(prefix,t);
2431eea217eSSatish Balay }
2441eea217eSSatish Balay 
245812c3f48SMatthew Knepley void PETSC_STDCALL matnullspaceremove_(MatNullSpace *sp,Vec *vec,Vec *out,PetscErrorCode *ierr)
246812c3f48SMatthew Knepley {
247812c3f48SMatthew Knepley   CHKFORTRANNULLOBJECT(out);
248812c3f48SMatthew Knepley   *ierr = MatNullSpaceRemove(*sp,*vec,out);
249812c3f48SMatthew Knepley }
2501eea217eSSatish Balay 
251*5ba43861SSatish Balay void PETSC_STDCALL   matgetinfo_(Mat *mat,MatInfoType *flag,MatInfo *info, int *__ierr )
252*5ba43861SSatish Balay {
253*5ba43861SSatish Balay   *__ierr = MatGetInfo(*mat,*flag,info);
254*5ba43861SSatish Balay }
255*5ba43861SSatish Balay 
256f4e70085SSatish Balay EXTERN_C_END
257