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