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