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