xref: /petsc/src/mat/interface/ftn-custom/zmatrixf.c (revision 7c54600c425fb8d0ad4ad4bb40f7fac17c980c51)
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
19*7c54600cSBarry Smith #define matgetvecs_                      MATGETVECS
20f4e70085SSatish Balay #elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
21*7c54600cSBarry Smith #define matgetvecs_                      matgetvecs
22f4e70085SSatish Balay #define matgetrowij_                     matgetrowij
23f4e70085SSatish Balay #define matrestorerowij_                 matrestorerowij
24f4e70085SSatish Balay #define matgetrow_                       matgetrow
25f4e70085SSatish Balay #define matrestorerow_                   matrestorerow
26f4e70085SSatish Balay #define matview_                         matview
27f4e70085SSatish Balay #define matgetarray_                     matgetarray
28f4e70085SSatish Balay #define matrestorearray_                 matrestorearray
29f4e70085SSatish Balay #define matconvert_                      matconvert
30f4e70085SSatish Balay #define matgetsubmatrices_               matgetsubmatrices
31f4e70085SSatish Balay #define matzerorows_                     matzerorows
32f4e70085SSatish Balay #define matzerorowsis_                   matzerorowsis
33f4e70085SSatish Balay #define matzerorowslocal_                matzerorowslocal
34f4e70085SSatish Balay #define matzerorowslocalis_              matzerorowslocalis
351eea217eSSatish Balay #define matsetoptionsprefix_             matsetoptionsprefix
36f4e70085SSatish Balay #endif
37f4e70085SSatish Balay 
38f4e70085SSatish Balay EXTERN_C_BEGIN
39f4e70085SSatish Balay 
40*7c54600cSBarry Smith void PETSC_STDCALL   matgetvecs_(Mat *mat,Vec *right,Vec *left, int *ierr )
41*7c54600cSBarry Smith {
42*7c54600cSBarry Smith   CHKFORTRANNULLOBJECT(right);
43*7c54600cSBarry Smith   CHKFORTRANNULLOBJECT(left);
44*7c54600cSBarry Smith   *ierr = MatGetVecs(*mat,right,left);
45*7c54600cSBarry Smith }
46*7c54600cSBarry Smith 
478f7157efSSatish Balay void PETSC_STDCALL matgetrowij_(Mat *B,PetscInt *shift,PetscTruth *sym,PetscTruth *blockcompressed,PetscInt *n,PetscInt *ia,size_t *iia,
488f7157efSSatish Balay                                 PetscInt *ja,size_t *jja,PetscTruth *done,PetscErrorCode *ierr)
49f4e70085SSatish Balay {
50f4e70085SSatish Balay   PetscInt *IA,*JA;
518f7157efSSatish Balay   *ierr = MatGetRowIJ(*B,*shift,*sym,*blockcompressed,n,&IA,&JA,done);if (*ierr) return;
52f4e70085SSatish Balay   *iia  = PetscIntAddressToFortran(ia,IA);
53f4e70085SSatish Balay   *jja  = PetscIntAddressToFortran(ja,JA);
54f4e70085SSatish Balay }
55f4e70085SSatish Balay 
568f7157efSSatish Balay void PETSC_STDCALL matrestorerowij_(Mat *B,PetscInt *shift,PetscTruth *sym,PetscTruth *blockcompressed, PetscInt *n,PetscInt *ia,size_t *iia,
578f7157efSSatish Balay                                     PetscInt *ja,size_t *jja,PetscTruth *done,PetscErrorCode *ierr)
58f4e70085SSatish Balay {
59f4e70085SSatish Balay   PetscInt *IA = PetscIntAddressFromFortran(ia,*iia),*JA = PetscIntAddressFromFortran(ja,*jja);
608f7157efSSatish Balay   *ierr = MatRestoreRowIJ(*B,*shift,*sym,*blockcompressed,n,&IA,&JA,done);
61f4e70085SSatish Balay }
62f4e70085SSatish Balay 
63f4e70085SSatish Balay /*
64f4e70085SSatish Balay    This is a poor way of storing the column and value pointers
65f4e70085SSatish Balay   generated by MatGetRow() to be returned with MatRestoreRow()
66f4e70085SSatish Balay   but there is not natural,good place else to store them. Hence
67f4e70085SSatish Balay   Fortran programmers can only have one outstanding MatGetRows()
68f4e70085SSatish Balay   at a time.
69f4e70085SSatish Balay */
70f4e70085SSatish Balay static PetscErrorCode    matgetrowactive = 0;
71f4e70085SSatish Balay static const PetscInt    *my_ocols = 0;
72f4e70085SSatish Balay static const PetscScalar *my_ovals = 0;
73f4e70085SSatish Balay 
74f4e70085SSatish Balay void PETSC_STDCALL matgetrow_(Mat *mat,PetscInt *row,PetscInt *ncols,PetscInt *cols,PetscScalar *vals,PetscErrorCode *ierr)
75f4e70085SSatish Balay {
76f4e70085SSatish Balay   const PetscInt    **oocols = &my_ocols;
77f4e70085SSatish Balay   const PetscScalar **oovals = &my_ovals;
78f4e70085SSatish Balay 
79f4e70085SSatish Balay   if (matgetrowactive) {
80f4e70085SSatish Balay      PetscError(__LINE__,"MatGetRow_Fortran",__FILE__,__SDIR__,1,0,
81f4e70085SSatish Balay                "Cannot have two MatGetRow() active simultaneously\n\
82f4e70085SSatish Balay                call MatRestoreRow() before calling MatGetRow() a second time");
83f4e70085SSatish Balay      *ierr = 1;
84f4e70085SSatish Balay      return;
85f4e70085SSatish Balay   }
86f4e70085SSatish Balay 
87f4e70085SSatish Balay   CHKFORTRANNULLINTEGER(cols); if (!cols) oocols = PETSC_NULL;
88f4e70085SSatish Balay   CHKFORTRANNULLSCALAR(vals);  if (!vals) oovals = PETSC_NULL;
89f4e70085SSatish Balay 
90f4e70085SSatish Balay   *ierr = MatGetRow(*mat,*row,ncols,oocols,oovals);
91f4e70085SSatish Balay   if (*ierr) return;
92f4e70085SSatish Balay 
93f4e70085SSatish Balay   if (oocols) { *ierr = PetscMemcpy(cols,my_ocols,(*ncols)*sizeof(PetscInt)); if (*ierr) return;}
94f4e70085SSatish Balay   if (oovals) { *ierr = PetscMemcpy(vals,my_ovals,(*ncols)*sizeof(PetscScalar)); if (*ierr) return; }
95f4e70085SSatish Balay   matgetrowactive = 1;
96f4e70085SSatish Balay }
97f4e70085SSatish Balay 
98f4e70085SSatish Balay void PETSC_STDCALL matrestorerow_(Mat *mat,PetscInt *row,PetscInt *ncols,PetscInt *cols,PetscScalar *vals,PetscErrorCode *ierr)
99f4e70085SSatish Balay {
100f4e70085SSatish Balay   const PetscInt         **oocols = &my_ocols;
101f4e70085SSatish Balay   const PetscScalar **oovals = &my_ovals;
102f4e70085SSatish Balay   if (!matgetrowactive) {
103f4e70085SSatish Balay      PetscError(__LINE__,"MatRestoreRow_Fortran",__FILE__,__SDIR__,1,0,
104f4e70085SSatish Balay                "Must call MatGetRow() first");
105f4e70085SSatish Balay      *ierr = 1;
106f4e70085SSatish Balay      return;
107f4e70085SSatish Balay   }
108f4e70085SSatish Balay   CHKFORTRANNULLINTEGER(cols); if (!cols) oocols = PETSC_NULL;
109f4e70085SSatish Balay   CHKFORTRANNULLSCALAR(vals);  if (!vals) oovals = PETSC_NULL;
110f4e70085SSatish Balay 
111f4e70085SSatish Balay   *ierr = MatRestoreRow(*mat,*row,ncols,oocols,oovals);
112f4e70085SSatish Balay   matgetrowactive = 0;
113f4e70085SSatish Balay }
114f4e70085SSatish Balay 
115f4e70085SSatish Balay void PETSC_STDCALL matview_(Mat *mat,PetscViewer *vin,PetscErrorCode *ierr)
116f4e70085SSatish Balay {
117f4e70085SSatish Balay   PetscViewer v;
118f4e70085SSatish Balay   PetscPatchDefaultViewers_Fortran(vin,v);
119f4e70085SSatish Balay   *ierr = MatView(*mat,v);
120f4e70085SSatish Balay }
121f4e70085SSatish Balay 
122f4e70085SSatish Balay void PETSC_STDCALL matgetarray_(Mat *mat,PetscScalar *fa,size_t *ia,PetscErrorCode *ierr)
123f4e70085SSatish Balay {
124f4e70085SSatish Balay   PetscScalar *mm;
125f4e70085SSatish Balay   PetscInt    m,n;
126f4e70085SSatish Balay 
127f4e70085SSatish Balay   *ierr = MatGetArray(*mat,&mm); if (*ierr) return;
128f4e70085SSatish Balay   *ierr = MatGetSize(*mat,&m,&n);  if (*ierr) return;
129f91d1997SBarry Smith   *ierr = PetscScalarAddressToFortran((PetscObject)*mat,1,fa,mm,m*n,ia); if (*ierr) return;
130f4e70085SSatish Balay }
131f4e70085SSatish Balay 
132f4e70085SSatish Balay void PETSC_STDCALL matrestorearray_(Mat *mat,PetscScalar *fa,size_t *ia,PetscErrorCode *ierr)
133f4e70085SSatish Balay {
134f4e70085SSatish Balay   PetscScalar          *lx;
135f4e70085SSatish Balay   PetscInt                  m,n;
136f4e70085SSatish Balay 
137f4e70085SSatish Balay   *ierr = MatGetSize(*mat,&m,&n); if (*ierr) return;
138f4e70085SSatish Balay   *ierr = PetscScalarAddressFromFortran((PetscObject)*mat,fa,*ia,m*n,&lx);if (*ierr) return;
139f4e70085SSatish Balay   *ierr = MatRestoreArray(*mat,&lx);if (*ierr) return;
140f4e70085SSatish Balay }
141f4e70085SSatish Balay 
142f4e70085SSatish Balay void PETSC_STDCALL matconvert_(Mat *mat,CHAR outtype PETSC_MIXED_LEN(len),MatReuse *reuse,Mat *M,PetscErrorCode *ierr PETSC_END_LEN(len))
143f4e70085SSatish Balay {
144f4e70085SSatish Balay   char *t;
145f4e70085SSatish Balay   FIXCHAR(outtype,len,t);
146f4e70085SSatish Balay   *ierr = MatConvert(*mat,t,*reuse,M);
147f4e70085SSatish Balay   FREECHAR(outtype,t);
148f4e70085SSatish Balay }
149f4e70085SSatish Balay 
150f4e70085SSatish Balay /*
151f4e70085SSatish Balay     MatGetSubmatrices() is slightly different from C since the
152f4e70085SSatish Balay     Fortran provides the array to hold the submatrix objects,while in C that
153f4e70085SSatish Balay     array is allocated by the MatGetSubmatrices()
154f4e70085SSatish Balay */
155f4e70085SSatish Balay void PETSC_STDCALL matgetsubmatrices_(Mat *mat,PetscInt *n,IS *isrow,IS *iscol,MatReuse *scall,Mat *smat,PetscErrorCode *ierr)
156f4e70085SSatish Balay {
157f4e70085SSatish Balay   Mat *lsmat;
158f4e70085SSatish Balay   PetscInt i;
159f4e70085SSatish Balay 
160f4e70085SSatish Balay   if (*scall == MAT_INITIAL_MATRIX) {
161f4e70085SSatish Balay     *ierr = MatGetSubMatrices(*mat,*n,isrow,iscol,*scall,&lsmat);
162f4e70085SSatish Balay     for (i=0; i<*n; i++) {
163f4e70085SSatish Balay       smat[i] = lsmat[i];
164f4e70085SSatish Balay     }
165f4e70085SSatish Balay     *ierr = PetscFree(lsmat);
166f4e70085SSatish Balay   } else {
167f4e70085SSatish Balay     *ierr = MatGetSubMatrices(*mat,*n,isrow,iscol,*scall,&smat);
168f4e70085SSatish Balay   }
169f4e70085SSatish Balay }
170f4e70085SSatish Balay 
171f4e70085SSatish Balay void PETSC_STDCALL matzerorows_(Mat *mat,PetscInt *numRows,PetscInt *rows,PetscScalar *diag,PetscErrorCode *ierr)
172f4e70085SSatish Balay {
173f4e70085SSatish Balay   *ierr = MatZeroRows(*mat,*numRows,rows,*diag);
174f4e70085SSatish Balay }
175f4e70085SSatish Balay 
176f4e70085SSatish Balay void PETSC_STDCALL matzerorowsis_(Mat *mat,IS *is,PetscScalar *diag,PetscErrorCode *ierr)
177f4e70085SSatish Balay {
178f4e70085SSatish Balay   *ierr = MatZeroRowsIS(*mat,*is,*diag);
179f4e70085SSatish Balay }
180f4e70085SSatish Balay 
181f4e70085SSatish Balay void PETSC_STDCALL matzerorowslocal_(Mat *mat,PetscInt *numRows,PetscInt *rows,PetscScalar *diag,PetscErrorCode *ierr)
182f4e70085SSatish Balay {
183f4e70085SSatish Balay   *ierr = MatZeroRowsLocal(*mat,*numRows,rows,*diag);
184f4e70085SSatish Balay }
185f4e70085SSatish Balay 
186f4e70085SSatish Balay void PETSC_STDCALL matzerorowslocalis_(Mat *mat,IS *is,PetscScalar *diag,PetscErrorCode *ierr)
187f4e70085SSatish Balay {
188f4e70085SSatish Balay   *ierr = MatZeroRowsLocalIS(*mat,*is,*diag);
189f4e70085SSatish Balay }
190f4e70085SSatish Balay 
1911eea217eSSatish Balay 
1921eea217eSSatish Balay void PETSC_STDCALL matsetoptionsprefix_(Mat *mat,CHAR prefix PETSC_MIXED_LEN(len),
1931eea217eSSatish Balay                                         PetscErrorCode *ierr PETSC_END_LEN(len))
1941eea217eSSatish Balay {
1951eea217eSSatish Balay   char *t;
1961eea217eSSatish Balay 
1971eea217eSSatish Balay   FIXCHAR(prefix,len,t);
1981eea217eSSatish Balay   *ierr = MatSetOptionsPrefix(*mat,t);
1991eea217eSSatish Balay   FREECHAR(prefix,t);
2001eea217eSSatish Balay }
2011eea217eSSatish Balay 
2021eea217eSSatish Balay 
203f4e70085SSatish Balay EXTERN_C_END
204