1cb512458SBarry Smith #ifndef lint 2*639f9d9dSBarry Smith static char vcid[] = "$Id: dense.c,v 1.112 1996/10/01 03:34:42 bsmith Exp bsmith $"; 3cb512458SBarry Smith #endif 467e560aaSBarry Smith /* 567e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 667e560aaSBarry Smith */ 7289bc588SBarry Smith 870f55243SBarry Smith #include "src/mat/impls/dense/seq/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 214e220ebcSLois Curfman McInnes static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 22289bc588SBarry Smith { 234e220ebcSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 244e220ebcSLois Curfman McInnes int i,N = a->m*a->n,count = 0; 254e220ebcSLois Curfman McInnes Scalar *v = a->v; 26289bc588SBarry Smith for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 274e220ebcSLois Curfman McInnes 284e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 294e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 304e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 314e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 324e220ebcSLois Curfman McInnes info->block_size = 1.0; 334e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 344e220ebcSLois Curfman McInnes info->nz_used = (double)count; 354e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(N-count); 364e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 374e220ebcSLois Curfman McInnes info->mallocs = 0; 384e220ebcSLois Curfman McInnes info->memory = A->mem; 394e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 404e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 414e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 424e220ebcSLois Curfman McInnes 43fa9ec3f1SLois Curfman McInnes return 0; 44289bc588SBarry Smith } 45289bc588SBarry Smith 4680cd9d93SLois Curfman McInnes static int MatScale_SeqDense(Scalar *alpha,Mat inA) 4780cd9d93SLois Curfman McInnes { 4880cd9d93SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) inA->data; 4980cd9d93SLois Curfman McInnes int one = 1, nz; 5080cd9d93SLois Curfman McInnes 5180cd9d93SLois Curfman McInnes nz = a->m*a->n; 5280cd9d93SLois Curfman McInnes BLscal_( &nz, alpha, a->v, &one ); 5380cd9d93SLois Curfman McInnes PLogFlops(nz); 5480cd9d93SLois Curfman McInnes return 0; 5580cd9d93SLois Curfman McInnes } 5680cd9d93SLois Curfman McInnes 57289bc588SBarry Smith /* ---------------------------------------------------------------*/ 58289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots, 59289bc588SBarry Smith rather than put it in the Mat->col slot.*/ 60c0bbcb79SLois Curfman McInnes static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f) 61289bc588SBarry Smith { 62c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 636abc6512SBarry Smith int info; 6455659b69SBarry Smith 65289bc588SBarry Smith if (!mat->pivots) { 660452661fSBarry Smith mat->pivots = (int *) PetscMalloc(mat->m*sizeof(int));CHKPTRQ(mat->pivots); 67c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,mat->m*sizeof(int)); 68289bc588SBarry Smith } 69289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 703d60af36SSatish Balay if (info<0) SETERRQ(1,"MatLUFactor_SeqDense:Bad argument to LU factorization"); 713d60af36SSatish Balay if (info>0) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization"); 72c0bbcb79SLois Curfman McInnes A->factor = FACTOR_LU; 7355659b69SBarry Smith PLogFlops((2*mat->n*mat->n*mat->n)/3); 74289bc588SBarry Smith return 0; 75289bc588SBarry Smith } 7602cad45dSBarry Smith static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues) 7702cad45dSBarry Smith { 7802cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l; 7902cad45dSBarry Smith int ierr; 8002cad45dSBarry Smith Mat newi; 8102cad45dSBarry Smith 8202cad45dSBarry Smith ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr); 8302cad45dSBarry Smith l = (Mat_SeqDense *) newi->data; 8402cad45dSBarry Smith if (cpvalues == COPY_VALUES) { 8502cad45dSBarry Smith PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 8602cad45dSBarry Smith } 8702cad45dSBarry Smith newi->assembled = PETSC_TRUE; 8802cad45dSBarry Smith *newmat = newi; 8902cad45dSBarry Smith return 0; 9002cad45dSBarry Smith } 9102cad45dSBarry Smith 92c0bbcb79SLois Curfman McInnes static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact) 93289bc588SBarry Smith { 9402cad45dSBarry Smith return MatConvertSameType_SeqDense(A,fact,PETSC_FALSE); 95289bc588SBarry Smith } 96c0bbcb79SLois Curfman McInnes static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 97289bc588SBarry Smith { 9802cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*) A->data, *l = (Mat_SeqDense*) (*fact)->data; 9902cad45dSBarry Smith /* copy the numerical values */ 10002cad45dSBarry Smith PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 10102cad45dSBarry Smith (*fact)->factor = 0; 10249d8b64dSBarry Smith return MatLUFactor(*fact,0,0,1.0); 103289bc588SBarry Smith } 104c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact) 105289bc588SBarry Smith { 106a57ad546SLois Curfman McInnes return MatConvert(A,MATSAME,fact); 107289bc588SBarry Smith } 108c0bbcb79SLois Curfman McInnes static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 109289bc588SBarry Smith { 11049d8b64dSBarry Smith return MatCholeskyFactor(*fact,0,1.0); 111289bc588SBarry Smith } 112c0bbcb79SLois Curfman McInnes static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f) 113289bc588SBarry Smith { 114c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1156abc6512SBarry Smith int info; 11655659b69SBarry Smith 117752f5567SLois Curfman McInnes if (mat->pivots) { 1180452661fSBarry Smith PetscFree(mat->pivots); 119c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,-mat->m*sizeof(int)); 120752f5567SLois Curfman McInnes mat->pivots = 0; 121752f5567SLois Curfman McInnes } 122289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 123ec8511deSBarry Smith if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization"); 124c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 12555659b69SBarry Smith PLogFlops((mat->n*mat->n*mat->n)/3); 126289bc588SBarry Smith return 0; 127289bc588SBarry Smith } 128289bc588SBarry Smith 129c0bbcb79SLois Curfman McInnes static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 130289bc588SBarry Smith { 131c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 132a57ad546SLois Curfman McInnes int one = 1, info, ierr; 1336abc6512SBarry Smith Scalar *x, *y; 13467e560aaSBarry Smith 135a57ad546SLois Curfman McInnes ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 136416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 137c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 13848d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 139289bc588SBarry Smith } 140c0bbcb79SLois Curfman McInnes else if (A->factor == FACTOR_CHOLESKY){ 14148d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 142289bc588SBarry Smith } 143ec8511deSBarry Smith else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve"); 144ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve"); 14555659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 146289bc588SBarry Smith return 0; 147289bc588SBarry Smith } 148c0bbcb79SLois Curfman McInnes static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy) 149da3a660dSBarry Smith { 150c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1516abc6512SBarry Smith int one = 1, info; 1526abc6512SBarry Smith Scalar *x, *y; 15367e560aaSBarry Smith 154da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 155416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 156752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 157da3a660dSBarry Smith if (mat->pivots) { 15848d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 159da3a660dSBarry Smith } 160da3a660dSBarry Smith else { 16148d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 162da3a660dSBarry Smith } 163ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve"); 16455659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 165da3a660dSBarry Smith return 0; 166da3a660dSBarry Smith } 167c0bbcb79SLois Curfman McInnes static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 168da3a660dSBarry Smith { 169c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1706abc6512SBarry Smith int one = 1, info,ierr; 1716abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 172da3a660dSBarry Smith Vec tmp = 0; 17367e560aaSBarry Smith 174da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 175da3a660dSBarry Smith if (yy == zz) { 17678b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 177c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 17878b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 179da3a660dSBarry Smith } 180416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 181752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 182da3a660dSBarry Smith if (mat->pivots) { 18348d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 184da3a660dSBarry Smith } 185da3a660dSBarry Smith else { 18648d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 187da3a660dSBarry Smith } 188ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve"); 189da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 190da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 19155659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 192da3a660dSBarry Smith return 0; 193da3a660dSBarry Smith } 19467e560aaSBarry Smith 195c0bbcb79SLois Curfman McInnes static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 196da3a660dSBarry Smith { 197c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1986abc6512SBarry Smith int one = 1, info,ierr; 1996abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 200da3a660dSBarry Smith Vec tmp; 20167e560aaSBarry Smith 202da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 203da3a660dSBarry Smith if (yy == zz) { 20478b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 205c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 20678b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 207da3a660dSBarry Smith } 208416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 209752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 210da3a660dSBarry Smith if (mat->pivots) { 21148d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 212da3a660dSBarry Smith } 213da3a660dSBarry Smith else { 21448d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 215da3a660dSBarry Smith } 216ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve"); 217da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 218da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 21955659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 220da3a660dSBarry Smith return 0; 221da3a660dSBarry Smith } 222289bc588SBarry Smith /* ------------------------------------------------------------------*/ 223c0bbcb79SLois Curfman McInnes static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 22420e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 225289bc588SBarry Smith { 226c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 227289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 2286abc6512SBarry Smith int o = 1,ierr, m = mat->m, i; 229289bc588SBarry Smith 230289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 231289bc588SBarry Smith /* this is a hack fix, should have another version without 232289bc588SBarry Smith the second BLdot */ 233bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx); CHKERRQ(ierr); 234289bc588SBarry Smith } 235289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 236289bc588SBarry Smith while (its--) { 237289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 238289bc588SBarry Smith for ( i=0; i<m; i++ ) { 239f1747703SBarry Smith #if defined(PETSC_COMPLEX) 240f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 241f1747703SBarry Smith not happy about returning a double complex */ 242f1747703SBarry Smith int _i; 243f1747703SBarry Smith Scalar sum = b[i]; 244f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 245f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 246f1747703SBarry Smith } 247f1747703SBarry Smith xt = sum; 248f1747703SBarry Smith #else 249289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 250f1747703SBarry Smith #endif 25155a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 252289bc588SBarry Smith } 253289bc588SBarry Smith } 254289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 255289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 256f1747703SBarry Smith #if defined(PETSC_COMPLEX) 257f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 258f1747703SBarry Smith not happy about returning a double complex */ 259f1747703SBarry Smith int _i; 260f1747703SBarry Smith Scalar sum = b[i]; 261f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 262f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 263f1747703SBarry Smith } 264f1747703SBarry Smith xt = sum; 265f1747703SBarry Smith #else 266289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 267f1747703SBarry Smith #endif 26855a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 269289bc588SBarry Smith } 270289bc588SBarry Smith } 271289bc588SBarry Smith } 272289bc588SBarry Smith return 0; 273289bc588SBarry Smith } 274289bc588SBarry Smith 275289bc588SBarry Smith /* -----------------------------------------------------------------*/ 27644cd7ae7SLois Curfman McInnes int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy) 277289bc588SBarry Smith { 278c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 279289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 280289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 281289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 28248d91487SBarry Smith LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 28355659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->n); 284289bc588SBarry Smith return 0; 285289bc588SBarry Smith } 28644cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 287289bc588SBarry Smith { 288c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 289289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 290289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 291289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 29248d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 29355659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->m); 294289bc588SBarry Smith return 0; 295289bc588SBarry Smith } 29644cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 297289bc588SBarry Smith { 298c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 299289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 3006abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 301289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 302416022c9SBarry Smith if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 30348d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 30455659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 305289bc588SBarry Smith return 0; 306289bc588SBarry Smith } 30744cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 308289bc588SBarry Smith { 309c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 310289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 311289bc588SBarry Smith int _One=1; 3126abc6512SBarry Smith Scalar _DOne=1.0; 31344cd7ae7SLois Curfman McInnes VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 314717eeb74SLois Curfman McInnes if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar)); 31548d91487SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 31655659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 317289bc588SBarry Smith return 0; 318289bc588SBarry Smith } 319289bc588SBarry Smith 320289bc588SBarry Smith /* -----------------------------------------------------------------*/ 321c0bbcb79SLois Curfman McInnes static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 322289bc588SBarry Smith { 323c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 324289bc588SBarry Smith Scalar *v; 325289bc588SBarry Smith int i; 32667e560aaSBarry Smith 327289bc588SBarry Smith *ncols = mat->n; 328289bc588SBarry Smith if (cols) { 3290452661fSBarry Smith *cols = (int *) PetscMalloc(mat->n*sizeof(int)); CHKPTRQ(*cols); 33073c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 331289bc588SBarry Smith } 332289bc588SBarry Smith if (vals) { 3330452661fSBarry Smith *vals = (Scalar *) PetscMalloc(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 334289bc588SBarry Smith v = mat->v + row; 33573c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 336289bc588SBarry Smith } 337289bc588SBarry Smith return 0; 338289bc588SBarry Smith } 339c0bbcb79SLois Curfman McInnes static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 340289bc588SBarry Smith { 3410452661fSBarry Smith if (cols) { PetscFree(*cols); } 3420452661fSBarry Smith if (vals) { PetscFree(*vals); } 343289bc588SBarry Smith return 0; 344289bc588SBarry Smith } 345289bc588SBarry Smith /* ----------------------------------------------------------------*/ 346ae80bb75SLois Curfman McInnes static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 347e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 348289bc588SBarry Smith { 349c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 350289bc588SBarry Smith int i,j; 351d6dfbf8fSBarry Smith 352289bc588SBarry Smith if (!mat->roworiented) { 353dbb450caSBarry Smith if (addv == INSERT_VALUES) { 354289bc588SBarry Smith for ( j=0; j<n; j++ ) { 355289bc588SBarry Smith for ( i=0; i<m; i++ ) { 356289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 357289bc588SBarry Smith } 358289bc588SBarry Smith } 359289bc588SBarry Smith } 360289bc588SBarry Smith else { 361289bc588SBarry Smith for ( j=0; j<n; j++ ) { 362289bc588SBarry Smith for ( i=0; i<m; i++ ) { 363289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 364289bc588SBarry Smith } 365289bc588SBarry Smith } 366289bc588SBarry Smith } 367e8d4e0b9SBarry Smith } 368e8d4e0b9SBarry Smith else { 369dbb450caSBarry Smith if (addv == INSERT_VALUES) { 370e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 371e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 372e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 373e8d4e0b9SBarry Smith } 374e8d4e0b9SBarry Smith } 375e8d4e0b9SBarry Smith } 376289bc588SBarry Smith else { 377289bc588SBarry Smith for ( i=0; i<m; i++ ) { 378289bc588SBarry Smith for ( j=0; j<n; j++ ) { 379289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 380289bc588SBarry Smith } 381289bc588SBarry Smith } 382289bc588SBarry Smith } 383e8d4e0b9SBarry Smith } 384289bc588SBarry Smith return 0; 385289bc588SBarry Smith } 386e8d4e0b9SBarry Smith 387ae80bb75SLois Curfman McInnes static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v) 388ae80bb75SLois Curfman McInnes { 389ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 390ae80bb75SLois Curfman McInnes int i, j; 391ae80bb75SLois Curfman McInnes Scalar *vpt = v; 392ae80bb75SLois Curfman McInnes 393ae80bb75SLois Curfman McInnes /* row-oriented output */ 394ae80bb75SLois Curfman McInnes for ( i=0; i<m; i++ ) { 395ae80bb75SLois Curfman McInnes for ( j=0; j<n; j++ ) { 396ae80bb75SLois Curfman McInnes *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]]; 397ae80bb75SLois Curfman McInnes } 398ae80bb75SLois Curfman McInnes } 399ae80bb75SLois Curfman McInnes return 0; 400ae80bb75SLois Curfman McInnes } 401ae80bb75SLois Curfman McInnes 402289bc588SBarry Smith /* -----------------------------------------------------------------*/ 403289bc588SBarry Smith 40477c4ece6SBarry Smith #include "sys.h" 40588e32edaSLois Curfman McInnes 40619bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A) 40788e32edaSLois Curfman McInnes { 40888e32edaSLois Curfman McInnes Mat_SeqDense *a; 40988e32edaSLois Curfman McInnes Mat B; 41088e32edaSLois Curfman McInnes int *scols, i, j, nz, ierr, fd, header[4], size; 41188e32edaSLois Curfman McInnes int *rowlengths = 0, M, N, *cols; 41288e32edaSLois Curfman McInnes Scalar *vals, *svals, *v; 41319bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 41488e32edaSLois Curfman McInnes 41588e32edaSLois Curfman McInnes MPI_Comm_size(comm,&size); 41688e32edaSLois Curfman McInnes if (size > 1) SETERRQ(1,"MatLoad_SeqDense: view must have one processor"); 41790ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 41877c4ece6SBarry Smith ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 41988e32edaSLois Curfman McInnes if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqDense:Not matrix object"); 42088e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 42188e32edaSLois Curfman McInnes 42290ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 42390ace30eSBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 42490ace30eSBarry Smith B = *A; 42590ace30eSBarry Smith a = (Mat_SeqDense *) B->data; 42690ace30eSBarry Smith 42790ace30eSBarry Smith /* read in nonzero values */ 42877c4ece6SBarry Smith ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr); 42990ace30eSBarry Smith 4306d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 4316d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 43290ace30eSBarry Smith ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr); 43390ace30eSBarry Smith } else { 43488e32edaSLois Curfman McInnes /* read row lengths */ 4350452661fSBarry Smith rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 43677c4ece6SBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 43788e32edaSLois Curfman McInnes 43888e32edaSLois Curfman McInnes /* create our matrix */ 439b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 44088e32edaSLois Curfman McInnes B = *A; 44188e32edaSLois Curfman McInnes a = (Mat_SeqDense *) B->data; 44288e32edaSLois Curfman McInnes v = a->v; 44388e32edaSLois Curfman McInnes 44488e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 4450452661fSBarry Smith cols = scols = (int *) PetscMalloc( nz*sizeof(int) ); CHKPTRQ(cols); 44677c4ece6SBarry Smith ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr); 4470452661fSBarry Smith vals = svals = (Scalar *) PetscMalloc( nz*sizeof(Scalar) ); CHKPTRQ(vals); 44877c4ece6SBarry Smith ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 44988e32edaSLois Curfman McInnes 45088e32edaSLois Curfman McInnes /* insert into matrix */ 45188e32edaSLois Curfman McInnes for ( i=0; i<M; i++ ) { 45288e32edaSLois Curfman McInnes for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j]; 45388e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 45488e32edaSLois Curfman McInnes } 4550452661fSBarry Smith PetscFree(vals); PetscFree(cols); PetscFree(rowlengths); 45688e32edaSLois Curfman McInnes 4576d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 4586d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 45990ace30eSBarry Smith } 46088e32edaSLois Curfman McInnes return 0; 46188e32edaSLois Curfman McInnes } 46288e32edaSLois Curfman McInnes 463932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h" 46477c4ece6SBarry Smith #include "sys.h" 465932b0c3eSLois Curfman McInnes 466932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 467289bc588SBarry Smith { 468932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 469932b0c3eSLois Curfman McInnes int ierr, i, j, format; 470d636dbe3SBarry Smith FILE *fd; 471932b0c3eSLois Curfman McInnes char *outputname; 472932b0c3eSLois Curfman McInnes Scalar *v; 473932b0c3eSLois Curfman McInnes 47490ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 475932b0c3eSLois Curfman McInnes ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 47690ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 477*639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 4787ddc982cSLois Curfman McInnes return 0; /* do nothing for now */ 479f72585f2SLois Curfman McInnes } 48002cad45dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_COMMON) { 48180cd9d93SLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 48244cd7ae7SLois Curfman McInnes v = a->v + i; 48380cd9d93SLois Curfman McInnes fprintf(fd,"row %d:",i); 48480cd9d93SLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 48580cd9d93SLois Curfman McInnes #if defined(PETSC_COMPLEX) 48680cd9d93SLois Curfman McInnes if (real(*v) != 0.0 && imag(*v) != 0.0) 48780cd9d93SLois Curfman McInnes fprintf(fd," %d %g + %g i",j,real(*v),imag(*v)); 48880cd9d93SLois Curfman McInnes else if (real(*v)) fprintf(fd," %d %g ",j,real(*v)); 48980cd9d93SLois Curfman McInnes v += a->m; 49080cd9d93SLois Curfman McInnes #else 49180cd9d93SLois Curfman McInnes if (*v) fprintf(fd," %d %g ",j,*v); 49280cd9d93SLois Curfman McInnes v += a->m; 49380cd9d93SLois Curfman McInnes #endif 49480cd9d93SLois Curfman McInnes } 49580cd9d93SLois Curfman McInnes fprintf(fd,"\n"); 49680cd9d93SLois Curfman McInnes } 49780cd9d93SLois Curfman McInnes } 498f72585f2SLois Curfman McInnes else { 49947989497SBarry Smith #if defined(PETSC_COMPLEX) 50047989497SBarry Smith int allreal = 1; 50147989497SBarry Smith /* determine if matrix has all real values */ 50247989497SBarry Smith v = a->v; 50347989497SBarry Smith for ( i=0; i<a->m*a->n; i++ ) { 50447989497SBarry Smith if (imag(v[i])) { allreal = 0; break ;} 50547989497SBarry Smith } 50647989497SBarry Smith #endif 507932b0c3eSLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 508932b0c3eSLois Curfman McInnes v = a->v + i; 509932b0c3eSLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 510289bc588SBarry Smith #if defined(PETSC_COMPLEX) 51147989497SBarry Smith if (allreal) { 51247989497SBarry Smith fprintf(fd,"%6.4e ",real(*v)); v += a->m; 51347989497SBarry Smith } else { 514932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m; 51547989497SBarry Smith } 516289bc588SBarry Smith #else 517932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e ",*v); v += a->m; 518289bc588SBarry Smith #endif 519289bc588SBarry Smith } 5208e37dc09SBarry Smith fprintf(fd,"\n"); 521289bc588SBarry Smith } 522da3a660dSBarry Smith } 523932b0c3eSLois Curfman McInnes fflush(fd); 524289bc588SBarry Smith return 0; 525289bc588SBarry Smith } 526289bc588SBarry Smith 527932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 528932b0c3eSLois Curfman McInnes { 529932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 530932b0c3eSLois Curfman McInnes int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 53190ace30eSBarry Smith int format; 53290ace30eSBarry Smith Scalar *v, *anonz,*vals; 533932b0c3eSLois Curfman McInnes 53490ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 53590ace30eSBarry Smith 53690ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 53702cad45dSBarry Smith if (format == VIEWER_FORMAT_BINARY_NATIVE) { 53890ace30eSBarry Smith /* store the matrix as a dense matrix */ 53990ace30eSBarry Smith col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens); 54090ace30eSBarry Smith col_lens[0] = MAT_COOKIE; 54190ace30eSBarry Smith col_lens[1] = m; 54290ace30eSBarry Smith col_lens[2] = n; 54390ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 54477c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr); 54590ace30eSBarry Smith PetscFree(col_lens); 54690ace30eSBarry Smith 54790ace30eSBarry Smith /* write out matrix, by rows */ 54890ace30eSBarry Smith vals = (Scalar *) PetscMalloc(m*n*sizeof(Scalar)); CHKPTRQ(vals); 54990ace30eSBarry Smith v = a->v; 55090ace30eSBarry Smith for ( i=0; i<m; i++ ) { 55190ace30eSBarry Smith for ( j=0; j<n; j++ ) { 55290ace30eSBarry Smith vals[i + j*m] = *v++; 55390ace30eSBarry Smith } 55490ace30eSBarry Smith } 55577c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr); 55690ace30eSBarry Smith PetscFree(vals); 55790ace30eSBarry Smith } else { 5580452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens); 559932b0c3eSLois Curfman McInnes col_lens[0] = MAT_COOKIE; 560932b0c3eSLois Curfman McInnes col_lens[1] = m; 561932b0c3eSLois Curfman McInnes col_lens[2] = n; 562932b0c3eSLois Curfman McInnes col_lens[3] = nz; 563932b0c3eSLois Curfman McInnes 564932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 565932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) col_lens[4+i] = n; 56677c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr); 567932b0c3eSLois Curfman McInnes 568932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 569932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 570932b0c3eSLois Curfman McInnes ict = 0; 571932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 572932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) col_lens[ict++] = j; 573932b0c3eSLois Curfman McInnes } 57477c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr); 5750452661fSBarry Smith PetscFree(col_lens); 576932b0c3eSLois Curfman McInnes 577932b0c3eSLois Curfman McInnes /* store nonzero values */ 5780452661fSBarry Smith anonz = (Scalar *) PetscMalloc((nz)*sizeof(Scalar)); CHKPTRQ(anonz); 579932b0c3eSLois Curfman McInnes ict = 0; 580932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 581932b0c3eSLois Curfman McInnes v = a->v + i; 582932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) { 583932b0c3eSLois Curfman McInnes anonz[ict++] = *v; v += a->m; 584932b0c3eSLois Curfman McInnes } 585932b0c3eSLois Curfman McInnes } 58677c4ece6SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr); 5870452661fSBarry Smith PetscFree(anonz); 58890ace30eSBarry Smith } 589932b0c3eSLois Curfman McInnes return 0; 590932b0c3eSLois Curfman McInnes } 591932b0c3eSLois Curfman McInnes 592932b0c3eSLois Curfman McInnes static int MatView_SeqDense(PetscObject obj,Viewer viewer) 593932b0c3eSLois Curfman McInnes { 594932b0c3eSLois Curfman McInnes Mat A = (Mat) obj; 595932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*) A->data; 596bcd2baecSBarry Smith ViewerType vtype; 597bcd2baecSBarry Smith int ierr; 598932b0c3eSLois Curfman McInnes 599bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 600bcd2baecSBarry Smith 601bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 602932b0c3eSLois Curfman McInnes return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); 603932b0c3eSLois Curfman McInnes } 60419bcc07fSBarry Smith else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 605932b0c3eSLois Curfman McInnes return MatView_SeqDense_ASCII(A,viewer); 606932b0c3eSLois Curfman McInnes } 607bcd2baecSBarry Smith else if (vtype == BINARY_FILE_VIEWER) { 608932b0c3eSLois Curfman McInnes return MatView_SeqDense_Binary(A,viewer); 609932b0c3eSLois Curfman McInnes } 610932b0c3eSLois Curfman McInnes return 0; 611932b0c3eSLois Curfman McInnes } 612289bc588SBarry Smith 613ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj) 614289bc588SBarry Smith { 615289bc588SBarry Smith Mat mat = (Mat) obj; 616ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 617a5a9c739SBarry Smith #if defined(PETSC_LOG) 618a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 619a5a9c739SBarry Smith #endif 6200452661fSBarry Smith if (l->pivots) PetscFree(l->pivots); 6213439631bSBarry Smith if (!l->user_alloc) PetscFree(l->v); 6220452661fSBarry Smith PetscFree(l); 623a5a9c739SBarry Smith PLogObjectDestroy(mat); 6240452661fSBarry Smith PetscHeaderDestroy(mat); 625289bc588SBarry Smith return 0; 626289bc588SBarry Smith } 627289bc588SBarry Smith 628c0bbcb79SLois Curfman McInnes static int MatTranspose_SeqDense(Mat A,Mat *matout) 629289bc588SBarry Smith { 630c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 631d3e5ee88SLois Curfman McInnes int k, j, m, n; 632d3e5ee88SLois Curfman McInnes Scalar *v, tmp; 63348b35521SBarry Smith 634d3e5ee88SLois Curfman McInnes v = mat->v; m = mat->m; n = mat->n; 6353638b69dSLois Curfman McInnes if (matout == PETSC_NULL) { /* in place transpose */ 636d9f96c7cSLois Curfman McInnes if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Supports square matrix only in-place"); 637d3e5ee88SLois Curfman McInnes for ( j=0; j<m; j++ ) { 638289bc588SBarry Smith for ( k=0; k<j; k++ ) { 639d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 640d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 641d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 642289bc588SBarry Smith } 643289bc588SBarry Smith } 64448b35521SBarry Smith } 645d3e5ee88SLois Curfman McInnes else { /* out-of-place transpose */ 646d3e5ee88SLois Curfman McInnes int ierr; 647d3e5ee88SLois Curfman McInnes Mat tmat; 648ec8511deSBarry Smith Mat_SeqDense *tmatd; 649d3e5ee88SLois Curfman McInnes Scalar *v2; 650b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr); 651ec8511deSBarry Smith tmatd = (Mat_SeqDense *) tmat->data; 6520de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 653d3e5ee88SLois Curfman McInnes for ( j=0; j<n; j++ ) { 6540de55854SLois Curfman McInnes for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 655d3e5ee88SLois Curfman McInnes } 6566d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 6576d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 658d3e5ee88SLois Curfman McInnes *matout = tmat; 65948b35521SBarry Smith } 660289bc588SBarry Smith return 0; 661289bc588SBarry Smith } 662289bc588SBarry Smith 66377c4ece6SBarry Smith static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg) 664289bc588SBarry Smith { 665c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 666c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 667289bc588SBarry Smith int i; 668289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 6699ea5d5aeSSatish Balay 6709ea5d5aeSSatish Balay if (A2->type != MATSEQDENSE) SETERRQ(1,"MatEqual_SeqDense:Both matrices should be of type MATSEQDENSE"); 67177c4ece6SBarry Smith if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;} 67277c4ece6SBarry Smith if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;} 673289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 67477c4ece6SBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;} 675289bc588SBarry Smith v1++; v2++; 676289bc588SBarry Smith } 67777c4ece6SBarry Smith *flg = PETSC_TRUE; 6789ea5d5aeSSatish Balay return 0; 679289bc588SBarry Smith } 680289bc588SBarry Smith 681c0bbcb79SLois Curfman McInnes static int MatGetDiagonal_SeqDense(Mat A,Vec v) 682289bc588SBarry Smith { 683c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 68444cd7ae7SLois Curfman McInnes int i, n, len; 68544cd7ae7SLois Curfman McInnes Scalar *x, zero = 0.0; 68644cd7ae7SLois Curfman McInnes 68744cd7ae7SLois Curfman McInnes VecSet(&zero,v); 688289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 68944cd7ae7SLois Curfman McInnes len = PetscMin(mat->m,mat->n); 690ec8511deSBarry Smith if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 69144cd7ae7SLois Curfman McInnes for ( i=0; i<len; i++ ) { 692289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 693289bc588SBarry Smith } 694289bc588SBarry Smith return 0; 695289bc588SBarry Smith } 696289bc588SBarry Smith 697052efed2SBarry Smith static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 698289bc588SBarry Smith { 699c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 700da3a660dSBarry Smith Scalar *l,*r,x,*v; 701da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 70255659b69SBarry Smith 70328988994SBarry Smith if (ll) { 704da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 705052efed2SBarry Smith if (m != mat->m) SETERRQ(1,"MatDiagonalScale_SeqDense:Left scaling vec wrong size"); 706da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 707da3a660dSBarry Smith x = l[i]; 708da3a660dSBarry Smith v = mat->v + i; 709da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 710da3a660dSBarry Smith } 71144cd7ae7SLois Curfman McInnes PLogFlops(n*m); 712da3a660dSBarry Smith } 71328988994SBarry Smith if (rr) { 714da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 715052efed2SBarry Smith if (n != mat->n) SETERRQ(1,"MatDiagonalScale_SeqDense:Right scaling vec wrong size"); 716da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 717da3a660dSBarry Smith x = r[i]; 718da3a660dSBarry Smith v = mat->v + i*m; 719da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 720da3a660dSBarry Smith } 72144cd7ae7SLois Curfman McInnes PLogFlops(n*m); 722da3a660dSBarry Smith } 723289bc588SBarry Smith return 0; 724289bc588SBarry Smith } 725289bc588SBarry Smith 726cddf8d76SBarry Smith static int MatNorm_SeqDense(Mat A,NormType type,double *norm) 727289bc588SBarry Smith { 728c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 729289bc588SBarry Smith Scalar *v = mat->v; 730289bc588SBarry Smith double sum = 0.0; 731289bc588SBarry Smith int i, j; 73255659b69SBarry Smith 733289bc588SBarry Smith if (type == NORM_FROBENIUS) { 734289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 735289bc588SBarry Smith #if defined(PETSC_COMPLEX) 736289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 737289bc588SBarry Smith #else 738289bc588SBarry Smith sum += (*v)*(*v); v++; 739289bc588SBarry Smith #endif 740289bc588SBarry Smith } 741289bc588SBarry Smith *norm = sqrt(sum); 74255659b69SBarry Smith PLogFlops(2*mat->n*mat->m); 743289bc588SBarry Smith } 744289bc588SBarry Smith else if (type == NORM_1) { 745289bc588SBarry Smith *norm = 0.0; 746289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 747289bc588SBarry Smith sum = 0.0; 748289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 74933a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 750289bc588SBarry Smith } 751289bc588SBarry Smith if (sum > *norm) *norm = sum; 752289bc588SBarry Smith } 75355659b69SBarry Smith PLogFlops(mat->n*mat->m); 754289bc588SBarry Smith } 755289bc588SBarry Smith else if (type == NORM_INFINITY) { 756289bc588SBarry Smith *norm = 0.0; 757289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 758289bc588SBarry Smith v = mat->v + j; 759289bc588SBarry Smith sum = 0.0; 760289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 761cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v += mat->m; 762289bc588SBarry Smith } 763289bc588SBarry Smith if (sum > *norm) *norm = sum; 764289bc588SBarry Smith } 76555659b69SBarry Smith PLogFlops(mat->n*mat->m); 766289bc588SBarry Smith } 767289bc588SBarry Smith else { 76848d91487SBarry Smith SETERRQ(1,"MatNorm_SeqDense:No two norm"); 769289bc588SBarry Smith } 770289bc588SBarry Smith return 0; 771289bc588SBarry Smith } 772289bc588SBarry Smith 773c0bbcb79SLois Curfman McInnes static int MatSetOption_SeqDense(Mat A,MatOption op) 774289bc588SBarry Smith { 775c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 77667e560aaSBarry Smith 7776d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) aij->roworiented = 1; 7786d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) aij->roworiented = 0; 7796d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 7806d4a8577SBarry Smith op == MAT_COLUMNS_SORTED || 7816d4a8577SBarry Smith op == MAT_SYMMETRIC || 7826d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 7836d4a8577SBarry Smith op == MAT_NO_NEW_NONZERO_LOCATIONS || 7846d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 7856d4a8577SBarry Smith op == MAT_NO_NEW_DIAGONALS || 7866d4a8577SBarry Smith op == MAT_YES_NEW_DIAGONALS) 78794a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n"); 788c0bbcb79SLois Curfman McInnes else 7894d39ef84SLois Curfman McInnes {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");} 790289bc588SBarry Smith return 0; 791289bc588SBarry Smith } 792289bc588SBarry Smith 793ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A) 7946f0a148fSBarry Smith { 795ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 796cddf8d76SBarry Smith PetscMemzero(l->v,l->m*l->n*sizeof(Scalar)); 7976f0a148fSBarry Smith return 0; 7986f0a148fSBarry Smith } 7996f0a148fSBarry Smith 8004e220ebcSLois Curfman McInnes static int MatGetBlockSize_SeqDense(Mat A,int *bs) 8014e220ebcSLois Curfman McInnes { 8024e220ebcSLois Curfman McInnes *bs = 1; 8034e220ebcSLois Curfman McInnes return 0; 8044e220ebcSLois Curfman McInnes } 8054e220ebcSLois Curfman McInnes 806ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 8076f0a148fSBarry Smith { 808ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 8096abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 8106f0a148fSBarry Smith Scalar *slot; 81155659b69SBarry Smith 81277c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 81378b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 8146f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 8156f0a148fSBarry Smith slot = l->v + rows[i]; 8166f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 8176f0a148fSBarry Smith } 8186f0a148fSBarry Smith if (diag) { 8196f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 8206f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 821260925b8SBarry Smith *slot = *diag; 8226f0a148fSBarry Smith } 8236f0a148fSBarry Smith } 824260925b8SBarry Smith ISRestoreIndices(is,&rows); 8256f0a148fSBarry Smith return 0; 8266f0a148fSBarry Smith } 827557bce09SLois Curfman McInnes 828c0bbcb79SLois Curfman McInnes static int MatGetSize_SeqDense(Mat A,int *m,int *n) 829557bce09SLois Curfman McInnes { 830c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 831557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 832557bce09SLois Curfman McInnes return 0; 833557bce09SLois Curfman McInnes } 834557bce09SLois Curfman McInnes 835c0bbcb79SLois Curfman McInnes static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 836d3e5ee88SLois Curfman McInnes { 837c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 838d3e5ee88SLois Curfman McInnes *m = 0; *n = mat->m; 839d3e5ee88SLois Curfman McInnes return 0; 840d3e5ee88SLois Curfman McInnes } 841d3e5ee88SLois Curfman McInnes 842c0bbcb79SLois Curfman McInnes static int MatGetArray_SeqDense(Mat A,Scalar **array) 84364e87e97SBarry Smith { 844c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 84564e87e97SBarry Smith *array = mat->v; 84664e87e97SBarry Smith return 0; 84764e87e97SBarry Smith } 8480754003eSLois Curfman McInnes 849ff14e315SSatish Balay static int MatRestoreArray_SeqDense(Mat A,Scalar **array) 850ff14e315SSatish Balay { 851ff14e315SSatish Balay return 0; 852ff14e315SSatish Balay } 8530754003eSLois Curfman McInnes 854cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall, 855cddf8d76SBarry Smith Mat *submat) 8560754003eSLois Curfman McInnes { 857c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 8580754003eSLois Curfman McInnes int nznew, *smap, i, j, ierr, oldcols = mat->n; 859160018dcSBarry Smith int *irow, *icol, nrows, ncols, *cwork; 8600754003eSLois Curfman McInnes Scalar *vwork, *val; 8610754003eSLois Curfman McInnes Mat newmat; 8620754003eSLois Curfman McInnes 86378b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 86478b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 86578b31e54SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 86678b31e54SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 8670754003eSLois Curfman McInnes 8680452661fSBarry Smith smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap); 8690452661fSBarry Smith cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork); 8700452661fSBarry Smith vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 871cddf8d76SBarry Smith PetscMemzero((char*)smap,oldcols*sizeof(int)); 8720754003eSLois Curfman McInnes for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 8730754003eSLois Curfman McInnes 8740754003eSLois Curfman McInnes /* Create and fill new matrix */ 875b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr); 8760754003eSLois Curfman McInnes for (i=0; i<nrows; i++) { 8770754003eSLois Curfman McInnes nznew = 0; 8780754003eSLois Curfman McInnes val = mat->v + irow[i]; 8790754003eSLois Curfman McInnes for (j=0; j<oldcols; j++) { 8800754003eSLois Curfman McInnes if (smap[j]) { 8810754003eSLois Curfman McInnes cwork[nznew] = smap[j] - 1; 8820754003eSLois Curfman McInnes vwork[nznew++] = val[j * mat->m]; 8830754003eSLois Curfman McInnes } 8840754003eSLois Curfman McInnes } 885dbb450caSBarry Smith ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES); 88678b31e54SBarry Smith CHKERRQ(ierr); 8870754003eSLois Curfman McInnes } 8886d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 8896d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 8900754003eSLois Curfman McInnes 8910754003eSLois Curfman McInnes /* Free work space */ 8920452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 89378b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 89478b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 8950754003eSLois Curfman McInnes *submat = newmat; 8960754003eSLois Curfman McInnes return 0; 8970754003eSLois Curfman McInnes } 8980754003eSLois Curfman McInnes 899905e6a2fSBarry Smith static int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 900905e6a2fSBarry Smith Mat **B) 901905e6a2fSBarry Smith { 902905e6a2fSBarry Smith int ierr,i; 903905e6a2fSBarry Smith 904905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 905905e6a2fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 906905e6a2fSBarry Smith } 907905e6a2fSBarry Smith 908905e6a2fSBarry Smith for ( i=0; i<n; i++ ) { 909905e6a2fSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 910905e6a2fSBarry Smith } 911905e6a2fSBarry Smith return 0; 912905e6a2fSBarry Smith } 913905e6a2fSBarry Smith 9144b0e389bSBarry Smith static int MatCopy_SeqDense(Mat A, Mat B) 9154b0e389bSBarry Smith { 9164b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data; 9174b0e389bSBarry Smith if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B); 9184b0e389bSBarry Smith if (a->m != b->m || a->n != b->n) SETERRQ(1,"MatCopy_SeqDense:size(B) != size(A)"); 9194b0e389bSBarry Smith PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar)); 9204b0e389bSBarry Smith return 0; 9214b0e389bSBarry Smith } 9224b0e389bSBarry Smith 923289bc588SBarry Smith /* -------------------------------------------------------------------*/ 924ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense, 925905e6a2fSBarry Smith MatGetRow_SeqDense, 926905e6a2fSBarry Smith MatRestoreRow_SeqDense, 927905e6a2fSBarry Smith MatMult_SeqDense, 928905e6a2fSBarry Smith MatMultAdd_SeqDense, 929905e6a2fSBarry Smith MatMultTrans_SeqDense, 930905e6a2fSBarry Smith MatMultTransAdd_SeqDense, 931905e6a2fSBarry Smith MatSolve_SeqDense, 932905e6a2fSBarry Smith MatSolveAdd_SeqDense, 933905e6a2fSBarry Smith MatSolveTrans_SeqDense, 934905e6a2fSBarry Smith MatSolveTransAdd_SeqDense, 935905e6a2fSBarry Smith MatLUFactor_SeqDense, 936905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 937ec8511deSBarry Smith MatRelax_SeqDense, 938ec8511deSBarry Smith MatTranspose_SeqDense, 939905e6a2fSBarry Smith MatGetInfo_SeqDense, 940905e6a2fSBarry Smith MatEqual_SeqDense, 941905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 942905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 943905e6a2fSBarry Smith MatNorm_SeqDense, 944905e6a2fSBarry Smith 0, 945905e6a2fSBarry Smith 0, 946905e6a2fSBarry Smith 0, 947905e6a2fSBarry Smith MatSetOption_SeqDense, 948905e6a2fSBarry Smith MatZeroEntries_SeqDense, 949905e6a2fSBarry Smith MatZeroRows_SeqDense, 950905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 951905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 952905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 953905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 954905e6a2fSBarry Smith MatGetSize_SeqDense, 955905e6a2fSBarry Smith MatGetSize_SeqDense, 956905e6a2fSBarry Smith MatGetOwnershipRange_SeqDense, 957905e6a2fSBarry Smith 0, 958905e6a2fSBarry Smith 0, 959905e6a2fSBarry Smith MatGetArray_SeqDense, 960905e6a2fSBarry Smith MatRestoreArray_SeqDense, 961905e6a2fSBarry Smith 0, 9623d1612f7SBarry Smith MatConvertSameType_SeqDense,0,0,0,0, 963905e6a2fSBarry Smith MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0, 9644b0e389bSBarry Smith MatGetValues_SeqDense, 9654e220ebcSLois Curfman McInnes MatCopy_SeqDense,0,MatScale_SeqDense, 9664e220ebcSLois Curfman McInnes 0,0,0,MatGetBlockSize_SeqDense}; 96790ace30eSBarry Smith 9684b828684SBarry Smith /*@C 969fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 970d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 971d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 972289bc588SBarry Smith 97320563c6bSBarry Smith Input Parameters: 9740c775827SLois Curfman McInnes . comm - MPI communicator, set to MPI_COMM_SELF 9750c775827SLois Curfman McInnes . m - number of rows 97618f449edSLois Curfman McInnes . n - number of columns 977b4fd4287SBarry Smith . data - optional location of matrix data. Set data=PETSC_NULL for PETSc 978dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 97920563c6bSBarry Smith 98020563c6bSBarry Smith Output Parameter: 98144cd7ae7SLois Curfman McInnes . A - the matrix 98220563c6bSBarry Smith 98318f449edSLois Curfman McInnes Notes: 98418f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 98518f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 986b4fd4287SBarry Smith set data=PETSC_NULL. 98718f449edSLois Curfman McInnes 988dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 989d65003e9SLois Curfman McInnes 990dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues() 99120563c6bSBarry Smith @*/ 99244cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A) 993289bc588SBarry Smith { 99444cd7ae7SLois Curfman McInnes Mat B; 99544cd7ae7SLois Curfman McInnes Mat_SeqDense *b; 9963b2fbd54SBarry Smith int ierr,flg,size; 9973b2fbd54SBarry Smith 9983b2fbd54SBarry Smith MPI_Comm_size(comm,&size); 9993b2fbd54SBarry Smith if (size > 1) SETERRQ(1,"MatCreateSeqDense:Comm must be of size 1"); 100055659b69SBarry Smith 100144cd7ae7SLois Curfman McInnes *A = 0; 100244cd7ae7SLois Curfman McInnes PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQDENSE,comm); 100344cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 100444cd7ae7SLois Curfman McInnes b = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b); 100544cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqDense)); 100644cd7ae7SLois Curfman McInnes PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 100744cd7ae7SLois Curfman McInnes B->destroy = MatDestroy_SeqDense; 100844cd7ae7SLois Curfman McInnes B->view = MatView_SeqDense; 100944cd7ae7SLois Curfman McInnes B->factor = 0; 101044cd7ae7SLois Curfman McInnes PLogObjectMemory(B,sizeof(struct _Mat)); 101144cd7ae7SLois Curfman McInnes B->data = (void *) b; 101218f449edSLois Curfman McInnes 101344cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 101444cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 101544cd7ae7SLois Curfman McInnes b->pivots = 0; 101644cd7ae7SLois Curfman McInnes b->roworiented = 1; 1017b4fd4287SBarry Smith if (data == PETSC_NULL) { 101844cd7ae7SLois Curfman McInnes b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v); 101944cd7ae7SLois Curfman McInnes PetscMemzero(b->v,m*n*sizeof(Scalar)); 102044cd7ae7SLois Curfman McInnes b->user_alloc = 0; 102144cd7ae7SLois Curfman McInnes PLogObjectMemory(B,n*m*sizeof(Scalar)); 102218f449edSLois Curfman McInnes } 10232dd2b441SLois Curfman McInnes else { /* user-allocated storage */ 102444cd7ae7SLois Curfman McInnes b->v = data; 102544cd7ae7SLois Curfman McInnes b->user_alloc = 1; 10262dd2b441SLois Curfman McInnes } 102725cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 102844cd7ae7SLois Curfman McInnes if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);} 10294e220ebcSLois Curfman McInnes 103044cd7ae7SLois Curfman McInnes *A = B; 1031289bc588SBarry Smith return 0; 1032289bc588SBarry Smith } 1033289bc588SBarry Smith 1034c0bbcb79SLois Curfman McInnes int MatCreate_SeqDense(Mat A,Mat *newmat) 1035289bc588SBarry Smith { 1036c0bbcb79SLois Curfman McInnes Mat_SeqDense *m = (Mat_SeqDense *) A->data; 1037b4fd4287SBarry Smith return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat); 1038289bc588SBarry Smith } 10393b2fbd54SBarry Smith 10403b2fbd54SBarry Smith 10413b2fbd54SBarry Smith 10423b2fbd54SBarry Smith 10433b2fbd54SBarry Smith 1044