1cb512458SBarry Smith #ifndef lint 2*932b0c3eSLois Curfman McInnes static char vcid[] = "$Id: dense.c,v 1.65 1995/10/13 02:05:06 curfman Exp bsmith $"; 3cb512458SBarry Smith #endif 4289bc588SBarry Smith 5289bc588SBarry Smith /* 6289bc588SBarry Smith Standard Fortran style matrices 7289bc588SBarry Smith */ 819b02663SBarry Smith #include "petsc.h" 9b16a3bb1SBarry Smith #include "pinclude/plapack.h" 10289bc588SBarry Smith #include "matimpl.h" 11289bc588SBarry Smith #include "math.h" 12289bc588SBarry Smith #include "vec/vecimpl.h" 13b16a3bb1SBarry Smith #include "pinclude/pviewer.h" 14289bc588SBarry Smith 15289bc588SBarry Smith typedef struct { 16c0bbcb79SLois Curfman McInnes Scalar *v; /* matrix elements */ 17c0bbcb79SLois Curfman McInnes int roworiented; /* if true, row oriented input (default) */ 18c0bbcb79SLois Curfman McInnes int m,n,pad; /* rows, columns, padding */ 19289bc588SBarry Smith int *pivots; /* pivots in LU factorization */ 20ec8511deSBarry Smith } Mat_SeqDense; 21289bc588SBarry Smith 221987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y) 231987afe7SBarry Smith { 241987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data; 251987afe7SBarry Smith int N = x->m*x->n, one = 1; 261987afe7SBarry Smith BLaxpy_( &N, alpha, x->v, &one, y->v, &one ); 271987afe7SBarry Smith PLogFlops(2*N-1); 281987afe7SBarry Smith return 0; 291987afe7SBarry Smith } 301987afe7SBarry Smith 31c0bbcb79SLois Curfman McInnes static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 32289bc588SBarry Smith { 33c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 34289bc588SBarry Smith int i,N = mat->m*mat->n,count = 0; 35289bc588SBarry Smith Scalar *v = mat->v; 36289bc588SBarry Smith for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 37c0bbcb79SLois Curfman McInnes *nz = count; *nzalloc = N; *mem = (int)A->mem; 38fa9ec3f1SLois Curfman McInnes return 0; 39289bc588SBarry Smith } 40289bc588SBarry Smith 41289bc588SBarry Smith /* ---------------------------------------------------------------*/ 42289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots, 43289bc588SBarry Smith rather than put it in the Mat->col slot.*/ 44c0bbcb79SLois Curfman McInnes static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f) 45289bc588SBarry Smith { 46c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 476abc6512SBarry Smith int info; 48289bc588SBarry Smith if (!mat->pivots) { 4948d91487SBarry Smith mat->pivots = (int *) PETSCMALLOC(mat->m*sizeof(int));CHKPTRQ(mat->pivots); 50c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,mat->m*sizeof(int)); 51289bc588SBarry Smith } 52289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 53ec8511deSBarry Smith if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization"); 54c0bbcb79SLois Curfman McInnes A->factor = FACTOR_LU; 55289bc588SBarry Smith return 0; 56289bc588SBarry Smith } 57c0bbcb79SLois Curfman McInnes static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact) 58289bc588SBarry Smith { 59289bc588SBarry Smith int ierr; 60c0bbcb79SLois Curfman McInnes ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr); 61289bc588SBarry Smith return 0; 62289bc588SBarry Smith } 63c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 64289bc588SBarry Smith { 6549d8b64dSBarry Smith return MatLUFactor(*fact,0,0,1.0); 66289bc588SBarry Smith } 67c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact) 68289bc588SBarry Smith { 69289bc588SBarry Smith int ierr; 70c0bbcb79SLois Curfman McInnes ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr); 71289bc588SBarry Smith return 0; 72289bc588SBarry Smith } 73c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 74289bc588SBarry Smith { 7549d8b64dSBarry Smith return MatCholeskyFactor(*fact,0,1.0); 76289bc588SBarry Smith } 77c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f) 78289bc588SBarry Smith { 79c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 806abc6512SBarry Smith int info; 81752f5567SLois Curfman McInnes if (mat->pivots) { 82752f5567SLois Curfman McInnes PETSCFREE(mat->pivots); 83c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,-mat->m*sizeof(int)); 84752f5567SLois Curfman McInnes mat->pivots = 0; 85752f5567SLois Curfman McInnes } 86289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 87ec8511deSBarry Smith if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization"); 88c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 89289bc588SBarry Smith return 0; 90289bc588SBarry Smith } 91289bc588SBarry Smith 92c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 93289bc588SBarry Smith { 94c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 956abc6512SBarry Smith int one = 1, info; 966abc6512SBarry Smith Scalar *x, *y; 97289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 98416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 99c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 10048d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 101289bc588SBarry Smith } 102c0bbcb79SLois Curfman McInnes else if (A->factor == FACTOR_CHOLESKY){ 10348d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 104289bc588SBarry Smith } 105ec8511deSBarry Smith else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve"); 106ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve"); 107289bc588SBarry Smith return 0; 108289bc588SBarry Smith } 109c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy) 110da3a660dSBarry Smith { 111c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1126abc6512SBarry Smith int one = 1, info; 1136abc6512SBarry Smith Scalar *x, *y; 114da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 115416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 116752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 117da3a660dSBarry Smith if (mat->pivots) { 11848d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 119da3a660dSBarry Smith } 120da3a660dSBarry Smith else { 12148d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 122da3a660dSBarry Smith } 123ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve"); 124da3a660dSBarry Smith return 0; 125da3a660dSBarry Smith } 126c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 127da3a660dSBarry Smith { 128c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1296abc6512SBarry Smith int one = 1, info,ierr; 1306abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 131da3a660dSBarry Smith Vec tmp = 0; 132da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 133da3a660dSBarry Smith if (yy == zz) { 13478b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 135c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 13678b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 137da3a660dSBarry Smith } 138416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 139752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 140da3a660dSBarry Smith if (mat->pivots) { 14148d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 142da3a660dSBarry Smith } 143da3a660dSBarry Smith else { 14448d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 145da3a660dSBarry Smith } 146ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve"); 147da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 148da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 149da3a660dSBarry Smith return 0; 150da3a660dSBarry Smith } 151c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 152da3a660dSBarry Smith { 153c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1546abc6512SBarry Smith int one = 1, info,ierr; 1556abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 156da3a660dSBarry Smith Vec tmp; 157da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 158da3a660dSBarry Smith if (yy == zz) { 15978b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 160c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 16178b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 162da3a660dSBarry Smith } 163416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 164752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 165da3a660dSBarry Smith if (mat->pivots) { 16648d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 167da3a660dSBarry Smith } 168da3a660dSBarry Smith else { 16948d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 170da3a660dSBarry Smith } 171ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve"); 172da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 173da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 174da3a660dSBarry Smith return 0; 175da3a660dSBarry Smith } 176289bc588SBarry Smith /* ------------------------------------------------------------------*/ 177c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 17820e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 179289bc588SBarry Smith { 180c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 181289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 1826abc6512SBarry Smith int o = 1,ierr, m = mat->m, i; 183289bc588SBarry Smith 184289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 185289bc588SBarry Smith /* this is a hack fix, should have another version without 186289bc588SBarry Smith the second BLdot */ 187bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx); CHKERRQ(ierr); 188289bc588SBarry Smith } 189289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 190289bc588SBarry Smith while (its--) { 191289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 192289bc588SBarry Smith for ( i=0; i<m; i++ ) { 193f1747703SBarry Smith #if defined(PETSC_COMPLEX) 194f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 195f1747703SBarry Smith not happy about returning a double complex */ 196f1747703SBarry Smith int _i; 197f1747703SBarry Smith Scalar sum = b[i]; 198f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 199f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 200f1747703SBarry Smith } 201f1747703SBarry Smith xt = sum; 202f1747703SBarry Smith #else 203289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 204f1747703SBarry Smith #endif 205d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 206289bc588SBarry Smith } 207289bc588SBarry Smith } 208289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 209289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 210f1747703SBarry Smith #if defined(PETSC_COMPLEX) 211f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 212f1747703SBarry Smith not happy about returning a double complex */ 213f1747703SBarry Smith int _i; 214f1747703SBarry Smith Scalar sum = b[i]; 215f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 216f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 217f1747703SBarry Smith } 218f1747703SBarry Smith xt = sum; 219f1747703SBarry Smith #else 220289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 221f1747703SBarry Smith #endif 222d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 223289bc588SBarry Smith } 224289bc588SBarry Smith } 225289bc588SBarry Smith } 226289bc588SBarry Smith return 0; 227289bc588SBarry Smith } 228289bc588SBarry Smith 229289bc588SBarry Smith /* -----------------------------------------------------------------*/ 230c0bbcb79SLois Curfman McInnes static int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy) 231289bc588SBarry Smith { 232c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 233289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 234289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 235289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 23648d91487SBarry Smith LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 237289bc588SBarry Smith return 0; 238289bc588SBarry Smith } 239c0bbcb79SLois Curfman McInnes static int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 240289bc588SBarry Smith { 241c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 242289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 243289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 244289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 24548d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 246289bc588SBarry Smith return 0; 247289bc588SBarry Smith } 248c0bbcb79SLois Curfman McInnes static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 249289bc588SBarry Smith { 250c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 251289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 2526abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 253289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 254416022c9SBarry Smith if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 25548d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 256289bc588SBarry Smith return 0; 257289bc588SBarry Smith } 258c0bbcb79SLois Curfman McInnes static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 259289bc588SBarry Smith { 260c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 261289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 262289bc588SBarry Smith int _One=1; 2636abc6512SBarry Smith Scalar _DOne=1.0; 264289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 265289bc588SBarry Smith VecGetArray(zz,&z); 266416022c9SBarry Smith if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 26748d91487SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 268289bc588SBarry Smith return 0; 269289bc588SBarry Smith } 270289bc588SBarry Smith 271289bc588SBarry Smith /* -----------------------------------------------------------------*/ 272c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 273289bc588SBarry Smith { 274c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 275289bc588SBarry Smith Scalar *v; 276289bc588SBarry Smith int i; 277289bc588SBarry Smith *ncols = mat->n; 278289bc588SBarry Smith if (cols) { 27978b31e54SBarry Smith *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols); 28073c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 281289bc588SBarry Smith } 282289bc588SBarry Smith if (vals) { 28378b31e54SBarry Smith *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 284289bc588SBarry Smith v = mat->v + row; 28573c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 286289bc588SBarry Smith } 287289bc588SBarry Smith return 0; 288289bc588SBarry Smith } 289c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 290289bc588SBarry Smith { 29178b31e54SBarry Smith if (cols) { PETSCFREE(*cols); } 29278b31e54SBarry Smith if (vals) { PETSCFREE(*vals); } 293289bc588SBarry Smith return 0; 294289bc588SBarry Smith } 295289bc588SBarry Smith /* ----------------------------------------------------------------*/ 296c0bbcb79SLois Curfman McInnes static int MatInsert_SeqDense(Mat A,int m,int *indexm,int n, 297e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 298289bc588SBarry Smith { 299c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 300289bc588SBarry Smith int i,j; 301d6dfbf8fSBarry Smith 302289bc588SBarry Smith if (!mat->roworiented) { 303dbb450caSBarry Smith if (addv == INSERT_VALUES) { 304289bc588SBarry Smith for ( j=0; j<n; j++ ) { 305289bc588SBarry Smith for ( i=0; i<m; i++ ) { 306289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 307289bc588SBarry Smith } 308289bc588SBarry Smith } 309289bc588SBarry Smith } 310289bc588SBarry Smith else { 311289bc588SBarry Smith for ( j=0; j<n; j++ ) { 312289bc588SBarry Smith for ( i=0; i<m; i++ ) { 313289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 314289bc588SBarry Smith } 315289bc588SBarry Smith } 316289bc588SBarry Smith } 317e8d4e0b9SBarry Smith } 318e8d4e0b9SBarry Smith else { 319dbb450caSBarry Smith if (addv == INSERT_VALUES) { 320e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 321e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 322e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 323e8d4e0b9SBarry Smith } 324e8d4e0b9SBarry Smith } 325e8d4e0b9SBarry Smith } 326289bc588SBarry Smith else { 327289bc588SBarry Smith for ( i=0; i<m; i++ ) { 328289bc588SBarry Smith for ( j=0; j<n; j++ ) { 329289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 330289bc588SBarry Smith } 331289bc588SBarry Smith } 332289bc588SBarry Smith } 333e8d4e0b9SBarry Smith } 334289bc588SBarry Smith return 0; 335289bc588SBarry Smith } 336e8d4e0b9SBarry Smith 337289bc588SBarry Smith /* -----------------------------------------------------------------*/ 338c0bbcb79SLois Curfman McInnes static int MatCopyPrivate_SeqDense(Mat A,Mat *newmat) 339289bc588SBarry Smith { 340c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l; 341289bc588SBarry Smith int ierr; 342289bc588SBarry Smith Mat newi; 34348d91487SBarry Smith 344c0bbcb79SLois Curfman McInnes ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,&newi); CHKERRQ(ierr); 345ec8511deSBarry Smith l = (Mat_SeqDense *) newi->data; 346416022c9SBarry Smith PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 347289bc588SBarry Smith *newmat = newi; 348289bc588SBarry Smith return 0; 349289bc588SBarry Smith } 350289bc588SBarry Smith 351*932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h" 352*932b0c3eSLois Curfman McInnes #include "sysio.h" 353*932b0c3eSLois Curfman McInnes 354*932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 355289bc588SBarry Smith { 356*932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 357*932b0c3eSLois Curfman McInnes int ierr, i, j, format; 358d636dbe3SBarry Smith FILE *fd; 359*932b0c3eSLois Curfman McInnes char *outputname; 360*932b0c3eSLois Curfman McInnes Scalar *v; 361*932b0c3eSLois Curfman McInnes 362*932b0c3eSLois Curfman McInnes ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 363*932b0c3eSLois Curfman McInnes ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 364*932b0c3eSLois Curfman McInnes ierr = ViewerFileGetFormat_Private(viewer,&format); 365f72585f2SLois Curfman McInnes if (format == FILE_FORMAT_INFO) { 366*932b0c3eSLois Curfman McInnes ; /* do nothing for now */ 367f72585f2SLois Curfman McInnes } 368f72585f2SLois Curfman McInnes else { 369*932b0c3eSLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 370*932b0c3eSLois Curfman McInnes v = a->v + i; 371*932b0c3eSLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 372289bc588SBarry Smith #if defined(PETSC_COMPLEX) 373*932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m; 374289bc588SBarry Smith #else 375*932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e ",*v); v += a->m; 376289bc588SBarry Smith #endif 377289bc588SBarry Smith } 3788e37dc09SBarry Smith fprintf(fd,"\n"); 379289bc588SBarry Smith } 380da3a660dSBarry Smith } 381*932b0c3eSLois Curfman McInnes fflush(fd); 382289bc588SBarry Smith return 0; 383289bc588SBarry Smith } 384289bc588SBarry Smith 385*932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 386*932b0c3eSLois Curfman McInnes { 387*932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 388*932b0c3eSLois Curfman McInnes int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 389*932b0c3eSLois Curfman McInnes Scalar *v, *anonz; 390*932b0c3eSLois Curfman McInnes 391*932b0c3eSLois Curfman McInnes ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 392*932b0c3eSLois Curfman McInnes col_lens = (int *) PETSCMALLOC( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens); 393*932b0c3eSLois Curfman McInnes col_lens[0] = MAT_COOKIE; 394*932b0c3eSLois Curfman McInnes col_lens[1] = m; 395*932b0c3eSLois Curfman McInnes col_lens[2] = n; 396*932b0c3eSLois Curfman McInnes col_lens[3] = nz; 397*932b0c3eSLois Curfman McInnes 398*932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 399*932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) col_lens[4+i] = n; 400*932b0c3eSLois Curfman McInnes ierr = SYWrite(fd,col_lens,4+m,SYINT,1); CHKERRQ(ierr); 401*932b0c3eSLois Curfman McInnes 402*932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 403*932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 404*932b0c3eSLois Curfman McInnes ict = 0; 405*932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 406*932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) col_lens[ict++] = j; 407*932b0c3eSLois Curfman McInnes } 408*932b0c3eSLois Curfman McInnes ierr = SYWrite(fd,col_lens,nz,SYINT,0); CHKERRQ(ierr); 409*932b0c3eSLois Curfman McInnes PETSCFREE(col_lens); 410*932b0c3eSLois Curfman McInnes 411*932b0c3eSLois Curfman McInnes /* store nonzero values */ 412*932b0c3eSLois Curfman McInnes anonz = (Scalar *) PETSCMALLOC((nz)*sizeof(Scalar)); CHKPTRQ(anonz); 413*932b0c3eSLois Curfman McInnes ict = 0; 414*932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 415*932b0c3eSLois Curfman McInnes v = a->v + i; 416*932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) { 417*932b0c3eSLois Curfman McInnes anonz[ict++] = *v; v += a->m; 418*932b0c3eSLois Curfman McInnes } 419*932b0c3eSLois Curfman McInnes } 420*932b0c3eSLois Curfman McInnes ierr = SYWrite(fd,anonz,nz,SYSCALAR,0); CHKERRQ(ierr); 421*932b0c3eSLois Curfman McInnes PETSCFREE(anonz); 422*932b0c3eSLois Curfman McInnes return 0; 423*932b0c3eSLois Curfman McInnes } 424*932b0c3eSLois Curfman McInnes 425*932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer) 426*932b0c3eSLois Curfman McInnes { 427*932b0c3eSLois Curfman McInnes Mat A = (Mat) obj; 428*932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*) A->data; 429*932b0c3eSLois Curfman McInnes PetscObject vobj = (PetscObject) viewer; 430*932b0c3eSLois Curfman McInnes 431*932b0c3eSLois Curfman McInnes if (!viewer) { 432*932b0c3eSLois Curfman McInnes viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 433*932b0c3eSLois Curfman McInnes } 434*932b0c3eSLois Curfman McInnes if (vobj->cookie == VIEWER_COOKIE) { 435*932b0c3eSLois Curfman McInnes if (vobj->type == MATLAB_VIEWER) { 436*932b0c3eSLois Curfman McInnes return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); 437*932b0c3eSLois Curfman McInnes } 438*932b0c3eSLois Curfman McInnes else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) { 439*932b0c3eSLois Curfman McInnes return MatView_SeqDense_ASCII(A,viewer); 440*932b0c3eSLois Curfman McInnes } 441*932b0c3eSLois Curfman McInnes else if (vobj->type == BINARY_FILE_VIEWER) { 442*932b0c3eSLois Curfman McInnes return MatView_SeqDense_Binary(A,viewer); 443*932b0c3eSLois Curfman McInnes } 444*932b0c3eSLois Curfman McInnes } 445*932b0c3eSLois Curfman McInnes return 0; 446*932b0c3eSLois Curfman McInnes } 447289bc588SBarry Smith 448ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj) 449289bc588SBarry Smith { 450289bc588SBarry Smith Mat mat = (Mat) obj; 451ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 452a5a9c739SBarry Smith #if defined(PETSC_LOG) 453a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 454a5a9c739SBarry Smith #endif 45578b31e54SBarry Smith if (l->pivots) PETSCFREE(l->pivots); 45678b31e54SBarry Smith PETSCFREE(l); 457a5a9c739SBarry Smith PLogObjectDestroy(mat); 4589e25ed09SBarry Smith PETSCHEADERDESTROY(mat); 459289bc588SBarry Smith return 0; 460289bc588SBarry Smith } 461289bc588SBarry Smith 462c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout) 463289bc588SBarry Smith { 464c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 465d3e5ee88SLois Curfman McInnes int k, j, m, n; 466d3e5ee88SLois Curfman McInnes Scalar *v, tmp; 46748b35521SBarry Smith 468d3e5ee88SLois Curfman McInnes v = mat->v; m = mat->m; n = mat->n; 46948b35521SBarry Smith if (!matout) { /* in place transpose */ 47048d91487SBarry Smith if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Not for rectangular matrix in place"); 471d3e5ee88SLois Curfman McInnes for ( j=0; j<m; j++ ) { 472289bc588SBarry Smith for ( k=0; k<j; k++ ) { 473d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 474d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 475d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 476289bc588SBarry Smith } 477289bc588SBarry Smith } 47848b35521SBarry Smith } 479d3e5ee88SLois Curfman McInnes else { /* out-of-place transpose */ 480d3e5ee88SLois Curfman McInnes int ierr; 481d3e5ee88SLois Curfman McInnes Mat tmat; 482ec8511deSBarry Smith Mat_SeqDense *tmatd; 483d3e5ee88SLois Curfman McInnes Scalar *v2; 484c0bbcb79SLois Curfman McInnes ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,&tmat); CHKERRQ(ierr); 485ec8511deSBarry Smith tmatd = (Mat_SeqDense *) tmat->data; 4860de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 487d3e5ee88SLois Curfman McInnes for ( j=0; j<n; j++ ) { 4880de55854SLois Curfman McInnes for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 489d3e5ee88SLois Curfman McInnes } 490d3e5ee88SLois Curfman McInnes ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 491d3e5ee88SLois Curfman McInnes ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 492d3e5ee88SLois Curfman McInnes *matout = tmat; 49348b35521SBarry Smith } 494289bc588SBarry Smith return 0; 495289bc588SBarry Smith } 496289bc588SBarry Smith 497c0bbcb79SLois Curfman McInnes static int MatEqual_SeqDense(Mat A1,Mat A2) 498289bc588SBarry Smith { 499c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 500c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 501289bc588SBarry Smith int i; 502289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 503289bc588SBarry Smith if (mat1->m != mat2->m) return 0; 504289bc588SBarry Smith if (mat1->n != mat2->n) return 0; 505289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 506289bc588SBarry Smith if (*v1 != *v2) return 0; 507289bc588SBarry Smith v1++; v2++; 508289bc588SBarry Smith } 509289bc588SBarry Smith return 1; 510289bc588SBarry Smith } 511289bc588SBarry Smith 512c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v) 513289bc588SBarry Smith { 514c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 5156abc6512SBarry Smith int i, n; 5166abc6512SBarry Smith Scalar *x; 517289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 518ec8511deSBarry Smith if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 519289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 520289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 521289bc588SBarry Smith } 522289bc588SBarry Smith return 0; 523289bc588SBarry Smith } 524289bc588SBarry Smith 525c0bbcb79SLois Curfman McInnes static int MatScale_SeqDense(Mat A,Vec ll,Vec rr) 526289bc588SBarry Smith { 527c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 528da3a660dSBarry Smith Scalar *l,*r,x,*v; 529da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 53028988994SBarry Smith if (ll) { 531da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 532ec8511deSBarry Smith if (m != mat->m) SETERRQ(1,"MatScale_SeqDense:Left scaling vec wrong size"); 533da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 534da3a660dSBarry Smith x = l[i]; 535da3a660dSBarry Smith v = mat->v + i; 536da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 537da3a660dSBarry Smith } 538da3a660dSBarry Smith } 53928988994SBarry Smith if (rr) { 540da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 541ec8511deSBarry Smith if (n != mat->n) SETERRQ(1,"MatScale_SeqDense:Right scaling vec wrong size"); 542da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 543da3a660dSBarry Smith x = r[i]; 544da3a660dSBarry Smith v = mat->v + i*m; 545da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 546da3a660dSBarry Smith } 547da3a660dSBarry Smith } 548289bc588SBarry Smith return 0; 549289bc588SBarry Smith } 550289bc588SBarry Smith 551c0bbcb79SLois Curfman McInnes static int MatNorm_SeqDense(Mat A,MatNormType type,double *norm) 552289bc588SBarry Smith { 553c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 554289bc588SBarry Smith Scalar *v = mat->v; 555289bc588SBarry Smith double sum = 0.0; 556289bc588SBarry Smith int i, j; 557289bc588SBarry Smith if (type == NORM_FROBENIUS) { 558289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 559289bc588SBarry Smith #if defined(PETSC_COMPLEX) 560289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 561289bc588SBarry Smith #else 562289bc588SBarry Smith sum += (*v)*(*v); v++; 563289bc588SBarry Smith #endif 564289bc588SBarry Smith } 565289bc588SBarry Smith *norm = sqrt(sum); 566289bc588SBarry Smith } 567289bc588SBarry Smith else if (type == NORM_1) { 568289bc588SBarry Smith *norm = 0.0; 569289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 570289bc588SBarry Smith sum = 0.0; 571289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 572289bc588SBarry Smith #if defined(PETSC_COMPLEX) 573289bc588SBarry Smith sum += abs(*v++); 574289bc588SBarry Smith #else 575289bc588SBarry Smith sum += fabs(*v++); 576289bc588SBarry Smith #endif 577289bc588SBarry Smith } 578289bc588SBarry Smith if (sum > *norm) *norm = sum; 579289bc588SBarry Smith } 580289bc588SBarry Smith } 581289bc588SBarry Smith else if (type == NORM_INFINITY) { 582289bc588SBarry Smith *norm = 0.0; 583289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 584289bc588SBarry Smith v = mat->v + j; 585289bc588SBarry Smith sum = 0.0; 586289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 587289bc588SBarry Smith #if defined(PETSC_COMPLEX) 588289bc588SBarry Smith sum += abs(*v); v += mat->m; 589289bc588SBarry Smith #else 590289bc588SBarry Smith sum += fabs(*v); v += mat->m; 591289bc588SBarry Smith #endif 592289bc588SBarry Smith } 593289bc588SBarry Smith if (sum > *norm) *norm = sum; 594289bc588SBarry Smith } 595289bc588SBarry Smith } 596289bc588SBarry Smith else { 59748d91487SBarry Smith SETERRQ(1,"MatNorm_SeqDense:No two norm"); 598289bc588SBarry Smith } 599289bc588SBarry Smith return 0; 600289bc588SBarry Smith } 601289bc588SBarry Smith 602c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op) 603289bc588SBarry Smith { 604c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 605289bc588SBarry Smith if (op == ROW_ORIENTED) aij->roworiented = 1; 606289bc588SBarry Smith else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 607c0bbcb79SLois Curfman McInnes else if (op == ROWS_SORTED || 608c0bbcb79SLois Curfman McInnes op == SYMMETRIC_MATRIX || 609c0bbcb79SLois Curfman McInnes op == STRUCTURALLY_SYMMETRIC_MATRIX || 610c0bbcb79SLois Curfman McInnes op == NO_NEW_NONZERO_LOCATIONS || 611c0bbcb79SLois Curfman McInnes op == YES_NEW_NONZERO_LOCATIONS || 612c0bbcb79SLois Curfman McInnes op == NO_NEW_DIAGONALS || 613c0bbcb79SLois Curfman McInnes op == YES_NEW_DIAGONALS) 614c0bbcb79SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_SeqDense:Option ignored\n"); 615c0bbcb79SLois Curfman McInnes else 6164d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");} 617289bc588SBarry Smith return 0; 618289bc588SBarry Smith } 619289bc588SBarry Smith 620ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A) 6216f0a148fSBarry Smith { 622ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 623416022c9SBarry Smith PetscZero(l->v,l->m*l->n*sizeof(Scalar)); 6246f0a148fSBarry Smith return 0; 6256f0a148fSBarry Smith } 6266f0a148fSBarry Smith 627ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 6286f0a148fSBarry Smith { 629ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 6306abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 6316f0a148fSBarry Smith Scalar *slot; 63278b31e54SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 63378b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 6346f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 6356f0a148fSBarry Smith slot = l->v + rows[i]; 6366f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 6376f0a148fSBarry Smith } 6386f0a148fSBarry Smith if (diag) { 6396f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 6406f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 641260925b8SBarry Smith *slot = *diag; 6426f0a148fSBarry Smith } 6436f0a148fSBarry Smith } 644260925b8SBarry Smith ISRestoreIndices(is,&rows); 6456f0a148fSBarry Smith return 0; 6466f0a148fSBarry Smith } 647557bce09SLois Curfman McInnes 648c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n) 649557bce09SLois Curfman McInnes { 650c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 651557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 652557bce09SLois Curfman McInnes return 0; 653557bce09SLois Curfman McInnes } 654557bce09SLois Curfman McInnes 655c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 656d3e5ee88SLois Curfman McInnes { 657c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 658d3e5ee88SLois Curfman McInnes *m = 0; *n = mat->m; 659d3e5ee88SLois Curfman McInnes return 0; 660d3e5ee88SLois Curfman McInnes } 661d3e5ee88SLois Curfman McInnes 662c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array) 66364e87e97SBarry Smith { 664c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 66564e87e97SBarry Smith *array = mat->v; 66664e87e97SBarry Smith return 0; 66764e87e97SBarry Smith } 6680754003eSLois Curfman McInnes 6690754003eSLois Curfman McInnes 670c0bbcb79SLois Curfman McInnes static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol) 6710754003eSLois Curfman McInnes { 672ec8511deSBarry Smith SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done"); 6730754003eSLois Curfman McInnes } 6740754003eSLois Curfman McInnes 675c0bbcb79SLois Curfman McInnes static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,Mat *submat) 6760754003eSLois Curfman McInnes { 677c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 6780754003eSLois Curfman McInnes int nznew, *smap, i, j, ierr, oldcols = mat->n; 679160018dcSBarry Smith int *irow, *icol, nrows, ncols, *cwork; 6800754003eSLois Curfman McInnes Scalar *vwork, *val; 6810754003eSLois Curfman McInnes Mat newmat; 6820754003eSLois Curfman McInnes 68378b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 68478b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 68578b31e54SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 68678b31e54SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 6870754003eSLois Curfman McInnes 68878b31e54SBarry Smith smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap); 68978b31e54SBarry Smith cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork); 69078b31e54SBarry Smith vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 691416022c9SBarry Smith PetscZero((char*)smap,oldcols*sizeof(int)); 6920754003eSLois Curfman McInnes for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 6930754003eSLois Curfman McInnes 6940754003eSLois Curfman McInnes /* Create and fill new matrix */ 695c0bbcb79SLois Curfman McInnes ierr = MatCreateSeqDense(A->comm,nrows,ncols,&newmat);CHKERRQ(ierr); 6960754003eSLois Curfman McInnes for (i=0; i<nrows; i++) { 6970754003eSLois Curfman McInnes nznew = 0; 6980754003eSLois Curfman McInnes val = mat->v + irow[i]; 6990754003eSLois Curfman McInnes for (j=0; j<oldcols; j++) { 7000754003eSLois Curfman McInnes if (smap[j]) { 7010754003eSLois Curfman McInnes cwork[nznew] = smap[j] - 1; 7020754003eSLois Curfman McInnes vwork[nznew++] = val[j * mat->m]; 7030754003eSLois Curfman McInnes } 7040754003eSLois Curfman McInnes } 705dbb450caSBarry Smith ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES); 70678b31e54SBarry Smith CHKERRQ(ierr); 7070754003eSLois Curfman McInnes } 70878b31e54SBarry Smith ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 70978b31e54SBarry Smith ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 7100754003eSLois Curfman McInnes 7110754003eSLois Curfman McInnes /* Free work space */ 71278b31e54SBarry Smith PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork); 71378b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 71478b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 7150754003eSLois Curfman McInnes *submat = newmat; 7160754003eSLois Curfman McInnes return 0; 7170754003eSLois Curfman McInnes } 7180754003eSLois Curfman McInnes 719289bc588SBarry Smith /* -------------------------------------------------------------------*/ 720ec8511deSBarry Smith static struct _MatOps MatOps = {MatInsert_SeqDense, 721ec8511deSBarry Smith MatGetRow_SeqDense, MatRestoreRow_SeqDense, 722ec8511deSBarry Smith MatMult_SeqDense, MatMultAdd_SeqDense, 723ec8511deSBarry Smith MatMultTrans_SeqDense, MatMultTransAdd_SeqDense, 724ec8511deSBarry Smith MatSolve_SeqDense,MatSolveAdd_SeqDense, 725ec8511deSBarry Smith MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense, 726ec8511deSBarry Smith MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense, 727ec8511deSBarry Smith MatRelax_SeqDense, 728ec8511deSBarry Smith MatTranspose_SeqDense, 729ec8511deSBarry Smith MatGetInfo_SeqDense,MatEqual_SeqDense, 730ec8511deSBarry Smith MatGetDiagonal_SeqDense,MatScale_SeqDense,MatNorm_SeqDense, 731289bc588SBarry Smith 0,0, 732ec8511deSBarry Smith 0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0, 733ec8511deSBarry Smith MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense, 734ec8511deSBarry Smith MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense, 735ec8511deSBarry Smith MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense, 736ec8511deSBarry Smith 0,0,MatGetArray_SeqDense,0,0, 737ec8511deSBarry Smith MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense, 7381987afe7SBarry Smith MatCopyPrivate_SeqDense,0,0,0,0, 7391987afe7SBarry Smith MatAXPY_SeqDense}; 7400754003eSLois Curfman McInnes 7414b828684SBarry Smith /*@C 742fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 743d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 744d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 745289bc588SBarry Smith 74620563c6bSBarry Smith Input Parameters: 7470c775827SLois Curfman McInnes . comm - MPI communicator, set to MPI_COMM_SELF 7480c775827SLois Curfman McInnes . m - number of rows 7490c775827SLois Curfman McInnes . n - number of column 75020563c6bSBarry Smith 75120563c6bSBarry Smith Output Parameter: 7520c775827SLois Curfman McInnes . newmat - the matrix 75320563c6bSBarry Smith 754dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 755d65003e9SLois Curfman McInnes 756dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues() 75720563c6bSBarry Smith @*/ 758fafbff53SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,Mat *newmat) 759289bc588SBarry Smith { 760ec8511deSBarry Smith int size = sizeof(Mat_SeqDense) + m*n*sizeof(Scalar); 761289bc588SBarry Smith Mat mat; 762ec8511deSBarry Smith Mat_SeqDense *l; 763289bc588SBarry Smith *newmat = 0; 764ec8511deSBarry Smith PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm); 765a5a9c739SBarry Smith PLogObjectCreate(mat); 766ec8511deSBarry Smith l = (Mat_SeqDense *) PETSCMALLOC(size); CHKPTRQ(l); 767416022c9SBarry Smith PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 768ec8511deSBarry Smith mat->destroy = MatDestroy_SeqDense; 769ec8511deSBarry Smith mat->view = MatView_SeqDense; 770289bc588SBarry Smith mat->data = (void *) l; 771289bc588SBarry Smith mat->factor = 0; 772752f5567SLois Curfman McInnes PLogObjectMemory(mat,sizeof(struct _Mat) + size); 773d6dfbf8fSBarry Smith 774289bc588SBarry Smith l->m = m; 775289bc588SBarry Smith l->n = n; 776289bc588SBarry Smith l->v = (Scalar *) (l + 1); 777289bc588SBarry Smith l->pivots = 0; 778289bc588SBarry Smith l->roworiented = 1; 779d6dfbf8fSBarry Smith 780416022c9SBarry Smith PetscZero(l->v,m*n*sizeof(Scalar)); 781289bc588SBarry Smith *newmat = mat; 782289bc588SBarry Smith return 0; 783289bc588SBarry Smith } 784289bc588SBarry Smith 785c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat) 786289bc588SBarry Smith { 787c0bbcb79SLois Curfman McInnes Mat_SeqDense *m = (Mat_SeqDense *) A->data; 788c0bbcb79SLois Curfman McInnes return MatCreateSeqDense(A->comm,m->m,m->n,newmat); 789289bc588SBarry Smith } 790