1cb512458SBarry Smith #ifndef lint 2*3d60af36SSatish Balay static char vcid[] = "$Id: dense.c,v 1.103 1996/04/09 02:23:22 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++;} 27bcd2baecSBarry Smith if (nz) *nz = count; 28bcd2baecSBarry Smith if (nzalloc) *nzalloc = N; 29bcd2baecSBarry Smith if (mem) *mem = (int)A->mem; 30fa9ec3f1SLois Curfman McInnes return 0; 31289bc588SBarry Smith } 32289bc588SBarry Smith 3380cd9d93SLois Curfman McInnes static int MatScale_SeqDense(Scalar *alpha,Mat inA) 3480cd9d93SLois Curfman McInnes { 3580cd9d93SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) inA->data; 3680cd9d93SLois Curfman McInnes int one = 1, nz; 3780cd9d93SLois Curfman McInnes 3880cd9d93SLois Curfman McInnes nz = a->m*a->n; 3980cd9d93SLois Curfman McInnes BLscal_( &nz, alpha, a->v, &one ); 4080cd9d93SLois Curfman McInnes PLogFlops(nz); 4180cd9d93SLois Curfman McInnes return 0; 4280cd9d93SLois Curfman McInnes } 4380cd9d93SLois Curfman McInnes 44289bc588SBarry Smith /* ---------------------------------------------------------------*/ 45289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots, 46289bc588SBarry Smith rather than put it in the Mat->col slot.*/ 47c0bbcb79SLois Curfman McInnes static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f) 48289bc588SBarry Smith { 49c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 506abc6512SBarry Smith int info; 5155659b69SBarry Smith 52289bc588SBarry Smith if (!mat->pivots) { 530452661fSBarry Smith mat->pivots = (int *) PetscMalloc(mat->m*sizeof(int));CHKPTRQ(mat->pivots); 54c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,mat->m*sizeof(int)); 55289bc588SBarry Smith } 56289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 57*3d60af36SSatish Balay if (info<0) SETERRQ(1,"MatLUFactor_SeqDense:Bad argument to LU factorization"); 58*3d60af36SSatish Balay if (info>0) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization"); 59c0bbcb79SLois Curfman McInnes A->factor = FACTOR_LU; 6055659b69SBarry Smith PLogFlops((2*mat->n*mat->n*mat->n)/3); 61289bc588SBarry Smith return 0; 62289bc588SBarry Smith } 63c0bbcb79SLois Curfman McInnes static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact) 64289bc588SBarry Smith { 65a57ad546SLois Curfman McInnes return MatConvert(A,MATSAME,fact); 66289bc588SBarry Smith } 67c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 68289bc588SBarry Smith { 6949d8b64dSBarry Smith return MatLUFactor(*fact,0,0,1.0); 70289bc588SBarry Smith } 71c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact) 72289bc588SBarry Smith { 73a57ad546SLois Curfman McInnes return MatConvert(A,MATSAME,fact); 74289bc588SBarry Smith } 75c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 76289bc588SBarry Smith { 7749d8b64dSBarry Smith return MatCholeskyFactor(*fact,0,1.0); 78289bc588SBarry Smith } 79c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f) 80289bc588SBarry Smith { 81c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 826abc6512SBarry Smith int info; 8355659b69SBarry Smith 84752f5567SLois Curfman McInnes if (mat->pivots) { 850452661fSBarry Smith PetscFree(mat->pivots); 86c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,-mat->m*sizeof(int)); 87752f5567SLois Curfman McInnes mat->pivots = 0; 88752f5567SLois Curfman McInnes } 89289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 90ec8511deSBarry Smith if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization"); 91c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 9255659b69SBarry Smith PLogFlops((mat->n*mat->n*mat->n)/3); 93289bc588SBarry Smith return 0; 94289bc588SBarry Smith } 95289bc588SBarry Smith 96c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 97289bc588SBarry Smith { 98c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 99a57ad546SLois Curfman McInnes int one = 1, info, ierr; 1006abc6512SBarry Smith Scalar *x, *y; 10167e560aaSBarry Smith 102a57ad546SLois Curfman McInnes ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 103416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 104c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 10548d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 106289bc588SBarry Smith } 107c0bbcb79SLois Curfman McInnes else if (A->factor == FACTOR_CHOLESKY){ 10848d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 109289bc588SBarry Smith } 110ec8511deSBarry Smith else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve"); 111ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve"); 11255659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 113289bc588SBarry Smith return 0; 114289bc588SBarry Smith } 115c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy) 116da3a660dSBarry Smith { 117c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1186abc6512SBarry Smith int one = 1, info; 1196abc6512SBarry Smith Scalar *x, *y; 12067e560aaSBarry Smith 121da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 122416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 123752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 124da3a660dSBarry Smith if (mat->pivots) { 12548d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 126da3a660dSBarry Smith } 127da3a660dSBarry Smith else { 12848d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 129da3a660dSBarry Smith } 130ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve"); 13155659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 132da3a660dSBarry Smith return 0; 133da3a660dSBarry Smith } 134c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 135da3a660dSBarry Smith { 136c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1376abc6512SBarry Smith int one = 1, info,ierr; 1386abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 139da3a660dSBarry Smith Vec tmp = 0; 14067e560aaSBarry Smith 141da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 142da3a660dSBarry Smith if (yy == zz) { 14378b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 144c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 14578b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 146da3a660dSBarry Smith } 147416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 148752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 149da3a660dSBarry Smith if (mat->pivots) { 15048d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 151da3a660dSBarry Smith } 152da3a660dSBarry Smith else { 15348d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 154da3a660dSBarry Smith } 155ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve"); 156da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 157da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 15855659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 159da3a660dSBarry Smith return 0; 160da3a660dSBarry Smith } 16167e560aaSBarry Smith 162c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 163da3a660dSBarry Smith { 164c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1656abc6512SBarry Smith int one = 1, info,ierr; 1666abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 167da3a660dSBarry Smith Vec tmp; 16867e560aaSBarry Smith 169da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 170da3a660dSBarry Smith if (yy == zz) { 17178b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 172c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 17378b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 174da3a660dSBarry Smith } 175416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 176752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 177da3a660dSBarry Smith if (mat->pivots) { 17848d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 179da3a660dSBarry Smith } 180da3a660dSBarry Smith else { 18148d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 182da3a660dSBarry Smith } 183ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve"); 184da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 185da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 18655659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 187da3a660dSBarry Smith return 0; 188da3a660dSBarry Smith } 189289bc588SBarry Smith /* ------------------------------------------------------------------*/ 190c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 19120e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 192289bc588SBarry Smith { 193c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 194289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 1956abc6512SBarry Smith int o = 1,ierr, m = mat->m, i; 196289bc588SBarry Smith 197289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 198289bc588SBarry Smith /* this is a hack fix, should have another version without 199289bc588SBarry Smith the second BLdot */ 200bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx); CHKERRQ(ierr); 201289bc588SBarry Smith } 202289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 203289bc588SBarry Smith while (its--) { 204289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 205289bc588SBarry Smith for ( i=0; i<m; i++ ) { 206f1747703SBarry Smith #if defined(PETSC_COMPLEX) 207f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 208f1747703SBarry Smith not happy about returning a double complex */ 209f1747703SBarry Smith int _i; 210f1747703SBarry Smith Scalar sum = b[i]; 211f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 212f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 213f1747703SBarry Smith } 214f1747703SBarry Smith xt = sum; 215f1747703SBarry Smith #else 216289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 217f1747703SBarry Smith #endif 218d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 219289bc588SBarry Smith } 220289bc588SBarry Smith } 221289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 222289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 223f1747703SBarry Smith #if defined(PETSC_COMPLEX) 224f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 225f1747703SBarry Smith not happy about returning a double complex */ 226f1747703SBarry Smith int _i; 227f1747703SBarry Smith Scalar sum = b[i]; 228f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 229f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 230f1747703SBarry Smith } 231f1747703SBarry Smith xt = sum; 232f1747703SBarry Smith #else 233289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 234f1747703SBarry Smith #endif 235d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 236289bc588SBarry Smith } 237289bc588SBarry Smith } 238289bc588SBarry Smith } 239289bc588SBarry Smith return 0; 240289bc588SBarry Smith } 241289bc588SBarry Smith 242289bc588SBarry Smith /* -----------------------------------------------------------------*/ 24344cd7ae7SLois Curfman McInnes int MatMultTrans_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_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 25055659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->n); 251289bc588SBarry Smith return 0; 252289bc588SBarry Smith } 25344cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 254289bc588SBarry Smith { 255c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 256289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 257289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 258289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 25948d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 26055659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->m); 261289bc588SBarry Smith return 0; 262289bc588SBarry Smith } 26344cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 264289bc588SBarry Smith { 265c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 266289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 2676abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 268289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 269416022c9SBarry Smith if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 27048d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 27155659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 272289bc588SBarry Smith return 0; 273289bc588SBarry Smith } 27444cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 275289bc588SBarry Smith { 276c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 277289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 278289bc588SBarry Smith int _One=1; 2796abc6512SBarry Smith Scalar _DOne=1.0; 28044cd7ae7SLois Curfman McInnes VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 281717eeb74SLois Curfman McInnes if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar)); 28248d91487SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 28355659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 284289bc588SBarry Smith return 0; 285289bc588SBarry Smith } 286289bc588SBarry Smith 287289bc588SBarry Smith /* -----------------------------------------------------------------*/ 288c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 289289bc588SBarry Smith { 290c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 291289bc588SBarry Smith Scalar *v; 292289bc588SBarry Smith int i; 29367e560aaSBarry Smith 294289bc588SBarry Smith *ncols = mat->n; 295289bc588SBarry Smith if (cols) { 2960452661fSBarry Smith *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols); 29773c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 298289bc588SBarry Smith } 299289bc588SBarry Smith if (vals) { 3000452661fSBarry Smith *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 301289bc588SBarry Smith v = mat->v + row; 30273c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 303289bc588SBarry Smith } 304289bc588SBarry Smith return 0; 305289bc588SBarry Smith } 306c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 307289bc588SBarry Smith { 3080452661fSBarry Smith if (cols) { PetscFree(*cols); } 3090452661fSBarry Smith if (vals) { PetscFree(*vals); } 310289bc588SBarry Smith return 0; 311289bc588SBarry Smith } 312289bc588SBarry Smith /* ----------------------------------------------------------------*/ 313ae80bb75SLois Curfman McInnes static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 314e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 315289bc588SBarry Smith { 316c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 317289bc588SBarry Smith int i,j; 318d6dfbf8fSBarry Smith 319289bc588SBarry Smith if (!mat->roworiented) { 320dbb450caSBarry Smith if (addv == INSERT_VALUES) { 321289bc588SBarry Smith for ( j=0; j<n; j++ ) { 322289bc588SBarry Smith for ( i=0; i<m; i++ ) { 323289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 324289bc588SBarry Smith } 325289bc588SBarry Smith } 326289bc588SBarry Smith } 327289bc588SBarry Smith else { 328289bc588SBarry Smith for ( j=0; j<n; j++ ) { 329289bc588SBarry Smith for ( i=0; i<m; i++ ) { 330289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 331289bc588SBarry Smith } 332289bc588SBarry Smith } 333289bc588SBarry Smith } 334e8d4e0b9SBarry Smith } 335e8d4e0b9SBarry Smith else { 336dbb450caSBarry Smith if (addv == INSERT_VALUES) { 337e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 338e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 339e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 340e8d4e0b9SBarry Smith } 341e8d4e0b9SBarry Smith } 342e8d4e0b9SBarry Smith } 343289bc588SBarry Smith else { 344289bc588SBarry Smith for ( i=0; i<m; i++ ) { 345289bc588SBarry Smith for ( j=0; j<n; j++ ) { 346289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 347289bc588SBarry Smith } 348289bc588SBarry Smith } 349289bc588SBarry Smith } 350e8d4e0b9SBarry Smith } 351289bc588SBarry Smith return 0; 352289bc588SBarry Smith } 353e8d4e0b9SBarry Smith 354ae80bb75SLois Curfman McInnes static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v) 355ae80bb75SLois Curfman McInnes { 356ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 357ae80bb75SLois Curfman McInnes int i, j; 358ae80bb75SLois Curfman McInnes Scalar *vpt = v; 359ae80bb75SLois Curfman McInnes 360ae80bb75SLois Curfman McInnes /* row-oriented output */ 361ae80bb75SLois Curfman McInnes for ( i=0; i<m; i++ ) { 362ae80bb75SLois Curfman McInnes for ( j=0; j<n; j++ ) { 363ae80bb75SLois Curfman McInnes *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]]; 364ae80bb75SLois Curfman McInnes } 365ae80bb75SLois Curfman McInnes } 366ae80bb75SLois Curfman McInnes return 0; 367ae80bb75SLois Curfman McInnes } 368ae80bb75SLois Curfman McInnes 369289bc588SBarry Smith /* -----------------------------------------------------------------*/ 3703d1612f7SBarry Smith static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues) 371289bc588SBarry Smith { 372c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l; 373289bc588SBarry Smith int ierr; 374289bc588SBarry Smith Mat newi; 37548d91487SBarry Smith 376b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr); 377ec8511deSBarry Smith l = (Mat_SeqDense *) newi->data; 37855659b69SBarry Smith if (cpvalues == COPY_VALUES) { 379416022c9SBarry Smith PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 38055659b69SBarry Smith } 381227d817aSBarry Smith newi->assembled = PETSC_TRUE; 382289bc588SBarry Smith *newmat = newi; 383289bc588SBarry Smith return 0; 384289bc588SBarry Smith } 385289bc588SBarry Smith 38677c4ece6SBarry Smith #include "sys.h" 38788e32edaSLois Curfman McInnes 38819bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A) 38988e32edaSLois Curfman McInnes { 39088e32edaSLois Curfman McInnes Mat_SeqDense *a; 39188e32edaSLois Curfman McInnes Mat B; 39288e32edaSLois Curfman McInnes int *scols, i, j, nz, ierr, fd, header[4], size; 39388e32edaSLois Curfman McInnes int *rowlengths = 0, M, N, *cols; 39488e32edaSLois Curfman McInnes Scalar *vals, *svals, *v; 39519bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 39688e32edaSLois Curfman McInnes 39788e32edaSLois Curfman McInnes MPI_Comm_size(comm,&size); 39888e32edaSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor"); 39990ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 40077c4ece6SBarry Smith ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 40188e32edaSLois Curfman McInnes if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object"); 40288e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 40388e32edaSLois Curfman McInnes 40490ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 40590ace30eSBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 40690ace30eSBarry Smith B = *A; 40790ace30eSBarry Smith a = (Mat_SeqDense *) B->data; 40890ace30eSBarry Smith 40990ace30eSBarry Smith /* read in nonzero values */ 41077c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr); 41190ace30eSBarry Smith 41290ace30eSBarry Smith ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 41390ace30eSBarry Smith ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 41490ace30eSBarry Smith ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr); 41590ace30eSBarry Smith } else { 41688e32edaSLois Curfman McInnes /* read row lengths */ 4170452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 41877c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 41988e32edaSLois Curfman McInnes 42088e32edaSLois Curfman McInnes /* create our matrix */ 421b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 42288e32edaSLois Curfman McInnes B = *A; 42388e32edaSLois Curfman McInnes a = (Mat_SeqDense *) B->data; 42488e32edaSLois Curfman McInnes v = a->v; 42588e32edaSLois Curfman McInnes 42688e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 4270452661fSBarry Smith cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols); 42877c4ece6SBarry Smith ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr); 4290452661fSBarry Smith vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 43077c4ece6SBarry Smith ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 43188e32edaSLois Curfman McInnes 43288e32edaSLois Curfman McInnes /* insert into matrix */ 43388e32edaSLois Curfman McInnes for ( i=0; i<M; i++ ) { 43488e32edaSLois Curfman McInnes for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j]; 43588e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 43688e32edaSLois Curfman McInnes } 4370452661fSBarry Smith PetscFree(vals); PetscFree(cols); PetscFree(rowlengths); 43888e32edaSLois Curfman McInnes 43988e32edaSLois Curfman McInnes ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 44088e32edaSLois Curfman McInnes ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 44190ace30eSBarry Smith } 44288e32edaSLois Curfman McInnes return 0; 44388e32edaSLois Curfman McInnes } 44488e32edaSLois Curfman McInnes 445932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h" 44677c4ece6SBarry Smith #include "sys.h" 447932b0c3eSLois Curfman McInnes 448932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 449289bc588SBarry Smith { 450932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 451932b0c3eSLois Curfman McInnes int ierr, i, j, format; 452d636dbe3SBarry Smith FILE *fd; 453932b0c3eSLois Curfman McInnes char *outputname; 454932b0c3eSLois Curfman McInnes Scalar *v; 455932b0c3eSLois Curfman McInnes 45690ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 457932b0c3eSLois Curfman McInnes ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 45890ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 4597ddc982cSLois Curfman McInnes if (format == ASCII_FORMAT_INFO || format == ASCII_FORMAT_INFO_DETAILED) { 4607ddc982cSLois Curfman McInnes return 0; /* do nothing for now */ 461f72585f2SLois Curfman McInnes } 46244cd7ae7SLois Curfman McInnes else if (format == ASCII_FORMAT_COMMON) { 46380cd9d93SLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 46444cd7ae7SLois Curfman McInnes v = a->v + i; 46580cd9d93SLois Curfman McInnes fprintf(fd,"row %d:",i); 46680cd9d93SLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 46780cd9d93SLois Curfman McInnes #if defined(PETSC_COMPLEX) 46880cd9d93SLois Curfman McInnes if (real(*v) != 0.0 && imag(*v) != 0.0) 46980cd9d93SLois Curfman McInnes fprintf(fd," %d %g + %g i",j,real(*v),imag(*v)); 47080cd9d93SLois Curfman McInnes else if (real(*v)) fprintf(fd," %d %g ",j,real(*v)); 47180cd9d93SLois Curfman McInnes v += a->m; 47280cd9d93SLois Curfman McInnes #else 47380cd9d93SLois Curfman McInnes if (*v) fprintf(fd," %d %g ",j,*v); 47480cd9d93SLois Curfman McInnes v += a->m; 47580cd9d93SLois Curfman McInnes #endif 47680cd9d93SLois Curfman McInnes } 47780cd9d93SLois Curfman McInnes fprintf(fd,"\n"); 47880cd9d93SLois Curfman McInnes } 47980cd9d93SLois Curfman McInnes } 480f72585f2SLois Curfman McInnes else { 48147989497SBarry Smith #if defined(PETSC_COMPLEX) 48247989497SBarry Smith int allreal = 1; 48347989497SBarry Smith /* determine if matrix has all real values */ 48447989497SBarry Smith v = a->v; 48547989497SBarry Smith for ( i=0; i<a->m*a->n; i++ ) { 48647989497SBarry Smith if (imag(v[i])) { allreal = 0; break ;} 48747989497SBarry Smith } 48847989497SBarry Smith #endif 489932b0c3eSLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 490932b0c3eSLois Curfman McInnes v = a->v + i; 491932b0c3eSLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 492289bc588SBarry Smith #if defined(PETSC_COMPLEX) 49347989497SBarry Smith if (allreal) { 49447989497SBarry Smith fprintf(fd,"%6.4e ",real(*v)); v += a->m; 49547989497SBarry Smith } else { 496932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m; 49747989497SBarry Smith } 498289bc588SBarry Smith #else 499932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e ",*v); v += a->m; 500289bc588SBarry Smith #endif 501289bc588SBarry Smith } 5028e37dc09SBarry Smith fprintf(fd,"\n"); 503289bc588SBarry Smith } 504da3a660dSBarry Smith } 505932b0c3eSLois Curfman McInnes fflush(fd); 506289bc588SBarry Smith return 0; 507289bc588SBarry Smith } 508289bc588SBarry Smith 509932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 510932b0c3eSLois Curfman McInnes { 511932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 512932b0c3eSLois Curfman McInnes int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 51390ace30eSBarry Smith int format; 51490ace30eSBarry Smith Scalar *v, *anonz,*vals; 515932b0c3eSLois Curfman McInnes 51690ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 51790ace30eSBarry Smith 51890ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 51990ace30eSBarry Smith if (format == BINARY_FORMAT_NATIVE) { 52090ace30eSBarry Smith /* store the matrix as a dense matrix */ 52190ace30eSBarry Smith col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens); 52290ace30eSBarry Smith col_lens[0] = MAT_COOKIE; 52390ace30eSBarry Smith col_lens[1] = m; 52490ace30eSBarry Smith col_lens[2] = n; 52590ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 52677c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr); 52790ace30eSBarry Smith PetscFree(col_lens); 52890ace30eSBarry Smith 52990ace30eSBarry Smith /* write out matrix, by rows */ 53090ace30eSBarry Smith vals = (Scalar *) PetscMalloc(m*n*sizeof(Scalar)); CHKPTRQ(vals); 53190ace30eSBarry Smith v = a->v; 53290ace30eSBarry Smith for ( i=0; i<m; i++ ) { 53390ace30eSBarry Smith for ( j=0; j<n; j++ ) { 53490ace30eSBarry Smith vals[i + j*m] = *v++; 53590ace30eSBarry Smith } 53690ace30eSBarry Smith } 53777c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr); 53890ace30eSBarry Smith PetscFree(vals); 53990ace30eSBarry Smith } else { 5400452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens); 541932b0c3eSLois Curfman McInnes col_lens[0] = MAT_COOKIE; 542932b0c3eSLois Curfman McInnes col_lens[1] = m; 543932b0c3eSLois Curfman McInnes col_lens[2] = n; 544932b0c3eSLois Curfman McInnes col_lens[3] = nz; 545932b0c3eSLois Curfman McInnes 546932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 547932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) col_lens[4+i] = n; 54877c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr); 549932b0c3eSLois Curfman McInnes 550932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 551932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 552932b0c3eSLois Curfman McInnes ict = 0; 553932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 554932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) col_lens[ict++] = j; 555932b0c3eSLois Curfman McInnes } 55677c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr); 5570452661fSBarry Smith PetscFree(col_lens); 558932b0c3eSLois Curfman McInnes 559932b0c3eSLois Curfman McInnes /* store nonzero values */ 5600452661fSBarry Smith anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz); 561932b0c3eSLois Curfman McInnes ict = 0; 562932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 563932b0c3eSLois Curfman McInnes v = a->v + i; 564932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) { 565932b0c3eSLois Curfman McInnes anonz[ict++] = *v; v += a->m; 566932b0c3eSLois Curfman McInnes } 567932b0c3eSLois Curfman McInnes } 56877c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr); 5690452661fSBarry Smith PetscFree(anonz); 57090ace30eSBarry Smith } 571932b0c3eSLois Curfman McInnes return 0; 572932b0c3eSLois Curfman McInnes } 573932b0c3eSLois Curfman McInnes 574932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer) 575932b0c3eSLois Curfman McInnes { 576932b0c3eSLois Curfman McInnes Mat A = (Mat) obj; 577932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*) A->data; 578bcd2baecSBarry Smith ViewerType vtype; 579bcd2baecSBarry Smith int ierr; 580932b0c3eSLois Curfman McInnes 581932b0c3eSLois Curfman McInnes if (!viewer) { 58219bcc07fSBarry Smith viewer = STDOUT_VIEWER_SELF; 583932b0c3eSLois Curfman McInnes } 584bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 585bcd2baecSBarry Smith 586bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 587932b0c3eSLois Curfman McInnes return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); 588932b0c3eSLois Curfman McInnes } 58919bcc07fSBarry Smith else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 590932b0c3eSLois Curfman McInnes return MatView_SeqDense_ASCII(A,viewer); 591932b0c3eSLois Curfman McInnes } 592bcd2baecSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 593932b0c3eSLois Curfman McInnes return MatView_SeqDense_Binary(A,viewer); 594932b0c3eSLois Curfman McInnes } 595932b0c3eSLois Curfman McInnes return 0; 596932b0c3eSLois Curfman McInnes } 597289bc588SBarry Smith 598ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj) 599289bc588SBarry Smith { 600289bc588SBarry Smith Mat mat = (Mat) obj; 601ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 602a5a9c739SBarry Smith #if defined(PETSC_LOG) 603a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 604a5a9c739SBarry Smith #endif 6050452661fSBarry Smith if (l->pivots) PetscFree(l->pivots); 6063439631bSBarry Smith if (!l->user_alloc) PetscFree(l->v); 6070452661fSBarry Smith PetscFree(l); 608a5a9c739SBarry Smith PLogObjectDestroy(mat); 6090452661fSBarry Smith PetscHeaderDestroy(mat); 610289bc588SBarry Smith return 0; 611289bc588SBarry Smith } 612289bc588SBarry Smith 613c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout) 614289bc588SBarry Smith { 615c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 616d3e5ee88SLois Curfman McInnes int k, j, m, n; 617d3e5ee88SLois Curfman McInnes Scalar *v, tmp; 61848b35521SBarry Smith 619d3e5ee88SLois Curfman McInnes v = mat->v; m = mat->m; n = mat->n; 6203638b69dSLois Curfman McInnes if (matout == PETSC_NULL) { /* in place transpose */ 621d9f96c7cSLois Curfman McInnes if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place"); 622d3e5ee88SLois Curfman McInnes for ( j=0; j<m; j++ ) { 623289bc588SBarry Smith for ( k=0; k<j; k++ ) { 624d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 625d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 626d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 627289bc588SBarry Smith } 628289bc588SBarry Smith } 62948b35521SBarry Smith } 630d3e5ee88SLois Curfman McInnes else { /* out-of-place transpose */ 631d3e5ee88SLois Curfman McInnes int ierr; 632d3e5ee88SLois Curfman McInnes Mat tmat; 633ec8511deSBarry Smith Mat_SeqDense *tmatd; 634d3e5ee88SLois Curfman McInnes Scalar *v2; 635b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr); 636ec8511deSBarry Smith tmatd = (Mat_SeqDense *) tmat->data; 6370de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 638d3e5ee88SLois Curfman McInnes for ( j=0; j<n; j++ ) { 6390de55854SLois Curfman McInnes for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 640d3e5ee88SLois Curfman McInnes } 641d3e5ee88SLois Curfman McInnes ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 642d3e5ee88SLois Curfman McInnes ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 643d3e5ee88SLois Curfman McInnes *matout = tmat; 64448b35521SBarry Smith } 645289bc588SBarry Smith return 0; 646289bc588SBarry Smith } 647289bc588SBarry Smith 64877c4ece6SBarry Smith static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg) 649289bc588SBarry Smith { 650c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 651c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 652289bc588SBarry Smith int i; 653289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 6549ea5d5aeSSatish Balay 6559ea5d5aeSSatish Balay if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type MATSEQDENSE"); 65677c4ece6SBarry Smith if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;} 65777c4ece6SBarry Smith if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;} 658289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 65977c4ece6SBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;} 660289bc588SBarry Smith v1++; v2++; 661289bc588SBarry Smith } 66277c4ece6SBarry Smith *flg = PETSC_TRUE; 6639ea5d5aeSSatish Balay return 0; 664289bc588SBarry Smith } 665289bc588SBarry Smith 666c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v) 667289bc588SBarry Smith { 668c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 66944cd7ae7SLois Curfman McInnes int i, n, len; 67044cd7ae7SLois Curfman McInnes Scalar *x, zero = 0.0; 67144cd7ae7SLois Curfman McInnes 67244cd7ae7SLois Curfman McInnes VecSet(&zero,v); 673289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 67444cd7ae7SLois Curfman McInnes len = PetscMin(mat->m,mat->n); 675ec8511deSBarry Smith if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 67644cd7ae7SLois Curfman McInnes for ( i=0; i<len; i++ ) { 677289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 678289bc588SBarry Smith } 679289bc588SBarry Smith return 0; 680289bc588SBarry Smith } 681289bc588SBarry Smith 682052efed2SBarry Smith static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 683289bc588SBarry Smith { 684c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 685da3a660dSBarry Smith Scalar *l,*r,x,*v; 686da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 68755659b69SBarry Smith 68828988994SBarry Smith if (ll) { 689da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 690052efed2SBarry Smith if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size"); 691da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 692da3a660dSBarry Smith x = l[i]; 693da3a660dSBarry Smith v = mat->v + i; 694da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 695da3a660dSBarry Smith } 69644cd7ae7SLois Curfman McInnes PLogFlops(n*m); 697da3a660dSBarry Smith } 69828988994SBarry Smith if (rr) { 699da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 700052efed2SBarry Smith if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size"); 701da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 702da3a660dSBarry Smith x = r[i]; 703da3a660dSBarry Smith v = mat->v + i*m; 704da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 705da3a660dSBarry Smith } 70644cd7ae7SLois Curfman McInnes PLogFlops(n*m); 707da3a660dSBarry Smith } 708289bc588SBarry Smith return 0; 709289bc588SBarry Smith } 710289bc588SBarry Smith 711cddf8d76SBarry Smith static int MatNorm_SeqDense(Mat A,NormType type,double *norm) 712289bc588SBarry Smith { 713c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 714289bc588SBarry Smith Scalar *v = mat->v; 715289bc588SBarry Smith double sum = 0.0; 716289bc588SBarry Smith int i, j; 71755659b69SBarry Smith 718289bc588SBarry Smith if (type == NORM_FROBENIUS) { 719289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 720289bc588SBarry Smith #if defined(PETSC_COMPLEX) 721289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 722289bc588SBarry Smith #else 723289bc588SBarry Smith sum += (*v)*(*v); v++; 724289bc588SBarry Smith #endif 725289bc588SBarry Smith } 726289bc588SBarry Smith *norm = sqrt(sum); 72755659b69SBarry Smith PLogFlops(2*mat->n*mat->m); 728289bc588SBarry Smith } 729289bc588SBarry Smith else if (type == NORM_1) { 730289bc588SBarry Smith *norm = 0.0; 731289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 732289bc588SBarry Smith sum = 0.0; 733289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 73433a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 735289bc588SBarry Smith } 736289bc588SBarry Smith if (sum > *norm) *norm = sum; 737289bc588SBarry Smith } 73855659b69SBarry Smith PLogFlops(mat->n*mat->m); 739289bc588SBarry Smith } 740289bc588SBarry Smith else if (type == NORM_INFINITY) { 741289bc588SBarry Smith *norm = 0.0; 742289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 743289bc588SBarry Smith v = mat->v + j; 744289bc588SBarry Smith sum = 0.0; 745289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 746cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v += mat->m; 747289bc588SBarry Smith } 748289bc588SBarry Smith if (sum > *norm) *norm = sum; 749289bc588SBarry Smith } 75055659b69SBarry Smith PLogFlops(mat->n*mat->m); 751289bc588SBarry Smith } 752289bc588SBarry Smith else { 75348d91487SBarry Smith SETERRQ(1,"MatNorm_SeqDense:No two norm"); 754289bc588SBarry Smith } 755289bc588SBarry Smith return 0; 756289bc588SBarry Smith } 757289bc588SBarry Smith 758c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op) 759289bc588SBarry Smith { 760c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 76167e560aaSBarry Smith 762289bc588SBarry Smith if (op == ROW_ORIENTED) aij->roworiented = 1; 763289bc588SBarry Smith else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 764c0bbcb79SLois Curfman McInnes else if (op == ROWS_SORTED || 76588e32edaSLois Curfman McInnes op == COLUMNS_SORTED || 766c0bbcb79SLois Curfman McInnes op == SYMMETRIC_MATRIX || 767c0bbcb79SLois Curfman McInnes op == STRUCTURALLY_SYMMETRIC_MATRIX || 768c0bbcb79SLois Curfman McInnes op == NO_NEW_NONZERO_LOCATIONS || 769c0bbcb79SLois Curfman McInnes op == YES_NEW_NONZERO_LOCATIONS || 770c0bbcb79SLois Curfman McInnes op == NO_NEW_DIAGONALS || 771c0bbcb79SLois Curfman McInnes op == YES_NEW_DIAGONALS) 77294a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n"); 773c0bbcb79SLois Curfman McInnes else 7744d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");} 775289bc588SBarry Smith return 0; 776289bc588SBarry Smith } 777289bc588SBarry Smith 778ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A) 7796f0a148fSBarry Smith { 780ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 781cddf8d76SBarry Smith PetscMemzero(l->v,l->m*l->n*sizeof(Scalar)); 7826f0a148fSBarry Smith return 0; 7836f0a148fSBarry Smith } 7846f0a148fSBarry Smith 785ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 7866f0a148fSBarry Smith { 787ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 7886abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 7896f0a148fSBarry Smith Scalar *slot; 79055659b69SBarry Smith 79177c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 79278b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 7936f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 7946f0a148fSBarry Smith slot = l->v + rows[i]; 7956f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 7966f0a148fSBarry Smith } 7976f0a148fSBarry Smith if (diag) { 7986f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 7996f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 800260925b8SBarry Smith *slot = *diag; 8016f0a148fSBarry Smith } 8026f0a148fSBarry Smith } 803260925b8SBarry Smith ISRestoreIndices(is,&rows); 8046f0a148fSBarry Smith return 0; 8056f0a148fSBarry Smith } 806557bce09SLois Curfman McInnes 807c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n) 808557bce09SLois Curfman McInnes { 809c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 810557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 811557bce09SLois Curfman McInnes return 0; 812557bce09SLois Curfman McInnes } 813557bce09SLois Curfman McInnes 814c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 815d3e5ee88SLois Curfman McInnes { 816c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 817d3e5ee88SLois Curfman McInnes *m = 0; *n = mat->m; 818d3e5ee88SLois Curfman McInnes return 0; 819d3e5ee88SLois Curfman McInnes } 820d3e5ee88SLois Curfman McInnes 821c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array) 82264e87e97SBarry Smith { 823c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 82464e87e97SBarry Smith *array = mat->v; 82564e87e97SBarry Smith return 0; 82664e87e97SBarry Smith } 8270754003eSLois Curfman McInnes 828ff14e315SSatish Balay static int MatRestoreArray_SeqDense(Mat A,Scalar **array) 829ff14e315SSatish Balay { 830ff14e315SSatish Balay return 0; 831ff14e315SSatish Balay } 8320754003eSLois Curfman McInnes 833c0bbcb79SLois Curfman McInnes static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol) 8340754003eSLois Curfman McInnes { 835ec8511deSBarry Smith SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done"); 8360754003eSLois Curfman McInnes } 8370754003eSLois Curfman McInnes 838cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall, 839cddf8d76SBarry Smith Mat *submat) 8400754003eSLois Curfman McInnes { 841c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 8420754003eSLois Curfman McInnes int nznew, *smap, i, j, ierr, oldcols = mat->n; 843160018dcSBarry Smith int *irow, *icol, nrows, ncols, *cwork; 8440754003eSLois Curfman McInnes Scalar *vwork, *val; 8450754003eSLois Curfman McInnes Mat newmat; 8460754003eSLois Curfman McInnes 84778b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 84878b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 84978b31e54SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 85078b31e54SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 8510754003eSLois Curfman McInnes 8520452661fSBarry Smith smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap); 8530452661fSBarry Smith cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork); 8540452661fSBarry Smith vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 855cddf8d76SBarry Smith PetscMemzero((char*)smap,oldcols*sizeof(int)); 8560754003eSLois Curfman McInnes for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 8570754003eSLois Curfman McInnes 8580754003eSLois Curfman McInnes /* Create and fill new matrix */ 859b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr); 8600754003eSLois Curfman McInnes for (i=0; i<nrows; i++) { 8610754003eSLois Curfman McInnes nznew = 0; 8620754003eSLois Curfman McInnes val = mat->v + irow[i]; 8630754003eSLois Curfman McInnes for (j=0; j<oldcols; j++) { 8640754003eSLois Curfman McInnes if (smap[j]) { 8650754003eSLois Curfman McInnes cwork[nznew] = smap[j] - 1; 8660754003eSLois Curfman McInnes vwork[nznew++] = val[j * mat->m]; 8670754003eSLois Curfman McInnes } 8680754003eSLois Curfman McInnes } 869dbb450caSBarry Smith ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES); 87078b31e54SBarry Smith CHKERRQ(ierr); 8710754003eSLois Curfman McInnes } 87278b31e54SBarry Smith ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 87378b31e54SBarry Smith ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 8740754003eSLois Curfman McInnes 8750754003eSLois Curfman McInnes /* Free work space */ 8760452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 87778b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 87878b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 8790754003eSLois Curfman McInnes *submat = newmat; 8800754003eSLois Curfman McInnes return 0; 8810754003eSLois Curfman McInnes } 8820754003eSLois Curfman McInnes 8834b0e389bSBarry Smith static int MatCopy_SeqDense(Mat A, Mat B) 8844b0e389bSBarry Smith { 8854b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data; 8864b0e389bSBarry Smith if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B); 8874b0e389bSBarry Smith if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)"); 8884b0e389bSBarry Smith PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar)); 8894b0e389bSBarry Smith return 0; 8904b0e389bSBarry Smith } 8914b0e389bSBarry Smith 892289bc588SBarry Smith /* -------------------------------------------------------------------*/ 893ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense, 894ec8511deSBarry Smith MatGetRow_SeqDense, MatRestoreRow_SeqDense, 895ec8511deSBarry Smith MatMult_SeqDense, MatMultAdd_SeqDense, 896ec8511deSBarry Smith MatMultTrans_SeqDense, MatMultTransAdd_SeqDense, 897ec8511deSBarry Smith MatSolve_SeqDense,MatSolveAdd_SeqDense, 898ec8511deSBarry Smith MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense, 899ec8511deSBarry Smith MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense, 900ec8511deSBarry Smith MatRelax_SeqDense, 901ec8511deSBarry Smith MatTranspose_SeqDense, 902ec8511deSBarry Smith MatGetInfo_SeqDense,MatEqual_SeqDense, 903052efed2SBarry Smith MatGetDiagonal_SeqDense,MatDiagonalScale_SeqDense,MatNorm_SeqDense, 904289bc588SBarry Smith 0,0, 905ec8511deSBarry Smith 0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0, 906ec8511deSBarry Smith MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense, 907ec8511deSBarry Smith MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense, 908ec8511deSBarry Smith MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense, 909ff14e315SSatish Balay 0,0, MatGetArray_SeqDense, MatRestoreArray_SeqDense,0, 910ec8511deSBarry Smith MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense, 9113d1612f7SBarry Smith MatConvertSameType_SeqDense,0,0,0,0, 912ae80bb75SLois Curfman McInnes MatAXPY_SeqDense,0,0, 9134b0e389bSBarry Smith MatGetValues_SeqDense, 91480cd9d93SLois Curfman McInnes MatCopy_SeqDense,0,MatScale_SeqDense}; 9150754003eSLois Curfman McInnes 91690ace30eSBarry Smith 9174b828684SBarry Smith /*@C 918fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 919d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 920d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 921289bc588SBarry Smith 92220563c6bSBarry Smith Input Parameters: 9230c775827SLois Curfman McInnes . comm - MPI communicator, set to MPI_COMM_SELF 9240c775827SLois Curfman McInnes . m - number of rows 92518f449edSLois Curfman McInnes . n - number of columns 926b4fd4287SBarry Smith . data - optional location of matrix data. Set data=PETSC_NULL for PETSc 927dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 92820563c6bSBarry Smith 92920563c6bSBarry Smith Output Parameter: 93044cd7ae7SLois Curfman McInnes . A - the matrix 93120563c6bSBarry Smith 93218f449edSLois Curfman McInnes Notes: 93318f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 93418f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 935b4fd4287SBarry Smith set data=PETSC_NULL. 93618f449edSLois Curfman McInnes 937dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 938d65003e9SLois Curfman McInnes 939dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues() 94020563c6bSBarry Smith @*/ 94144cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A) 942289bc588SBarry Smith { 94344cd7ae7SLois Curfman McInnes Mat B; 94444cd7ae7SLois Curfman McInnes Mat_SeqDense *b; 94525cdf11fSBarry Smith int ierr,flg; 94655659b69SBarry Smith 94744cd7ae7SLois Curfman McInnes *A = 0; 94844cd7ae7SLois Curfman McInnes PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQDENSE,comm); 94944cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 95044cd7ae7SLois Curfman McInnes b = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b); 95144cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqDense)); 95244cd7ae7SLois Curfman McInnes PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 95344cd7ae7SLois Curfman McInnes B->destroy = MatDestroy_SeqDense; 95444cd7ae7SLois Curfman McInnes B->view = MatView_SeqDense; 95544cd7ae7SLois Curfman McInnes B->factor = 0; 95644cd7ae7SLois Curfman McInnes PLogObjectMemory(B,sizeof(struct _Mat)); 95744cd7ae7SLois Curfman McInnes B->data = (void *) b; 95818f449edSLois Curfman McInnes 95944cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 96044cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 96144cd7ae7SLois Curfman McInnes b->pivots = 0; 96244cd7ae7SLois Curfman McInnes b->roworiented = 1; 963b4fd4287SBarry Smith if (data == PETSC_NULL) { 96444cd7ae7SLois Curfman McInnes b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v); 96544cd7ae7SLois Curfman McInnes PetscMemzero(b->v,m*n*sizeof(Scalar)); 96644cd7ae7SLois Curfman McInnes b->user_alloc = 0; 96744cd7ae7SLois Curfman McInnes PLogObjectMemory(B,n*m*sizeof(Scalar)); 96818f449edSLois Curfman McInnes } 9692dd2b441SLois Curfman McInnes else { /* user-allocated storage */ 97044cd7ae7SLois Curfman McInnes b->v = data; 97144cd7ae7SLois Curfman McInnes b->user_alloc = 1; 9722dd2b441SLois Curfman McInnes } 97325cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 97444cd7ae7SLois Curfman McInnes if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);} 97544cd7ae7SLois Curfman McInnes *A = B; 976289bc588SBarry Smith return 0; 977289bc588SBarry Smith } 978289bc588SBarry Smith 979c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat) 980289bc588SBarry Smith { 981c0bbcb79SLois Curfman McInnes Mat_SeqDense *m = (Mat_SeqDense *) A->data; 982b4fd4287SBarry Smith return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat); 983289bc588SBarry Smith } 984