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