1cb512458SBarry Smith #ifndef lint 2*18f449edSLois Curfman McInnes static char vcid[] = "$Id: dense.c,v 1.74 1995/11/03 02:47:26 bsmith Exp curfman $"; 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 { 51289bc588SBarry Smith int ierr; 52c0bbcb79SLois Curfman McInnes ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr); 53289bc588SBarry Smith return 0; 54289bc588SBarry Smith } 55c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 56289bc588SBarry Smith { 5749d8b64dSBarry Smith return MatLUFactor(*fact,0,0,1.0); 58289bc588SBarry Smith } 59c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact) 60289bc588SBarry Smith { 61289bc588SBarry Smith int ierr; 62c0bbcb79SLois Curfman McInnes ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr); 63289bc588SBarry Smith return 0; 64289bc588SBarry Smith } 65c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 66289bc588SBarry Smith { 6749d8b64dSBarry Smith return MatCholeskyFactor(*fact,0,1.0); 68289bc588SBarry Smith } 69c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f) 70289bc588SBarry Smith { 71c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 726abc6512SBarry Smith int info; 7355659b69SBarry Smith 74752f5567SLois Curfman McInnes if (mat->pivots) { 750452661fSBarry Smith PetscFree(mat->pivots); 76c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,-mat->m*sizeof(int)); 77752f5567SLois Curfman McInnes mat->pivots = 0; 78752f5567SLois Curfman McInnes } 79289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 80ec8511deSBarry Smith if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization"); 81c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 8255659b69SBarry Smith PLogFlops((mat->n*mat->n*mat->n)/3); 83289bc588SBarry Smith return 0; 84289bc588SBarry Smith } 85289bc588SBarry Smith 86c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 87289bc588SBarry Smith { 88c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 896abc6512SBarry Smith int one = 1, info; 906abc6512SBarry Smith Scalar *x, *y; 9167e560aaSBarry Smith 92289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 93416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 94c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 9548d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 96289bc588SBarry Smith } 97c0bbcb79SLois Curfman McInnes else if (A->factor == FACTOR_CHOLESKY){ 9848d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 99289bc588SBarry Smith } 100ec8511deSBarry Smith else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve"); 101ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve"); 10255659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 103289bc588SBarry Smith return 0; 104289bc588SBarry Smith } 105c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy) 106da3a660dSBarry Smith { 107c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1086abc6512SBarry Smith int one = 1, info; 1096abc6512SBarry Smith Scalar *x, *y; 11067e560aaSBarry Smith 111da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 112416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 113752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 114da3a660dSBarry Smith if (mat->pivots) { 11548d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 116da3a660dSBarry Smith } 117da3a660dSBarry Smith else { 11848d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 119da3a660dSBarry Smith } 120ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve"); 12155659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 122da3a660dSBarry Smith return 0; 123da3a660dSBarry Smith } 124c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 125da3a660dSBarry Smith { 126c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1276abc6512SBarry Smith int one = 1, info,ierr; 1286abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 129da3a660dSBarry Smith Vec tmp = 0; 13067e560aaSBarry Smith 131da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 132da3a660dSBarry Smith if (yy == zz) { 13378b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 134c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 13578b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 136da3a660dSBarry Smith } 137416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 138752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 139da3a660dSBarry Smith if (mat->pivots) { 14048d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 141da3a660dSBarry Smith } 142da3a660dSBarry Smith else { 14348d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 144da3a660dSBarry Smith } 145ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve"); 146da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 147da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 14855659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 149da3a660dSBarry Smith return 0; 150da3a660dSBarry Smith } 15167e560aaSBarry Smith 152c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 153da3a660dSBarry Smith { 154c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1556abc6512SBarry Smith int one = 1, info,ierr; 1566abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 157da3a660dSBarry Smith Vec tmp; 15867e560aaSBarry Smith 159da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 160da3a660dSBarry Smith if (yy == zz) { 16178b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 162c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 16378b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 164da3a660dSBarry Smith } 165416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 166752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 167da3a660dSBarry Smith if (mat->pivots) { 16848d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 169da3a660dSBarry Smith } 170da3a660dSBarry Smith else { 17148d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 172da3a660dSBarry Smith } 173ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve"); 174da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 175da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 17655659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 177da3a660dSBarry Smith return 0; 178da3a660dSBarry Smith } 179289bc588SBarry Smith /* ------------------------------------------------------------------*/ 180c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 18120e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 182289bc588SBarry Smith { 183c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 184289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 1856abc6512SBarry Smith int o = 1,ierr, m = mat->m, i; 186289bc588SBarry Smith 187289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 188289bc588SBarry Smith /* this is a hack fix, should have another version without 189289bc588SBarry Smith the second BLdot */ 190bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx); CHKERRQ(ierr); 191289bc588SBarry Smith } 192289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 193289bc588SBarry Smith while (its--) { 194289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 195289bc588SBarry Smith for ( i=0; i<m; i++ ) { 196f1747703SBarry Smith #if defined(PETSC_COMPLEX) 197f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 198f1747703SBarry Smith not happy about returning a double complex */ 199f1747703SBarry Smith int _i; 200f1747703SBarry Smith Scalar sum = b[i]; 201f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 202f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 203f1747703SBarry Smith } 204f1747703SBarry Smith xt = sum; 205f1747703SBarry Smith #else 206289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 207f1747703SBarry Smith #endif 208d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 209289bc588SBarry Smith } 210289bc588SBarry Smith } 211289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 212289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 213f1747703SBarry Smith #if defined(PETSC_COMPLEX) 214f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 215f1747703SBarry Smith not happy about returning a double complex */ 216f1747703SBarry Smith int _i; 217f1747703SBarry Smith Scalar sum = b[i]; 218f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 219f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 220f1747703SBarry Smith } 221f1747703SBarry Smith xt = sum; 222f1747703SBarry Smith #else 223289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 224f1747703SBarry Smith #endif 225d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 226289bc588SBarry Smith } 227289bc588SBarry Smith } 228289bc588SBarry Smith } 229289bc588SBarry Smith return 0; 230289bc588SBarry Smith } 231289bc588SBarry Smith 232289bc588SBarry Smith /* -----------------------------------------------------------------*/ 233c0bbcb79SLois Curfman McInnes static int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy) 234289bc588SBarry Smith { 235c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 236289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 237289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 238289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 23948d91487SBarry Smith LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 24055659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->n); 241289bc588SBarry Smith return 0; 242289bc588SBarry Smith } 243c0bbcb79SLois Curfman McInnes static int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 244289bc588SBarry Smith { 245c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 246289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 247289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 248289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 24948d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 25055659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->m); 251289bc588SBarry Smith return 0; 252289bc588SBarry Smith } 253c0bbcb79SLois Curfman McInnes static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 254289bc588SBarry Smith { 255c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 256289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 2576abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 258289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 259416022c9SBarry Smith if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 26048d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 26155659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 262289bc588SBarry Smith return 0; 263289bc588SBarry Smith } 264c0bbcb79SLois Curfman McInnes static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 265289bc588SBarry Smith { 266c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 267289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 268289bc588SBarry Smith int _One=1; 2696abc6512SBarry Smith Scalar _DOne=1.0; 270289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 271289bc588SBarry Smith VecGetArray(zz,&z); 272416022c9SBarry Smith if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 27348d91487SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 27455659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 275289bc588SBarry Smith return 0; 276289bc588SBarry Smith } 277289bc588SBarry Smith 278289bc588SBarry Smith /* -----------------------------------------------------------------*/ 279c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 280289bc588SBarry Smith { 281c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 282289bc588SBarry Smith Scalar *v; 283289bc588SBarry Smith int i; 28467e560aaSBarry Smith 285289bc588SBarry Smith *ncols = mat->n; 286289bc588SBarry Smith if (cols) { 2870452661fSBarry Smith *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols); 28873c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 289289bc588SBarry Smith } 290289bc588SBarry Smith if (vals) { 2910452661fSBarry Smith *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 292289bc588SBarry Smith v = mat->v + row; 29373c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 294289bc588SBarry Smith } 295289bc588SBarry Smith return 0; 296289bc588SBarry Smith } 297c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 298289bc588SBarry Smith { 2990452661fSBarry Smith if (cols) { PetscFree(*cols); } 3000452661fSBarry Smith if (vals) { PetscFree(*vals); } 301289bc588SBarry Smith return 0; 302289bc588SBarry Smith } 303289bc588SBarry Smith /* ----------------------------------------------------------------*/ 304c0bbcb79SLois Curfman McInnes static int MatInsert_SeqDense(Mat A,int m,int *indexm,int n, 305e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 306289bc588SBarry Smith { 307c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 308289bc588SBarry Smith int i,j; 309d6dfbf8fSBarry Smith 310289bc588SBarry Smith if (!mat->roworiented) { 311dbb450caSBarry Smith if (addv == INSERT_VALUES) { 312289bc588SBarry Smith for ( j=0; j<n; j++ ) { 313289bc588SBarry Smith for ( i=0; i<m; i++ ) { 314289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 315289bc588SBarry Smith } 316289bc588SBarry Smith } 317289bc588SBarry Smith } 318289bc588SBarry Smith else { 319289bc588SBarry Smith for ( j=0; j<n; j++ ) { 320289bc588SBarry Smith for ( i=0; i<m; i++ ) { 321289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 322289bc588SBarry Smith } 323289bc588SBarry Smith } 324289bc588SBarry Smith } 325e8d4e0b9SBarry Smith } 326e8d4e0b9SBarry Smith else { 327dbb450caSBarry Smith if (addv == INSERT_VALUES) { 328e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 329e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 330e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 331e8d4e0b9SBarry Smith } 332e8d4e0b9SBarry Smith } 333e8d4e0b9SBarry Smith } 334289bc588SBarry Smith else { 335289bc588SBarry Smith for ( i=0; i<m; i++ ) { 336289bc588SBarry Smith for ( j=0; j<n; j++ ) { 337289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 338289bc588SBarry Smith } 339289bc588SBarry Smith } 340289bc588SBarry Smith } 341e8d4e0b9SBarry Smith } 342289bc588SBarry Smith return 0; 343289bc588SBarry Smith } 344e8d4e0b9SBarry Smith 345289bc588SBarry Smith /* -----------------------------------------------------------------*/ 34655659b69SBarry Smith static int MatCopyPrivate_SeqDense(Mat A,Mat *newmat,int cpvalues) 347289bc588SBarry Smith { 348c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l; 349289bc588SBarry Smith int ierr; 350289bc588SBarry Smith Mat newi; 35148d91487SBarry Smith 352*18f449edSLois Curfman McInnes ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,0,&newi); CHKERRQ(ierr); 353ec8511deSBarry Smith l = (Mat_SeqDense *) newi->data; 35455659b69SBarry Smith if (cpvalues == COPY_VALUES) { 355416022c9SBarry Smith PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 35655659b69SBarry Smith } 357289bc588SBarry Smith *newmat = newi; 358289bc588SBarry Smith return 0; 359289bc588SBarry Smith } 360289bc588SBarry Smith 36188e32edaSLois Curfman McInnes #include "sysio.h" 36288e32edaSLois Curfman McInnes 36388e32edaSLois Curfman McInnes int MatLoad_SeqDense(Viewer bview,MatType type,Mat *A) 36488e32edaSLois Curfman McInnes { 36588e32edaSLois Curfman McInnes Mat_SeqDense *a; 36688e32edaSLois Curfman McInnes Mat B; 36788e32edaSLois Curfman McInnes int *scols, i, j, nz, ierr, fd, header[4], size; 36888e32edaSLois Curfman McInnes int *rowlengths = 0, M, N, *cols; 36988e32edaSLois Curfman McInnes Scalar *vals, *svals, *v; 37088e32edaSLois Curfman McInnes PetscObject vobj = (PetscObject) bview; 37188e32edaSLois Curfman McInnes MPI_Comm comm = vobj->comm; 37288e32edaSLois Curfman McInnes 37388e32edaSLois Curfman McInnes MPI_Comm_size(comm,&size); 37488e32edaSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor"); 37588e32edaSLois Curfman McInnes ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 37688e32edaSLois Curfman McInnes ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 37788e32edaSLois Curfman McInnes if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object"); 37888e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 37988e32edaSLois Curfman McInnes 38088e32edaSLois Curfman McInnes /* read row lengths */ 3810452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 38288e32edaSLois Curfman McInnes ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 38388e32edaSLois Curfman McInnes 38488e32edaSLois Curfman McInnes /* create our matrix */ 385*18f449edSLois Curfman McInnes ierr = MatCreateSeqDense(comm,M,N,0,A); CHKERRQ(ierr); 38688e32edaSLois Curfman McInnes B = *A; 38788e32edaSLois Curfman McInnes a = (Mat_SeqDense *) B->data; 38888e32edaSLois Curfman McInnes v = a->v; 38988e32edaSLois Curfman McInnes 39088e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 3910452661fSBarry Smith cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols); 39288e32edaSLois Curfman McInnes ierr = SYRead(fd,cols,nz,SYINT); CHKERRQ(ierr); 3930452661fSBarry Smith vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 39488e32edaSLois Curfman McInnes ierr = SYRead(fd,vals,nz,SYSCALAR); CHKERRQ(ierr); 39588e32edaSLois Curfman McInnes 39688e32edaSLois Curfman McInnes /* insert into matrix */ 39788e32edaSLois Curfman McInnes for ( i=0; i<M; i++ ) { 39888e32edaSLois Curfman McInnes for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j]; 39988e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 40088e32edaSLois Curfman McInnes } 4010452661fSBarry Smith PetscFree(vals); PetscFree(cols); PetscFree(rowlengths); 40288e32edaSLois Curfman McInnes 40388e32edaSLois Curfman McInnes ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 40488e32edaSLois Curfman McInnes ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 40588e32edaSLois Curfman McInnes return 0; 40688e32edaSLois Curfman McInnes } 40788e32edaSLois Curfman McInnes 408932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h" 409932b0c3eSLois Curfman McInnes #include "sysio.h" 410932b0c3eSLois Curfman McInnes 411932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 412289bc588SBarry Smith { 413932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 414932b0c3eSLois Curfman McInnes int ierr, i, j, format; 415d636dbe3SBarry Smith FILE *fd; 416932b0c3eSLois Curfman McInnes char *outputname; 417932b0c3eSLois Curfman McInnes Scalar *v; 418932b0c3eSLois Curfman McInnes 419932b0c3eSLois Curfman McInnes ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 420932b0c3eSLois Curfman McInnes ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 421932b0c3eSLois Curfman McInnes ierr = ViewerFileGetFormat_Private(viewer,&format); 422f72585f2SLois Curfman McInnes if (format == FILE_FORMAT_INFO) { 423932b0c3eSLois Curfman McInnes ; /* do nothing for now */ 424f72585f2SLois Curfman McInnes } 425f72585f2SLois Curfman McInnes else { 426932b0c3eSLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 427932b0c3eSLois Curfman McInnes v = a->v + i; 428932b0c3eSLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 429289bc588SBarry Smith #if defined(PETSC_COMPLEX) 430932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m; 431289bc588SBarry Smith #else 432932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e ",*v); v += a->m; 433289bc588SBarry Smith #endif 434289bc588SBarry Smith } 4358e37dc09SBarry Smith fprintf(fd,"\n"); 436289bc588SBarry Smith } 437da3a660dSBarry Smith } 438932b0c3eSLois Curfman McInnes fflush(fd); 439289bc588SBarry Smith return 0; 440289bc588SBarry Smith } 441289bc588SBarry Smith 442932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 443932b0c3eSLois Curfman McInnes { 444932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 445932b0c3eSLois Curfman McInnes int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 446932b0c3eSLois Curfman McInnes Scalar *v, *anonz; 447932b0c3eSLois Curfman McInnes 448932b0c3eSLois Curfman McInnes ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 4490452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens); 450932b0c3eSLois Curfman McInnes col_lens[0] = MAT_COOKIE; 451932b0c3eSLois Curfman McInnes col_lens[1] = m; 452932b0c3eSLois Curfman McInnes col_lens[2] = n; 453932b0c3eSLois Curfman McInnes col_lens[3] = nz; 454932b0c3eSLois Curfman McInnes 455932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 456932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) col_lens[4+i] = n; 457932b0c3eSLois Curfman McInnes ierr = SYWrite(fd,col_lens,4+m,SYINT,1); CHKERRQ(ierr); 458932b0c3eSLois Curfman McInnes 459932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 460932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 461932b0c3eSLois Curfman McInnes ict = 0; 462932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 463932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) col_lens[ict++] = j; 464932b0c3eSLois Curfman McInnes } 465932b0c3eSLois Curfman McInnes ierr = SYWrite(fd,col_lens,nz,SYINT,0); CHKERRQ(ierr); 4660452661fSBarry Smith PetscFree(col_lens); 467932b0c3eSLois Curfman McInnes 468932b0c3eSLois Curfman McInnes /* store nonzero values */ 4690452661fSBarry Smith anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz); 470932b0c3eSLois Curfman McInnes ict = 0; 471932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 472932b0c3eSLois Curfman McInnes v = a->v + i; 473932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) { 474932b0c3eSLois Curfman McInnes anonz[ict++] = *v; v += a->m; 475932b0c3eSLois Curfman McInnes } 476932b0c3eSLois Curfman McInnes } 477932b0c3eSLois Curfman McInnes ierr = SYWrite(fd,anonz,nz,SYSCALAR,0); CHKERRQ(ierr); 4780452661fSBarry Smith PetscFree(anonz); 479932b0c3eSLois Curfman McInnes return 0; 480932b0c3eSLois Curfman McInnes } 481932b0c3eSLois Curfman McInnes 482932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer) 483932b0c3eSLois Curfman McInnes { 484932b0c3eSLois Curfman McInnes Mat A = (Mat) obj; 485932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*) A->data; 486932b0c3eSLois Curfman McInnes PetscObject vobj = (PetscObject) viewer; 487932b0c3eSLois Curfman McInnes 488932b0c3eSLois Curfman McInnes if (!viewer) { 489932b0c3eSLois Curfman McInnes viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 490932b0c3eSLois Curfman McInnes } 491932b0c3eSLois Curfman McInnes if (vobj->cookie == VIEWER_COOKIE) { 492932b0c3eSLois Curfman McInnes if (vobj->type == MATLAB_VIEWER) { 493932b0c3eSLois Curfman McInnes return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); 494932b0c3eSLois Curfman McInnes } 495932b0c3eSLois Curfman McInnes else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) { 496932b0c3eSLois Curfman McInnes return MatView_SeqDense_ASCII(A,viewer); 497932b0c3eSLois Curfman McInnes } 498932b0c3eSLois Curfman McInnes else if (vobj->type == BINARY_FILE_VIEWER) { 499932b0c3eSLois Curfman McInnes return MatView_SeqDense_Binary(A,viewer); 500932b0c3eSLois Curfman McInnes } 501932b0c3eSLois Curfman McInnes } 502932b0c3eSLois Curfman McInnes return 0; 503932b0c3eSLois Curfman McInnes } 504289bc588SBarry Smith 505ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj) 506289bc588SBarry Smith { 507289bc588SBarry Smith Mat mat = (Mat) obj; 508ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 509a5a9c739SBarry Smith #if defined(PETSC_LOG) 510a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 511a5a9c739SBarry Smith #endif 5120452661fSBarry Smith if (l->pivots) PetscFree(l->pivots); 5130452661fSBarry Smith PetscFree(l); 514a5a9c739SBarry Smith PLogObjectDestroy(mat); 5150452661fSBarry Smith PetscHeaderDestroy(mat); 516289bc588SBarry Smith return 0; 517289bc588SBarry Smith } 518289bc588SBarry Smith 519c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout) 520289bc588SBarry Smith { 521c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 522d3e5ee88SLois Curfman McInnes int k, j, m, n; 523d3e5ee88SLois Curfman McInnes Scalar *v, tmp; 52448b35521SBarry Smith 525d3e5ee88SLois Curfman McInnes v = mat->v; m = mat->m; n = mat->n; 52648b35521SBarry Smith if (!matout) { /* in place transpose */ 527d9f96c7cSLois Curfman McInnes if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place"); 528d3e5ee88SLois Curfman McInnes for ( j=0; j<m; j++ ) { 529289bc588SBarry Smith for ( k=0; k<j; k++ ) { 530d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 531d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 532d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 533289bc588SBarry Smith } 534289bc588SBarry Smith } 53548b35521SBarry Smith } 536d3e5ee88SLois Curfman McInnes else { /* out-of-place transpose */ 537d3e5ee88SLois Curfman McInnes int ierr; 538d3e5ee88SLois Curfman McInnes Mat tmat; 539ec8511deSBarry Smith Mat_SeqDense *tmatd; 540d3e5ee88SLois Curfman McInnes Scalar *v2; 541*18f449edSLois Curfman McInnes ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,0,&tmat); CHKERRQ(ierr); 542ec8511deSBarry Smith tmatd = (Mat_SeqDense *) tmat->data; 5430de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 544d3e5ee88SLois Curfman McInnes for ( j=0; j<n; j++ ) { 5450de55854SLois Curfman McInnes for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 546d3e5ee88SLois Curfman McInnes } 547d3e5ee88SLois Curfman McInnes ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 548d3e5ee88SLois Curfman McInnes ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 549d3e5ee88SLois Curfman McInnes *matout = tmat; 55048b35521SBarry Smith } 551289bc588SBarry Smith return 0; 552289bc588SBarry Smith } 553289bc588SBarry Smith 554c0bbcb79SLois Curfman McInnes static int MatEqual_SeqDense(Mat A1,Mat A2) 555289bc588SBarry Smith { 556c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 557c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 558289bc588SBarry Smith int i; 559289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 560289bc588SBarry Smith if (mat1->m != mat2->m) return 0; 561289bc588SBarry Smith if (mat1->n != mat2->n) return 0; 562289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 563289bc588SBarry Smith if (*v1 != *v2) return 0; 564289bc588SBarry Smith v1++; v2++; 565289bc588SBarry Smith } 566289bc588SBarry Smith return 1; 567289bc588SBarry Smith } 568289bc588SBarry Smith 569c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v) 570289bc588SBarry Smith { 571c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 5726abc6512SBarry Smith int i, n; 5736abc6512SBarry Smith Scalar *x; 574289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 575ec8511deSBarry Smith if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 576289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 577289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 578289bc588SBarry Smith } 579289bc588SBarry Smith return 0; 580289bc588SBarry Smith } 581289bc588SBarry Smith 582c0bbcb79SLois Curfman McInnes static int MatScale_SeqDense(Mat A,Vec ll,Vec rr) 583289bc588SBarry Smith { 584c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 585da3a660dSBarry Smith Scalar *l,*r,x,*v; 586da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 58755659b69SBarry Smith 58828988994SBarry Smith if (ll) { 589da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 590ec8511deSBarry Smith if (m != mat->m) SETERRQ(1,"MatScale_SeqDense:Left scaling vec wrong size"); 59155659b69SBarry Smith PLogFlops(n*m); 592da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 593da3a660dSBarry Smith x = l[i]; 594da3a660dSBarry Smith v = mat->v + i; 595da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 596da3a660dSBarry Smith } 597da3a660dSBarry Smith } 59828988994SBarry Smith if (rr) { 599da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 600ec8511deSBarry Smith if (n != mat->n) SETERRQ(1,"MatScale_SeqDense:Right scaling vec wrong size"); 60155659b69SBarry Smith PLogFlops(n*m); 602da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 603da3a660dSBarry Smith x = r[i]; 604da3a660dSBarry Smith v = mat->v + i*m; 605da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 606da3a660dSBarry Smith } 607da3a660dSBarry Smith } 608289bc588SBarry Smith return 0; 609289bc588SBarry Smith } 610289bc588SBarry Smith 611cddf8d76SBarry Smith static int MatNorm_SeqDense(Mat A,NormType type,double *norm) 612289bc588SBarry Smith { 613c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 614289bc588SBarry Smith Scalar *v = mat->v; 615289bc588SBarry Smith double sum = 0.0; 616289bc588SBarry Smith int i, j; 61755659b69SBarry Smith 618289bc588SBarry Smith if (type == NORM_FROBENIUS) { 619289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 620289bc588SBarry Smith #if defined(PETSC_COMPLEX) 621289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 622289bc588SBarry Smith #else 623289bc588SBarry Smith sum += (*v)*(*v); v++; 624289bc588SBarry Smith #endif 625289bc588SBarry Smith } 626289bc588SBarry Smith *norm = sqrt(sum); 62755659b69SBarry Smith PLogFlops(2*mat->n*mat->m); 628289bc588SBarry Smith } 629289bc588SBarry Smith else if (type == NORM_1) { 630289bc588SBarry Smith *norm = 0.0; 631289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 632289bc588SBarry Smith sum = 0.0; 633289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 63433a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 635289bc588SBarry Smith } 636289bc588SBarry Smith if (sum > *norm) *norm = sum; 637289bc588SBarry Smith } 63855659b69SBarry Smith PLogFlops(mat->n*mat->m); 639289bc588SBarry Smith } 640289bc588SBarry Smith else if (type == NORM_INFINITY) { 641289bc588SBarry Smith *norm = 0.0; 642289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 643289bc588SBarry Smith v = mat->v + j; 644289bc588SBarry Smith sum = 0.0; 645289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 646cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v += mat->m; 647289bc588SBarry Smith } 648289bc588SBarry Smith if (sum > *norm) *norm = sum; 649289bc588SBarry Smith } 65055659b69SBarry Smith PLogFlops(mat->n*mat->m); 651289bc588SBarry Smith } 652289bc588SBarry Smith else { 65348d91487SBarry Smith SETERRQ(1,"MatNorm_SeqDense:No two norm"); 654289bc588SBarry Smith } 655289bc588SBarry Smith return 0; 656289bc588SBarry Smith } 657289bc588SBarry Smith 658c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op) 659289bc588SBarry Smith { 660c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 66167e560aaSBarry Smith 662289bc588SBarry Smith if (op == ROW_ORIENTED) aij->roworiented = 1; 663289bc588SBarry Smith else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 664c0bbcb79SLois Curfman McInnes else if (op == ROWS_SORTED || 66588e32edaSLois Curfman McInnes op == COLUMNS_SORTED || 666c0bbcb79SLois Curfman McInnes op == SYMMETRIC_MATRIX || 667c0bbcb79SLois Curfman McInnes op == STRUCTURALLY_SYMMETRIC_MATRIX || 668c0bbcb79SLois Curfman McInnes op == NO_NEW_NONZERO_LOCATIONS || 669c0bbcb79SLois Curfman McInnes op == YES_NEW_NONZERO_LOCATIONS || 670c0bbcb79SLois Curfman McInnes op == NO_NEW_DIAGONALS || 671c0bbcb79SLois Curfman McInnes op == YES_NEW_DIAGONALS) 672c0bbcb79SLois Curfman McInnes PLogInfo((PetscObject)A,"Info:MatSetOption_SeqDense:Option ignored\n"); 673c0bbcb79SLois Curfman McInnes else 6744d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");} 675289bc588SBarry Smith return 0; 676289bc588SBarry Smith } 677289bc588SBarry Smith 678ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A) 6796f0a148fSBarry Smith { 680ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 681cddf8d76SBarry Smith PetscMemzero(l->v,l->m*l->n*sizeof(Scalar)); 6826f0a148fSBarry Smith return 0; 6836f0a148fSBarry Smith } 6846f0a148fSBarry Smith 685ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 6866f0a148fSBarry Smith { 687ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 6886abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 6896f0a148fSBarry Smith Scalar *slot; 69055659b69SBarry Smith 69178b31e54SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 69278b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 6936f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 6946f0a148fSBarry Smith slot = l->v + rows[i]; 6956f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 6966f0a148fSBarry Smith } 6976f0a148fSBarry Smith if (diag) { 6986f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 6996f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 700260925b8SBarry Smith *slot = *diag; 7016f0a148fSBarry Smith } 7026f0a148fSBarry Smith } 703260925b8SBarry Smith ISRestoreIndices(is,&rows); 7046f0a148fSBarry Smith return 0; 7056f0a148fSBarry Smith } 706557bce09SLois Curfman McInnes 707c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n) 708557bce09SLois Curfman McInnes { 709c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 710557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 711557bce09SLois Curfman McInnes return 0; 712557bce09SLois Curfman McInnes } 713557bce09SLois Curfman McInnes 714c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 715d3e5ee88SLois Curfman McInnes { 716c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 717d3e5ee88SLois Curfman McInnes *m = 0; *n = mat->m; 718d3e5ee88SLois Curfman McInnes return 0; 719d3e5ee88SLois Curfman McInnes } 720d3e5ee88SLois Curfman McInnes 721c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array) 72264e87e97SBarry Smith { 723c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 72464e87e97SBarry Smith *array = mat->v; 72564e87e97SBarry Smith return 0; 72664e87e97SBarry Smith } 7270754003eSLois Curfman McInnes 7280754003eSLois Curfman McInnes 729c0bbcb79SLois Curfman McInnes static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol) 7300754003eSLois Curfman McInnes { 731ec8511deSBarry Smith SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done"); 7320754003eSLois Curfman McInnes } 7330754003eSLois Curfman McInnes 734cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall, 735cddf8d76SBarry Smith Mat *submat) 7360754003eSLois Curfman McInnes { 737c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 7380754003eSLois Curfman McInnes int nznew, *smap, i, j, ierr, oldcols = mat->n; 739160018dcSBarry Smith int *irow, *icol, nrows, ncols, *cwork; 7400754003eSLois Curfman McInnes Scalar *vwork, *val; 7410754003eSLois Curfman McInnes Mat newmat; 7420754003eSLois Curfman McInnes 74378b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 74478b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 74578b31e54SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 74678b31e54SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 7470754003eSLois Curfman McInnes 7480452661fSBarry Smith smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap); 7490452661fSBarry Smith cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork); 7500452661fSBarry Smith vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 751cddf8d76SBarry Smith PetscMemzero((char*)smap,oldcols*sizeof(int)); 7520754003eSLois Curfman McInnes for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 7530754003eSLois Curfman McInnes 7540754003eSLois Curfman McInnes /* Create and fill new matrix */ 755*18f449edSLois Curfman McInnes ierr = MatCreateSeqDense(A->comm,nrows,ncols,0,&newmat);CHKERRQ(ierr); 7560754003eSLois Curfman McInnes for (i=0; i<nrows; i++) { 7570754003eSLois Curfman McInnes nznew = 0; 7580754003eSLois Curfman McInnes val = mat->v + irow[i]; 7590754003eSLois Curfman McInnes for (j=0; j<oldcols; j++) { 7600754003eSLois Curfman McInnes if (smap[j]) { 7610754003eSLois Curfman McInnes cwork[nznew] = smap[j] - 1; 7620754003eSLois Curfman McInnes vwork[nznew++] = val[j * mat->m]; 7630754003eSLois Curfman McInnes } 7640754003eSLois Curfman McInnes } 765dbb450caSBarry Smith ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES); 76678b31e54SBarry Smith CHKERRQ(ierr); 7670754003eSLois Curfman McInnes } 76878b31e54SBarry Smith ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 76978b31e54SBarry Smith ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 7700754003eSLois Curfman McInnes 7710754003eSLois Curfman McInnes /* Free work space */ 7720452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 77378b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 77478b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 7750754003eSLois Curfman McInnes *submat = newmat; 7760754003eSLois Curfman McInnes return 0; 7770754003eSLois Curfman McInnes } 7780754003eSLois Curfman McInnes 779289bc588SBarry Smith /* -------------------------------------------------------------------*/ 780ec8511deSBarry Smith static struct _MatOps MatOps = {MatInsert_SeqDense, 781ec8511deSBarry Smith MatGetRow_SeqDense, MatRestoreRow_SeqDense, 782ec8511deSBarry Smith MatMult_SeqDense, MatMultAdd_SeqDense, 783ec8511deSBarry Smith MatMultTrans_SeqDense, MatMultTransAdd_SeqDense, 784ec8511deSBarry Smith MatSolve_SeqDense,MatSolveAdd_SeqDense, 785ec8511deSBarry Smith MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense, 786ec8511deSBarry Smith MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense, 787ec8511deSBarry Smith MatRelax_SeqDense, 788ec8511deSBarry Smith MatTranspose_SeqDense, 789ec8511deSBarry Smith MatGetInfo_SeqDense,MatEqual_SeqDense, 790ec8511deSBarry Smith MatGetDiagonal_SeqDense,MatScale_SeqDense,MatNorm_SeqDense, 791289bc588SBarry Smith 0,0, 792ec8511deSBarry Smith 0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0, 793ec8511deSBarry Smith MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense, 794ec8511deSBarry Smith MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense, 795ec8511deSBarry Smith MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense, 796ec8511deSBarry Smith 0,0,MatGetArray_SeqDense,0,0, 797ec8511deSBarry Smith MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense, 7981987afe7SBarry Smith MatCopyPrivate_SeqDense,0,0,0,0, 7991987afe7SBarry Smith MatAXPY_SeqDense}; 8000754003eSLois Curfman McInnes 8014b828684SBarry Smith /*@C 802fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 803d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 804d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 805289bc588SBarry Smith 80620563c6bSBarry Smith Input Parameters: 8070c775827SLois Curfman McInnes . comm - MPI communicator, set to MPI_COMM_SELF 8080c775827SLois Curfman McInnes . m - number of rows 809*18f449edSLois Curfman McInnes . n - number of columns 810*18f449edSLois Curfman McInnes . data - optional location of matrix data. Set data=0 for PETSc to 811*18f449edSLois Curfman McInnes control all matrix memory allocation. 81220563c6bSBarry Smith 81320563c6bSBarry Smith Output Parameter: 8140c775827SLois Curfman McInnes . newmat - the matrix 81520563c6bSBarry Smith 816*18f449edSLois Curfman McInnes Notes: 817*18f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 818*18f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 819*18f449edSLois Curfman McInnes set data=0. 820*18f449edSLois Curfman McInnes 821dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 822d65003e9SLois Curfman McInnes 823dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues() 82420563c6bSBarry Smith @*/ 825*18f449edSLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *newmat) 826289bc588SBarry Smith { 827ec8511deSBarry Smith int size = sizeof(Mat_SeqDense) + m*n*sizeof(Scalar); 828289bc588SBarry Smith Mat mat; 829ec8511deSBarry Smith Mat_SeqDense *l; 83055659b69SBarry Smith 831289bc588SBarry Smith *newmat = 0; 832*18f449edSLois Curfman McInnes 833*18f449edSLois Curfman McInnes if (!data) size = sizeof(Mat_SeqDense) + m*n*sizeof(Scalar); 834*18f449edSLois Curfman McInnes else size = sizeof(Mat_SeqDense); 8350452661fSBarry Smith PetscHeaderCreate(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm); 836a5a9c739SBarry Smith PLogObjectCreate(mat); 8370452661fSBarry Smith l = (Mat_SeqDense *) PetscMalloc(size); CHKPTRQ(l); 838416022c9SBarry Smith PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 839ec8511deSBarry Smith mat->destroy = MatDestroy_SeqDense; 840ec8511deSBarry Smith mat->view = MatView_SeqDense; 841289bc588SBarry Smith mat->factor = 0; 842752f5567SLois Curfman McInnes PLogObjectMemory(mat,sizeof(struct _Mat) + size); 843*18f449edSLois Curfman McInnes mat->data = (void *) l; 844d6dfbf8fSBarry Smith 845289bc588SBarry Smith l->m = m; 846289bc588SBarry Smith l->n = n; 847289bc588SBarry Smith l->pivots = 0; 848289bc588SBarry Smith l->roworiented = 1; 849*18f449edSLois Curfman McInnes if (!data) { 850*18f449edSLois Curfman McInnes l->v = (Scalar *) (l + 1); 851cddf8d76SBarry Smith PetscMemzero(l->v,m*n*sizeof(Scalar)); 852*18f449edSLois Curfman McInnes } 853*18f449edSLois Curfman McInnes else l->v = data; /* user-allocated storage */ 854*18f449edSLois Curfman McInnes 855289bc588SBarry Smith *newmat = mat; 856289bc588SBarry Smith return 0; 857289bc588SBarry Smith } 858289bc588SBarry Smith 859c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat) 860289bc588SBarry Smith { 861c0bbcb79SLois Curfman McInnes Mat_SeqDense *m = (Mat_SeqDense *) A->data; 862*18f449edSLois Curfman McInnes return MatCreateSeqDense(A->comm,m->m,m->n,0,newmat); 863289bc588SBarry Smith } 864