1cb512458SBarry Smith #ifndef lint 2*9ea5d5aeSSatish Balay static char vcid[] = "$Id: dense.c,v 1.90 1996/02/01 18:52:28 curfman Exp balay $"; 3cb512458SBarry Smith #endif 467e560aaSBarry Smith /* 567e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 667e560aaSBarry Smith */ 7289bc588SBarry Smith 817699dbbSLois Curfman McInnes #include "dense.h" 9b16a3bb1SBarry Smith #include "pinclude/plapack.h" 10b16a3bb1SBarry Smith #include "pinclude/pviewer.h" 11289bc588SBarry Smith 121987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y) 131987afe7SBarry Smith { 141987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data; 151987afe7SBarry Smith int N = x->m*x->n, one = 1; 161987afe7SBarry Smith BLaxpy_( &N, alpha, x->v, &one, y->v, &one ); 171987afe7SBarry Smith PLogFlops(2*N-1); 181987afe7SBarry Smith return 0; 191987afe7SBarry Smith } 201987afe7SBarry Smith 21c0bbcb79SLois Curfman McInnes static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 22289bc588SBarry Smith { 23c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 24289bc588SBarry Smith int i,N = mat->m*mat->n,count = 0; 25289bc588SBarry Smith Scalar *v = mat->v; 26289bc588SBarry Smith for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 27c0bbcb79SLois Curfman McInnes *nz = count; *nzalloc = N; *mem = (int)A->mem; 28fa9ec3f1SLois Curfman McInnes return 0; 29289bc588SBarry Smith } 30289bc588SBarry Smith 31289bc588SBarry Smith /* ---------------------------------------------------------------*/ 32289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots, 33289bc588SBarry Smith rather than put it in the Mat->col slot.*/ 34c0bbcb79SLois Curfman McInnes static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f) 35289bc588SBarry Smith { 36c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 376abc6512SBarry Smith int info; 3855659b69SBarry Smith 39289bc588SBarry Smith if (!mat->pivots) { 400452661fSBarry Smith mat->pivots = (int *) PetscMalloc(mat->m*sizeof(int));CHKPTRQ(mat->pivots); 41c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,mat->m*sizeof(int)); 42289bc588SBarry Smith } 43289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 44ec8511deSBarry Smith if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization"); 45c0bbcb79SLois Curfman McInnes A->factor = FACTOR_LU; 4655659b69SBarry Smith PLogFlops((2*mat->n*mat->n*mat->n)/3); 47289bc588SBarry Smith return 0; 48289bc588SBarry Smith } 49c0bbcb79SLois Curfman McInnes static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact) 50289bc588SBarry Smith { 51a57ad546SLois Curfman McInnes return MatConvert(A,MATSAME,fact); 52289bc588SBarry Smith } 53c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 54289bc588SBarry Smith { 5549d8b64dSBarry Smith return MatLUFactor(*fact,0,0,1.0); 56289bc588SBarry Smith } 57c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact) 58289bc588SBarry Smith { 59a57ad546SLois Curfman McInnes return MatConvert(A,MATSAME,fact); 60289bc588SBarry Smith } 61c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 62289bc588SBarry Smith { 6349d8b64dSBarry Smith return MatCholeskyFactor(*fact,0,1.0); 64289bc588SBarry Smith } 65c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f) 66289bc588SBarry Smith { 67c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 686abc6512SBarry Smith int info; 6955659b69SBarry Smith 70752f5567SLois Curfman McInnes if (mat->pivots) { 710452661fSBarry Smith PetscFree(mat->pivots); 72c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,-mat->m*sizeof(int)); 73752f5567SLois Curfman McInnes mat->pivots = 0; 74752f5567SLois Curfman McInnes } 75289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 76ec8511deSBarry Smith if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization"); 77c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 7855659b69SBarry Smith PLogFlops((mat->n*mat->n*mat->n)/3); 79289bc588SBarry Smith return 0; 80289bc588SBarry Smith } 81289bc588SBarry Smith 82c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 83289bc588SBarry Smith { 84c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 85a57ad546SLois Curfman McInnes int one = 1, info, ierr; 866abc6512SBarry Smith Scalar *x, *y; 8767e560aaSBarry Smith 88a57ad546SLois Curfman McInnes ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 89416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 90c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 9148d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 92289bc588SBarry Smith } 93c0bbcb79SLois Curfman McInnes else if (A->factor == FACTOR_CHOLESKY){ 9448d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 95289bc588SBarry Smith } 96ec8511deSBarry Smith else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve"); 97ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve"); 9855659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 99289bc588SBarry Smith return 0; 100289bc588SBarry Smith } 101c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy) 102da3a660dSBarry Smith { 103c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1046abc6512SBarry Smith int one = 1, info; 1056abc6512SBarry Smith Scalar *x, *y; 10667e560aaSBarry Smith 107da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 108416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 109752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 110da3a660dSBarry Smith if (mat->pivots) { 11148d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 112da3a660dSBarry Smith } 113da3a660dSBarry Smith else { 11448d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 115da3a660dSBarry Smith } 116ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve"); 11755659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 118da3a660dSBarry Smith return 0; 119da3a660dSBarry Smith } 120c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 121da3a660dSBarry Smith { 122c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1236abc6512SBarry Smith int one = 1, info,ierr; 1246abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 125da3a660dSBarry Smith Vec tmp = 0; 12667e560aaSBarry Smith 127da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 128da3a660dSBarry Smith if (yy == zz) { 12978b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 130c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 13178b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 132da3a660dSBarry Smith } 133416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 134752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 135da3a660dSBarry Smith if (mat->pivots) { 13648d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 137da3a660dSBarry Smith } 138da3a660dSBarry Smith else { 13948d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 140da3a660dSBarry Smith } 141ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve"); 142da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 143da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 14455659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 145da3a660dSBarry Smith return 0; 146da3a660dSBarry Smith } 14767e560aaSBarry Smith 148c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 149da3a660dSBarry Smith { 150c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1516abc6512SBarry Smith int one = 1, info,ierr; 1526abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 153da3a660dSBarry Smith Vec tmp; 15467e560aaSBarry Smith 155da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 156da3a660dSBarry Smith if (yy == zz) { 15778b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 158c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 15978b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 160da3a660dSBarry Smith } 161416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 162752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 163da3a660dSBarry Smith if (mat->pivots) { 16448d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 165da3a660dSBarry Smith } 166da3a660dSBarry Smith else { 16748d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 168da3a660dSBarry Smith } 169ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve"); 170da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 171da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 17255659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 173da3a660dSBarry Smith return 0; 174da3a660dSBarry Smith } 175289bc588SBarry Smith /* ------------------------------------------------------------------*/ 176c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 17720e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 178289bc588SBarry Smith { 179c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 180289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 1816abc6512SBarry Smith int o = 1,ierr, m = mat->m, i; 182289bc588SBarry Smith 183289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 184289bc588SBarry Smith /* this is a hack fix, should have another version without 185289bc588SBarry Smith the second BLdot */ 186bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx); CHKERRQ(ierr); 187289bc588SBarry Smith } 188289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 189289bc588SBarry Smith while (its--) { 190289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 191289bc588SBarry Smith for ( i=0; i<m; i++ ) { 192f1747703SBarry Smith #if defined(PETSC_COMPLEX) 193f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 194f1747703SBarry Smith not happy about returning a double complex */ 195f1747703SBarry Smith int _i; 196f1747703SBarry Smith Scalar sum = b[i]; 197f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 198f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 199f1747703SBarry Smith } 200f1747703SBarry Smith xt = sum; 201f1747703SBarry Smith #else 202289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 203f1747703SBarry Smith #endif 204d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 205289bc588SBarry Smith } 206289bc588SBarry Smith } 207289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 208289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 209f1747703SBarry Smith #if defined(PETSC_COMPLEX) 210f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 211f1747703SBarry Smith not happy about returning a double complex */ 212f1747703SBarry Smith int _i; 213f1747703SBarry Smith Scalar sum = b[i]; 214f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 215f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 216f1747703SBarry Smith } 217f1747703SBarry Smith xt = sum; 218f1747703SBarry Smith #else 219289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 220f1747703SBarry Smith #endif 221d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 222289bc588SBarry Smith } 223289bc588SBarry Smith } 224289bc588SBarry Smith } 225289bc588SBarry Smith return 0; 226289bc588SBarry Smith } 227289bc588SBarry Smith 228289bc588SBarry Smith /* -----------------------------------------------------------------*/ 229c0bbcb79SLois Curfman McInnes static int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy) 230289bc588SBarry Smith { 231c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 232289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 233289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 234289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 23548d91487SBarry Smith LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 23655659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->n); 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); 24655659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->m); 247289bc588SBarry Smith return 0; 248289bc588SBarry Smith } 249c0bbcb79SLois Curfman McInnes static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 250289bc588SBarry Smith { 251c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 252289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 2536abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 254289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 255416022c9SBarry Smith if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 25648d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 25755659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 258289bc588SBarry Smith return 0; 259289bc588SBarry Smith } 260c0bbcb79SLois Curfman McInnes static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 261289bc588SBarry Smith { 262c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 263289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 264289bc588SBarry Smith int _One=1; 2656abc6512SBarry Smith Scalar _DOne=1.0; 266289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 267289bc588SBarry Smith VecGetArray(zz,&z); 268416022c9SBarry Smith if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 26948d91487SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 27055659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 271289bc588SBarry Smith return 0; 272289bc588SBarry Smith } 273289bc588SBarry Smith 274289bc588SBarry Smith /* -----------------------------------------------------------------*/ 275c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 276289bc588SBarry Smith { 277c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 278289bc588SBarry Smith Scalar *v; 279289bc588SBarry Smith int i; 28067e560aaSBarry Smith 281289bc588SBarry Smith *ncols = mat->n; 282289bc588SBarry Smith if (cols) { 2830452661fSBarry Smith *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols); 28473c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 285289bc588SBarry Smith } 286289bc588SBarry Smith if (vals) { 2870452661fSBarry Smith *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 288289bc588SBarry Smith v = mat->v + row; 28973c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 290289bc588SBarry Smith } 291289bc588SBarry Smith return 0; 292289bc588SBarry Smith } 293c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 294289bc588SBarry Smith { 2950452661fSBarry Smith if (cols) { PetscFree(*cols); } 2960452661fSBarry Smith if (vals) { PetscFree(*vals); } 297289bc588SBarry Smith return 0; 298289bc588SBarry Smith } 299289bc588SBarry Smith /* ----------------------------------------------------------------*/ 300ae80bb75SLois Curfman McInnes static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 301e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 302289bc588SBarry Smith { 303c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 304289bc588SBarry Smith int i,j; 305d6dfbf8fSBarry Smith 306289bc588SBarry Smith if (!mat->roworiented) { 307dbb450caSBarry Smith if (addv == INSERT_VALUES) { 308289bc588SBarry Smith for ( j=0; j<n; j++ ) { 309289bc588SBarry Smith for ( i=0; i<m; i++ ) { 310289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 311289bc588SBarry Smith } 312289bc588SBarry Smith } 313289bc588SBarry Smith } 314289bc588SBarry Smith else { 315289bc588SBarry Smith for ( j=0; j<n; j++ ) { 316289bc588SBarry Smith for ( i=0; i<m; i++ ) { 317289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 318289bc588SBarry Smith } 319289bc588SBarry Smith } 320289bc588SBarry Smith } 321e8d4e0b9SBarry Smith } 322e8d4e0b9SBarry Smith else { 323dbb450caSBarry Smith if (addv == INSERT_VALUES) { 324e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 325e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 326e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 327e8d4e0b9SBarry Smith } 328e8d4e0b9SBarry Smith } 329e8d4e0b9SBarry Smith } 330289bc588SBarry Smith else { 331289bc588SBarry Smith for ( i=0; i<m; i++ ) { 332289bc588SBarry Smith for ( j=0; j<n; j++ ) { 333289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 334289bc588SBarry Smith } 335289bc588SBarry Smith } 336289bc588SBarry Smith } 337e8d4e0b9SBarry Smith } 338289bc588SBarry Smith return 0; 339289bc588SBarry Smith } 340e8d4e0b9SBarry Smith 341ae80bb75SLois Curfman McInnes static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v) 342ae80bb75SLois Curfman McInnes { 343ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 344ae80bb75SLois Curfman McInnes int i, j; 345ae80bb75SLois Curfman McInnes Scalar *vpt = v; 346ae80bb75SLois Curfman McInnes 347ae80bb75SLois Curfman McInnes /* row-oriented output */ 348ae80bb75SLois Curfman McInnes for ( i=0; i<m; i++ ) { 349ae80bb75SLois Curfman McInnes for ( j=0; j<n; j++ ) { 350ae80bb75SLois Curfman McInnes *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]]; 351ae80bb75SLois Curfman McInnes } 352ae80bb75SLois Curfman McInnes } 353ae80bb75SLois Curfman McInnes return 0; 354ae80bb75SLois Curfman McInnes } 355ae80bb75SLois Curfman McInnes 356289bc588SBarry Smith /* -----------------------------------------------------------------*/ 3573d1612f7SBarry Smith static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues) 358289bc588SBarry Smith { 359c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l; 360289bc588SBarry Smith int ierr; 361289bc588SBarry Smith Mat newi; 36248d91487SBarry Smith 363b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr); 364ec8511deSBarry Smith l = (Mat_SeqDense *) newi->data; 36555659b69SBarry Smith if (cpvalues == COPY_VALUES) { 366416022c9SBarry Smith PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 36755659b69SBarry Smith } 368227d817aSBarry Smith newi->assembled = PETSC_TRUE; 369289bc588SBarry Smith *newmat = newi; 370289bc588SBarry Smith return 0; 371289bc588SBarry Smith } 372289bc588SBarry Smith 37388e32edaSLois Curfman McInnes #include "sysio.h" 37488e32edaSLois Curfman McInnes 37588e32edaSLois Curfman McInnes int MatLoad_SeqDense(Viewer bview,MatType type,Mat *A) 37688e32edaSLois Curfman McInnes { 37788e32edaSLois Curfman McInnes Mat_SeqDense *a; 37888e32edaSLois Curfman McInnes Mat B; 37988e32edaSLois Curfman McInnes int *scols, i, j, nz, ierr, fd, header[4], size; 38088e32edaSLois Curfman McInnes int *rowlengths = 0, M, N, *cols; 38188e32edaSLois Curfman McInnes Scalar *vals, *svals, *v; 38288e32edaSLois Curfman McInnes PetscObject vobj = (PetscObject) bview; 38388e32edaSLois Curfman McInnes MPI_Comm comm = vobj->comm; 38488e32edaSLois Curfman McInnes 38588e32edaSLois Curfman McInnes MPI_Comm_size(comm,&size); 38688e32edaSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor"); 38788e32edaSLois Curfman McInnes ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 38888e32edaSLois Curfman McInnes ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 38988e32edaSLois Curfman McInnes if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object"); 39088e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 39188e32edaSLois Curfman McInnes 39288e32edaSLois Curfman McInnes /* read row lengths */ 3930452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 39488e32edaSLois Curfman McInnes ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 39588e32edaSLois Curfman McInnes 39688e32edaSLois Curfman McInnes /* create our matrix */ 397b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 39888e32edaSLois Curfman McInnes B = *A; 39988e32edaSLois Curfman McInnes a = (Mat_SeqDense *) B->data; 40088e32edaSLois Curfman McInnes v = a->v; 40188e32edaSLois Curfman McInnes 40288e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 4030452661fSBarry Smith cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols); 40488e32edaSLois Curfman McInnes ierr = SYRead(fd,cols,nz,SYINT); CHKERRQ(ierr); 4050452661fSBarry Smith vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 40688e32edaSLois Curfman McInnes ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr); 40788e32edaSLois Curfman McInnes 40888e32edaSLois Curfman McInnes /* insert into matrix */ 40988e32edaSLois Curfman McInnes for ( i=0; i<M; i++ ) { 41088e32edaSLois Curfman McInnes for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j]; 41188e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 41288e32edaSLois Curfman McInnes } 4130452661fSBarry Smith PetscFree(vals); PetscFree(cols); PetscFree(rowlengths); 41488e32edaSLois Curfman McInnes 41588e32edaSLois Curfman McInnes ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 41688e32edaSLois Curfman McInnes ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 41788e32edaSLois Curfman McInnes return 0; 41888e32edaSLois Curfman McInnes } 41988e32edaSLois Curfman McInnes 420932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h" 421932b0c3eSLois Curfman McInnes #include "sysio.h" 422932b0c3eSLois Curfman McInnes 423932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 424289bc588SBarry Smith { 425932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 426932b0c3eSLois Curfman McInnes int ierr, i, j, format; 427d636dbe3SBarry Smith FILE *fd; 428932b0c3eSLois Curfman McInnes char *outputname; 429932b0c3eSLois Curfman McInnes Scalar *v; 430932b0c3eSLois Curfman McInnes 431227d817aSBarry Smith ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr); 432932b0c3eSLois Curfman McInnes ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 433932b0c3eSLois Curfman McInnes ierr = ViewerFileGetFormat_Private(viewer,&format); 434f72585f2SLois Curfman McInnes if (format == FILE_FORMAT_INFO) { 435932b0c3eSLois Curfman McInnes ; /* do nothing for now */ 436f72585f2SLois Curfman McInnes } 437f72585f2SLois Curfman McInnes else { 438932b0c3eSLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 439932b0c3eSLois Curfman McInnes v = a->v + i; 440932b0c3eSLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 441289bc588SBarry Smith #if defined(PETSC_COMPLEX) 442932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m; 443289bc588SBarry Smith #else 444932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e ",*v); v += a->m; 445289bc588SBarry Smith #endif 446289bc588SBarry Smith } 4478e37dc09SBarry Smith fprintf(fd,"\n"); 448289bc588SBarry Smith } 449da3a660dSBarry Smith } 450932b0c3eSLois Curfman McInnes fflush(fd); 451289bc588SBarry Smith return 0; 452289bc588SBarry Smith } 453289bc588SBarry Smith 454932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 455932b0c3eSLois Curfman McInnes { 456932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 457932b0c3eSLois Curfman McInnes int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 458932b0c3eSLois Curfman McInnes Scalar *v, *anonz; 459932b0c3eSLois Curfman McInnes 460932b0c3eSLois Curfman McInnes ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 4610452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens); 462932b0c3eSLois Curfman McInnes col_lens[0] = MAT_COOKIE; 463932b0c3eSLois Curfman McInnes col_lens[1] = m; 464932b0c3eSLois Curfman McInnes col_lens[2] = n; 465932b0c3eSLois Curfman McInnes col_lens[3] = nz; 466932b0c3eSLois Curfman McInnes 467932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 468932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) col_lens[4+i] = n; 469932b0c3eSLois Curfman McInnes ierr = SYWrite(fd,col_lens,4+m,SYINT,1); CHKERRQ(ierr); 470932b0c3eSLois Curfman McInnes 471932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 472932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 473932b0c3eSLois Curfman McInnes ict = 0; 474932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 475932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) col_lens[ict++] = j; 476932b0c3eSLois Curfman McInnes } 477932b0c3eSLois Curfman McInnes ierr = SYWrite(fd,col_lens,nz,SYINT,0); CHKERRQ(ierr); 4780452661fSBarry Smith PetscFree(col_lens); 479932b0c3eSLois Curfman McInnes 480932b0c3eSLois Curfman McInnes /* store nonzero values */ 4810452661fSBarry Smith anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz); 482932b0c3eSLois Curfman McInnes ict = 0; 483932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 484932b0c3eSLois Curfman McInnes v = a->v + i; 485932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) { 486932b0c3eSLois Curfman McInnes anonz[ict++] = *v; v += a->m; 487932b0c3eSLois Curfman McInnes } 488932b0c3eSLois Curfman McInnes } 489932b0c3eSLois Curfman McInnes ierr = SYWrite(fd,anonz,nz,SYSCALAR,0); CHKERRQ(ierr); 4900452661fSBarry Smith PetscFree(anonz); 491932b0c3eSLois Curfman McInnes return 0; 492932b0c3eSLois Curfman McInnes } 493932b0c3eSLois Curfman McInnes 494932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer) 495932b0c3eSLois Curfman McInnes { 496932b0c3eSLois Curfman McInnes Mat A = (Mat) obj; 497932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*) A->data; 498932b0c3eSLois Curfman McInnes PetscObject vobj = (PetscObject) viewer; 499932b0c3eSLois Curfman McInnes 500932b0c3eSLois Curfman McInnes if (!viewer) { 501932b0c3eSLois Curfman McInnes viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 502932b0c3eSLois Curfman McInnes } 503932b0c3eSLois Curfman McInnes if (vobj->cookie == VIEWER_COOKIE) { 504932b0c3eSLois Curfman McInnes if (vobj->type == MATLAB_VIEWER) { 505932b0c3eSLois Curfman McInnes return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); 506932b0c3eSLois Curfman McInnes } 507932b0c3eSLois Curfman McInnes else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) { 508932b0c3eSLois Curfman McInnes return MatView_SeqDense_ASCII(A,viewer); 509932b0c3eSLois Curfman McInnes } 510932b0c3eSLois Curfman McInnes else if (vobj->type == BINARY_FILE_VIEWER) { 511932b0c3eSLois Curfman McInnes return MatView_SeqDense_Binary(A,viewer); 512932b0c3eSLois Curfman McInnes } 513932b0c3eSLois Curfman McInnes } 514932b0c3eSLois Curfman McInnes return 0; 515932b0c3eSLois Curfman McInnes } 516289bc588SBarry Smith 517ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj) 518289bc588SBarry Smith { 519289bc588SBarry Smith Mat mat = (Mat) obj; 520ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 521a5a9c739SBarry Smith #if defined(PETSC_LOG) 522a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 523a5a9c739SBarry Smith #endif 5240452661fSBarry Smith if (l->pivots) PetscFree(l->pivots); 5253439631bSBarry Smith if (!l->user_alloc) PetscFree(l->v); 5260452661fSBarry Smith PetscFree(l); 527a5a9c739SBarry Smith PLogObjectDestroy(mat); 5280452661fSBarry Smith PetscHeaderDestroy(mat); 529289bc588SBarry Smith return 0; 530289bc588SBarry Smith } 531289bc588SBarry Smith 532c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout) 533289bc588SBarry Smith { 534c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 535d3e5ee88SLois Curfman McInnes int k, j, m, n; 536d3e5ee88SLois Curfman McInnes Scalar *v, tmp; 53748b35521SBarry Smith 538d3e5ee88SLois Curfman McInnes v = mat->v; m = mat->m; n = mat->n; 5393638b69dSLois Curfman McInnes if (matout == PETSC_NULL) { /* in place transpose */ 540d9f96c7cSLois Curfman McInnes if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place"); 541d3e5ee88SLois Curfman McInnes for ( j=0; j<m; j++ ) { 542289bc588SBarry Smith for ( k=0; k<j; k++ ) { 543d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 544d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 545d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 546289bc588SBarry Smith } 547289bc588SBarry Smith } 54848b35521SBarry Smith } 549d3e5ee88SLois Curfman McInnes else { /* out-of-place transpose */ 550d3e5ee88SLois Curfman McInnes int ierr; 551d3e5ee88SLois Curfman McInnes Mat tmat; 552ec8511deSBarry Smith Mat_SeqDense *tmatd; 553d3e5ee88SLois Curfman McInnes Scalar *v2; 554b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr); 555ec8511deSBarry Smith tmatd = (Mat_SeqDense *) tmat->data; 5560de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 557d3e5ee88SLois Curfman McInnes for ( j=0; j<n; j++ ) { 5580de55854SLois Curfman McInnes for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 559d3e5ee88SLois Curfman McInnes } 560d3e5ee88SLois Curfman McInnes ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 561d3e5ee88SLois Curfman McInnes ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 562d3e5ee88SLois Curfman McInnes *matout = tmat; 56348b35521SBarry Smith } 564289bc588SBarry Smith return 0; 565289bc588SBarry Smith } 566289bc588SBarry Smith 567*9ea5d5aeSSatish Balay static int MatEqual_SeqDense(Mat A1,Mat A2, int *flg) 568289bc588SBarry Smith { 569c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 570c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 571289bc588SBarry Smith int i; 572289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 573*9ea5d5aeSSatish Balay 574*9ea5d5aeSSatish Balay if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type MATSEQDENSE"); 575*9ea5d5aeSSatish Balay if (mat1->m != mat2->m) { *flg = 0; return 0;} 576*9ea5d5aeSSatish Balay if (mat1->n != mat2->n) {*flg =0; return 0;} 577289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 578*9ea5d5aeSSatish Balay if (*v1 != *v2) {*flg =0; return 0;} 579289bc588SBarry Smith v1++; v2++; 580289bc588SBarry Smith } 581*9ea5d5aeSSatish Balay *flg = 1; 582*9ea5d5aeSSatish Balay return 0; 583289bc588SBarry Smith } 584289bc588SBarry Smith 585c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v) 586289bc588SBarry Smith { 587c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 5886abc6512SBarry Smith int i, n; 5896abc6512SBarry Smith Scalar *x; 590289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 591ec8511deSBarry Smith if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 592289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 593289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 594289bc588SBarry Smith } 595289bc588SBarry Smith return 0; 596289bc588SBarry Smith } 597289bc588SBarry Smith 598052efed2SBarry Smith static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 599289bc588SBarry Smith { 600c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 601da3a660dSBarry Smith Scalar *l,*r,x,*v; 602da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 60355659b69SBarry Smith 60428988994SBarry Smith if (ll) { 605da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 606052efed2SBarry Smith if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size"); 60755659b69SBarry Smith PLogFlops(n*m); 608da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 609da3a660dSBarry Smith x = l[i]; 610da3a660dSBarry Smith v = mat->v + i; 611da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 612da3a660dSBarry Smith } 613da3a660dSBarry Smith } 61428988994SBarry Smith if (rr) { 615da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 616052efed2SBarry Smith if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size"); 61755659b69SBarry Smith PLogFlops(n*m); 618da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 619da3a660dSBarry Smith x = r[i]; 620da3a660dSBarry Smith v = mat->v + i*m; 621da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 622da3a660dSBarry Smith } 623da3a660dSBarry Smith } 624289bc588SBarry Smith return 0; 625289bc588SBarry Smith } 626289bc588SBarry Smith 627cddf8d76SBarry Smith static int MatNorm_SeqDense(Mat A,NormType type,double *norm) 628289bc588SBarry Smith { 629c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 630289bc588SBarry Smith Scalar *v = mat->v; 631289bc588SBarry Smith double sum = 0.0; 632289bc588SBarry Smith int i, j; 63355659b69SBarry Smith 634289bc588SBarry Smith if (type == NORM_FROBENIUS) { 635289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 636289bc588SBarry Smith #if defined(PETSC_COMPLEX) 637289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 638289bc588SBarry Smith #else 639289bc588SBarry Smith sum += (*v)*(*v); v++; 640289bc588SBarry Smith #endif 641289bc588SBarry Smith } 642289bc588SBarry Smith *norm = sqrt(sum); 64355659b69SBarry Smith PLogFlops(2*mat->n*mat->m); 644289bc588SBarry Smith } 645289bc588SBarry Smith else if (type == NORM_1) { 646289bc588SBarry Smith *norm = 0.0; 647289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 648289bc588SBarry Smith sum = 0.0; 649289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 65033a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 651289bc588SBarry Smith } 652289bc588SBarry Smith if (sum > *norm) *norm = sum; 653289bc588SBarry Smith } 65455659b69SBarry Smith PLogFlops(mat->n*mat->m); 655289bc588SBarry Smith } 656289bc588SBarry Smith else if (type == NORM_INFINITY) { 657289bc588SBarry Smith *norm = 0.0; 658289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 659289bc588SBarry Smith v = mat->v + j; 660289bc588SBarry Smith sum = 0.0; 661289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 662cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v += mat->m; 663289bc588SBarry Smith } 664289bc588SBarry Smith if (sum > *norm) *norm = sum; 665289bc588SBarry Smith } 66655659b69SBarry Smith PLogFlops(mat->n*mat->m); 667289bc588SBarry Smith } 668289bc588SBarry Smith else { 66948d91487SBarry Smith SETERRQ(1,"MatNorm_SeqDense:No two norm"); 670289bc588SBarry Smith } 671289bc588SBarry Smith return 0; 672289bc588SBarry Smith } 673289bc588SBarry Smith 674c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op) 675289bc588SBarry Smith { 676c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 67767e560aaSBarry Smith 678289bc588SBarry Smith if (op == ROW_ORIENTED) aij->roworiented = 1; 679289bc588SBarry Smith else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 680c0bbcb79SLois Curfman McInnes else if (op == ROWS_SORTED || 68188e32edaSLois Curfman McInnes op == COLUMNS_SORTED || 682c0bbcb79SLois Curfman McInnes op == SYMMETRIC_MATRIX || 683c0bbcb79SLois Curfman McInnes op == STRUCTURALLY_SYMMETRIC_MATRIX || 684c0bbcb79SLois Curfman McInnes op == NO_NEW_NONZERO_LOCATIONS || 685c0bbcb79SLois Curfman McInnes op == YES_NEW_NONZERO_LOCATIONS || 686c0bbcb79SLois Curfman McInnes op == NO_NEW_DIAGONALS || 687c0bbcb79SLois Curfman McInnes op == YES_NEW_DIAGONALS) 688c0bbcb79SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_SeqDense:Option ignored\n"); 689c0bbcb79SLois Curfman McInnes else 6904d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");} 691289bc588SBarry Smith return 0; 692289bc588SBarry Smith } 693289bc588SBarry Smith 694ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A) 6956f0a148fSBarry Smith { 696ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 697cddf8d76SBarry Smith PetscMemzero(l->v,l->m*l->n*sizeof(Scalar)); 6986f0a148fSBarry Smith return 0; 6996f0a148fSBarry Smith } 7006f0a148fSBarry Smith 701ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 7026f0a148fSBarry Smith { 703ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 7046abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 7056f0a148fSBarry Smith Scalar *slot; 70655659b69SBarry Smith 70778b31e54SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 70878b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 7096f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 7106f0a148fSBarry Smith slot = l->v + rows[i]; 7116f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 7126f0a148fSBarry Smith } 7136f0a148fSBarry Smith if (diag) { 7146f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 7156f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 716260925b8SBarry Smith *slot = *diag; 7176f0a148fSBarry Smith } 7186f0a148fSBarry Smith } 719260925b8SBarry Smith ISRestoreIndices(is,&rows); 7206f0a148fSBarry Smith return 0; 7216f0a148fSBarry Smith } 722557bce09SLois Curfman McInnes 723c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n) 724557bce09SLois Curfman McInnes { 725c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 726557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 727557bce09SLois Curfman McInnes return 0; 728557bce09SLois Curfman McInnes } 729557bce09SLois Curfman McInnes 730c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 731d3e5ee88SLois Curfman McInnes { 732c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 733d3e5ee88SLois Curfman McInnes *m = 0; *n = mat->m; 734d3e5ee88SLois Curfman McInnes return 0; 735d3e5ee88SLois Curfman McInnes } 736d3e5ee88SLois Curfman McInnes 737c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array) 73864e87e97SBarry Smith { 739c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 74064e87e97SBarry Smith *array = mat->v; 74164e87e97SBarry Smith return 0; 74264e87e97SBarry Smith } 7430754003eSLois Curfman McInnes 7440754003eSLois Curfman McInnes 745c0bbcb79SLois Curfman McInnes static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol) 7460754003eSLois Curfman McInnes { 747ec8511deSBarry Smith SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done"); 7480754003eSLois Curfman McInnes } 7490754003eSLois Curfman McInnes 750cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall, 751cddf8d76SBarry Smith Mat *submat) 7520754003eSLois Curfman McInnes { 753c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 7540754003eSLois Curfman McInnes int nznew, *smap, i, j, ierr, oldcols = mat->n; 755160018dcSBarry Smith int *irow, *icol, nrows, ncols, *cwork; 7560754003eSLois Curfman McInnes Scalar *vwork, *val; 7570754003eSLois Curfman McInnes Mat newmat; 7580754003eSLois Curfman McInnes 75978b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 76078b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 76178b31e54SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 76278b31e54SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 7630754003eSLois Curfman McInnes 7640452661fSBarry Smith smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap); 7650452661fSBarry Smith cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork); 7660452661fSBarry Smith vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 767cddf8d76SBarry Smith PetscMemzero((char*)smap,oldcols*sizeof(int)); 7680754003eSLois Curfman McInnes for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 7690754003eSLois Curfman McInnes 7700754003eSLois Curfman McInnes /* Create and fill new matrix */ 771b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr); 7720754003eSLois Curfman McInnes for (i=0; i<nrows; i++) { 7730754003eSLois Curfman McInnes nznew = 0; 7740754003eSLois Curfman McInnes val = mat->v + irow[i]; 7750754003eSLois Curfman McInnes for (j=0; j<oldcols; j++) { 7760754003eSLois Curfman McInnes if (smap[j]) { 7770754003eSLois Curfman McInnes cwork[nznew] = smap[j] - 1; 7780754003eSLois Curfman McInnes vwork[nznew++] = val[j * mat->m]; 7790754003eSLois Curfman McInnes } 7800754003eSLois Curfman McInnes } 781dbb450caSBarry Smith ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES); 78278b31e54SBarry Smith CHKERRQ(ierr); 7830754003eSLois Curfman McInnes } 78478b31e54SBarry Smith ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 78578b31e54SBarry Smith ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 7860754003eSLois Curfman McInnes 7870754003eSLois Curfman McInnes /* Free work space */ 7880452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 78978b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 79078b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 7910754003eSLois Curfman McInnes *submat = newmat; 7920754003eSLois Curfman McInnes return 0; 7930754003eSLois Curfman McInnes } 7940754003eSLois Curfman McInnes 7954b0e389bSBarry Smith static int MatCopy_SeqDense(Mat A, Mat B) 7964b0e389bSBarry Smith { 7974b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data; 7984b0e389bSBarry Smith if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B); 7994b0e389bSBarry Smith if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)"); 8004b0e389bSBarry Smith PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar)); 8014b0e389bSBarry Smith return 0; 8024b0e389bSBarry Smith } 8034b0e389bSBarry Smith 804289bc588SBarry Smith /* -------------------------------------------------------------------*/ 805ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense, 806ec8511deSBarry Smith MatGetRow_SeqDense, MatRestoreRow_SeqDense, 807ec8511deSBarry Smith MatMult_SeqDense, MatMultAdd_SeqDense, 808ec8511deSBarry Smith MatMultTrans_SeqDense, MatMultTransAdd_SeqDense, 809ec8511deSBarry Smith MatSolve_SeqDense,MatSolveAdd_SeqDense, 810ec8511deSBarry Smith MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense, 811ec8511deSBarry Smith MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense, 812ec8511deSBarry Smith MatRelax_SeqDense, 813ec8511deSBarry Smith MatTranspose_SeqDense, 814ec8511deSBarry Smith MatGetInfo_SeqDense,MatEqual_SeqDense, 815052efed2SBarry Smith MatGetDiagonal_SeqDense,MatDiagonalScale_SeqDense,MatNorm_SeqDense, 816289bc588SBarry Smith 0,0, 817ec8511deSBarry Smith 0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0, 818ec8511deSBarry Smith MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense, 819ec8511deSBarry Smith MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense, 820ec8511deSBarry Smith MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense, 821ec8511deSBarry Smith 0,0,MatGetArray_SeqDense,0,0, 822ec8511deSBarry Smith MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense, 8233d1612f7SBarry Smith MatConvertSameType_SeqDense,0,0,0,0, 824ae80bb75SLois Curfman McInnes MatAXPY_SeqDense,0,0, 8254b0e389bSBarry Smith MatGetValues_SeqDense, 8264b0e389bSBarry Smith MatCopy_SeqDense}; 8270754003eSLois Curfman McInnes 8284b828684SBarry Smith /*@C 829fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 830d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 831d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 832289bc588SBarry Smith 83320563c6bSBarry Smith Input Parameters: 8340c775827SLois Curfman McInnes . comm - MPI communicator, set to MPI_COMM_SELF 8350c775827SLois Curfman McInnes . m - number of rows 83618f449edSLois Curfman McInnes . n - number of columns 837b4fd4287SBarry Smith . data - optional location of matrix data. Set data=PETSC_NULL for PETSc 838dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 83920563c6bSBarry Smith 84020563c6bSBarry Smith Output Parameter: 8410c775827SLois Curfman McInnes . newmat - the matrix 84220563c6bSBarry Smith 84318f449edSLois Curfman McInnes Notes: 84418f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 84518f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 846b4fd4287SBarry Smith set data=PETSC_NULL. 84718f449edSLois Curfman McInnes 848dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 849d65003e9SLois Curfman McInnes 850dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues() 85120563c6bSBarry Smith @*/ 85218f449edSLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *newmat) 853289bc588SBarry Smith { 854289bc588SBarry Smith Mat mat; 855ec8511deSBarry Smith Mat_SeqDense *l; 85625cdf11fSBarry Smith int ierr,flg; 85755659b69SBarry Smith 858289bc588SBarry Smith *newmat = 0; 85918f449edSLois Curfman McInnes 8600452661fSBarry Smith PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm); 861a5a9c739SBarry Smith PLogObjectCreate(mat); 8623439631bSBarry Smith l = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(l); 863416022c9SBarry Smith PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 864ec8511deSBarry Smith mat->destroy = MatDestroy_SeqDense; 865ec8511deSBarry Smith mat->view = MatView_SeqDense; 866289bc588SBarry Smith mat->factor = 0; 8673439631bSBarry Smith PLogObjectMemory(mat,sizeof(struct _Mat)); 86818f449edSLois Curfman McInnes mat->data = (void *) l; 869d6dfbf8fSBarry Smith 870289bc588SBarry Smith l->m = m; 871289bc588SBarry Smith l->n = n; 872289bc588SBarry Smith l->pivots = 0; 873289bc588SBarry Smith l->roworiented = 1; 874b4fd4287SBarry Smith if (data == PETSC_NULL) { 875ba7af9b1SBarry Smith l->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(l->v); 876cddf8d76SBarry Smith PetscMemzero(l->v,m*n*sizeof(Scalar)); 8772dd2b441SLois Curfman McInnes l->user_alloc = 0; 8783439631bSBarry Smith PLogObjectMemory(mat,n*m*sizeof(Scalar)); 87918f449edSLois Curfman McInnes } 8802dd2b441SLois Curfman McInnes else { /* user-allocated storage */ 8812dd2b441SLois Curfman McInnes l->v = data; 8822dd2b441SLois Curfman McInnes l->user_alloc = 1; 8832dd2b441SLois Curfman McInnes } 88418f449edSLois Curfman McInnes 88525cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 88625cdf11fSBarry Smith if (flg) { 887682d7d0cSBarry Smith ierr = MatPrintHelp(mat); CHKERRQ(ierr); 888682d7d0cSBarry Smith } 889289bc588SBarry Smith *newmat = mat; 890289bc588SBarry Smith return 0; 891289bc588SBarry Smith } 892289bc588SBarry Smith 893c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat) 894289bc588SBarry Smith { 895c0bbcb79SLois Curfman McInnes Mat_SeqDense *m = (Mat_SeqDense *) A->data; 896b4fd4287SBarry Smith return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat); 897289bc588SBarry Smith } 898