1cb512458SBarry Smith #ifndef lint 2*48d91487SBarry Smith static char vcid[] = "$Id: dense.c,v 1.60 1995/09/30 19:28:39 bsmith Exp bsmith $"; 3cb512458SBarry Smith #endif 4289bc588SBarry Smith 5289bc588SBarry Smith /* 6289bc588SBarry Smith Standard Fortran style matrices 7289bc588SBarry Smith */ 819b02663SBarry Smith #include "petsc.h" 9b16a3bb1SBarry Smith #include "pinclude/plapack.h" 10289bc588SBarry Smith #include "matimpl.h" 11289bc588SBarry Smith #include "math.h" 12289bc588SBarry Smith #include "vec/vecimpl.h" 13b16a3bb1SBarry Smith #include "pinclude/pviewer.h" 14289bc588SBarry Smith 15289bc588SBarry Smith typedef struct { 16289bc588SBarry Smith Scalar *v; 17289bc588SBarry Smith int roworiented; 18289bc588SBarry Smith int m,n,pad; 19289bc588SBarry Smith int *pivots; /* pivots in LU factorization */ 20ec8511deSBarry Smith } Mat_SeqDense; 21289bc588SBarry Smith 22*48d91487SBarry Smith static int MatGetInfo_SeqDense(Mat matin,MatInfoType flag,int *nz,int *nzalloc,int *mem) 23289bc588SBarry Smith { 24ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 25289bc588SBarry Smith int i,N = mat->m*mat->n,count = 0; 26289bc588SBarry Smith Scalar *v = mat->v; 27289bc588SBarry Smith for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 28752f5567SLois Curfman McInnes *nz = count; *nzalloc = N; *mem = (int)matin->mem; 29fa9ec3f1SLois Curfman McInnes return 0; 30289bc588SBarry Smith } 31289bc588SBarry Smith 32289bc588SBarry Smith /* ---------------------------------------------------------------*/ 33289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots, 34289bc588SBarry Smith rather than put it in the Mat->col slot.*/ 35ec8511deSBarry Smith static int MatLUFactor_SeqDense(Mat matin,IS row,IS col,double f) 36289bc588SBarry Smith { 37ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 386abc6512SBarry Smith int info; 39289bc588SBarry Smith if (!mat->pivots) { 40*48d91487SBarry Smith mat->pivots = (int *) PETSCMALLOC(mat->m*sizeof(int));CHKPTRQ(mat->pivots); 41752f5567SLois Curfman McInnes PLogObjectMemory(matin,mat->m*sizeof(int)); 42289bc588SBarry Smith } 43289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 44ec8511deSBarry Smith if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization"); 45289bc588SBarry Smith matin->factor = FACTOR_LU; 46289bc588SBarry Smith return 0; 47289bc588SBarry Smith } 48*48d91487SBarry Smith static int MatLUFactorSymbolic_SeqDense(Mat matin,IS row,IS col,double f,Mat *fact) 49289bc588SBarry Smith { 50289bc588SBarry Smith int ierr; 51bbb6d6a8SBarry Smith ierr = MatConvert(matin,MATSAME,fact); CHKERRQ(ierr); 52289bc588SBarry Smith return 0; 53289bc588SBarry Smith } 54ec8511deSBarry Smith static int MatLUFactorNumeric_SeqDense(Mat matin,Mat *fact) 55289bc588SBarry Smith { 5649d8b64dSBarry Smith return MatLUFactor(*fact,0,0,1.0); 57289bc588SBarry Smith } 58ec8511deSBarry Smith static int MatCholeskyFactorSymbolic_SeqDense(Mat matin,IS row,double f,Mat *fact) 59289bc588SBarry Smith { 60289bc588SBarry Smith int ierr; 61bbb6d6a8SBarry Smith ierr = MatConvert(matin,MATSAME,fact); CHKERRQ(ierr); 62289bc588SBarry Smith return 0; 63289bc588SBarry Smith } 64ec8511deSBarry Smith static int MatCholeskyFactorNumeric_SeqDense(Mat matin,Mat *fact) 65289bc588SBarry Smith { 6649d8b64dSBarry Smith return MatCholeskyFactor(*fact,0,1.0); 67289bc588SBarry Smith } 68ec8511deSBarry Smith static int MatCholeskyFactor_SeqDense(Mat matin,IS perm,double f) 69289bc588SBarry Smith { 70ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 716abc6512SBarry Smith int info; 72752f5567SLois Curfman McInnes if (mat->pivots) { 73752f5567SLois Curfman McInnes PETSCFREE(mat->pivots); 74752f5567SLois Curfman McInnes PLogObjectMemory(matin,-mat->m*sizeof(int)); 75752f5567SLois Curfman McInnes mat->pivots = 0; 76752f5567SLois Curfman McInnes } 77289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 78ec8511deSBarry Smith if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization"); 79289bc588SBarry Smith matin->factor = FACTOR_CHOLESKY; 80289bc588SBarry Smith return 0; 81289bc588SBarry Smith } 82289bc588SBarry Smith 83ec8511deSBarry Smith static int MatSolve_SeqDense(Mat matin,Vec xx,Vec yy) 84289bc588SBarry Smith { 85ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 866abc6512SBarry Smith int one = 1, info; 876abc6512SBarry Smith Scalar *x, *y; 88289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 89416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 909e25ed09SBarry Smith if (matin->factor == FACTOR_LU) { 91*48d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 92289bc588SBarry Smith } 939e25ed09SBarry Smith else if (matin->factor == FACTOR_CHOLESKY){ 94*48d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 95289bc588SBarry Smith } 96ec8511deSBarry Smith else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve"); 97ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve"); 98289bc588SBarry Smith return 0; 99289bc588SBarry Smith } 100ec8511deSBarry Smith static int MatSolveTrans_SeqDense(Mat matin,Vec xx,Vec yy) 101da3a660dSBarry Smith { 102ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 1036abc6512SBarry Smith int one = 1, info; 1046abc6512SBarry Smith Scalar *x, *y; 105da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 106416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 107752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 108da3a660dSBarry Smith if (mat->pivots) { 109*48d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 110da3a660dSBarry Smith } 111da3a660dSBarry Smith else { 112*48d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 113da3a660dSBarry Smith } 114ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve"); 115da3a660dSBarry Smith return 0; 116da3a660dSBarry Smith } 117ec8511deSBarry Smith static int MatSolveAdd_SeqDense(Mat matin,Vec xx,Vec zz,Vec yy) 118da3a660dSBarry Smith { 119ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 1206abc6512SBarry Smith int one = 1, info,ierr; 1216abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 122da3a660dSBarry Smith Vec tmp = 0; 123da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 124da3a660dSBarry Smith if (yy == zz) { 12578b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 126464493b3SBarry Smith PLogObjectParent(matin,tmp); 12778b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 128da3a660dSBarry Smith } 129416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 130752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 131da3a660dSBarry Smith if (mat->pivots) { 132*48d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 133da3a660dSBarry Smith } 134da3a660dSBarry Smith else { 135*48d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 136da3a660dSBarry Smith } 137ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve"); 138da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 139da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 140da3a660dSBarry Smith return 0; 141da3a660dSBarry Smith } 142ec8511deSBarry Smith static int MatSolveTransAdd_SeqDense(Mat matin,Vec xx,Vec zz, Vec yy) 143da3a660dSBarry Smith { 144ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 1456abc6512SBarry Smith int one = 1, info,ierr; 1466abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 147da3a660dSBarry Smith Vec tmp; 148da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 149da3a660dSBarry Smith if (yy == zz) { 15078b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 151464493b3SBarry Smith PLogObjectParent(matin,tmp); 15278b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 153da3a660dSBarry Smith } 154416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 155752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 156da3a660dSBarry Smith if (mat->pivots) { 157*48d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 158da3a660dSBarry Smith } 159da3a660dSBarry Smith else { 160*48d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 161da3a660dSBarry Smith } 162ec8511deSBarry Smith if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve"); 163da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 164da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 165da3a660dSBarry Smith return 0; 166da3a660dSBarry Smith } 167289bc588SBarry Smith /* ------------------------------------------------------------------*/ 168ec8511deSBarry Smith static int MatRelax_SeqDense(Mat matin,Vec bb,double omega,MatSORType flag, 16920e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 170289bc588SBarry Smith { 171ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 172289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 1736abc6512SBarry Smith int o = 1,ierr, m = mat->m, i; 174289bc588SBarry Smith 175289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 176289bc588SBarry Smith /* this is a hack fix, should have another version without 177289bc588SBarry Smith the second BLdot */ 178bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx); CHKERRQ(ierr); 179289bc588SBarry Smith } 180289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 181289bc588SBarry Smith while (its--) { 182289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 183289bc588SBarry Smith for ( i=0; i<m; i++ ) { 184f1747703SBarry Smith #if defined(PETSC_COMPLEX) 185f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 186f1747703SBarry Smith not happy about returning a double complex */ 187f1747703SBarry Smith int _i; 188f1747703SBarry Smith Scalar sum = b[i]; 189f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 190f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 191f1747703SBarry Smith } 192f1747703SBarry Smith xt = sum; 193f1747703SBarry Smith #else 194289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 195f1747703SBarry Smith #endif 196d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 197289bc588SBarry Smith } 198289bc588SBarry Smith } 199289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 200289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 201f1747703SBarry Smith #if defined(PETSC_COMPLEX) 202f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 203f1747703SBarry Smith not happy about returning a double complex */ 204f1747703SBarry Smith int _i; 205f1747703SBarry Smith Scalar sum = b[i]; 206f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 207f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 208f1747703SBarry Smith } 209f1747703SBarry Smith xt = sum; 210f1747703SBarry Smith #else 211289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 212f1747703SBarry Smith #endif 213d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 214289bc588SBarry Smith } 215289bc588SBarry Smith } 216289bc588SBarry Smith } 217289bc588SBarry Smith return 0; 218289bc588SBarry Smith } 219289bc588SBarry Smith 220289bc588SBarry Smith /* -----------------------------------------------------------------*/ 221ec8511deSBarry Smith static int MatMultTrans_SeqDense(Mat matin,Vec xx,Vec yy) 222289bc588SBarry Smith { 223ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 224289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 225289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 226289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 227*48d91487SBarry Smith LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 228289bc588SBarry Smith return 0; 229289bc588SBarry Smith } 230ec8511deSBarry Smith static int MatMult_SeqDense(Mat matin,Vec xx,Vec yy) 231289bc588SBarry Smith { 232ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 233289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 234289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 235289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 236*48d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 237289bc588SBarry Smith return 0; 238289bc588SBarry Smith } 239ec8511deSBarry Smith static int MatMultAdd_SeqDense(Mat matin,Vec xx,Vec zz,Vec yy) 240289bc588SBarry Smith { 241ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 242289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 2436abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 244289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 245416022c9SBarry Smith if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 246*48d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 247289bc588SBarry Smith return 0; 248289bc588SBarry Smith } 249ec8511deSBarry Smith static int MatMultTransAdd_SeqDense(Mat matin,Vec xx,Vec zz,Vec yy) 250289bc588SBarry Smith { 251ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 252289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 253289bc588SBarry Smith int _One=1; 2546abc6512SBarry Smith Scalar _DOne=1.0; 255289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 256289bc588SBarry Smith VecGetArray(zz,&z); 257416022c9SBarry Smith if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 258*48d91487SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 259289bc588SBarry Smith return 0; 260289bc588SBarry Smith } 261289bc588SBarry Smith 262289bc588SBarry Smith /* -----------------------------------------------------------------*/ 263*48d91487SBarry Smith static int MatGetRow_SeqDense(Mat matin,int row,int *ncols,int **cols,Scalar **vals) 264289bc588SBarry Smith { 265ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 266289bc588SBarry Smith Scalar *v; 267289bc588SBarry Smith int i; 268289bc588SBarry Smith *ncols = mat->n; 269289bc588SBarry Smith if (cols) { 27078b31e54SBarry Smith *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols); 271289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) *cols[i] = i; 272289bc588SBarry Smith } 273289bc588SBarry Smith if (vals) { 27478b31e54SBarry Smith *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 275289bc588SBarry Smith v = mat->v + row; 276289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 277289bc588SBarry Smith } 278289bc588SBarry Smith return 0; 279289bc588SBarry Smith } 280*48d91487SBarry Smith static int MatRestoreRow_SeqDense(Mat matin,int row,int *ncols,int **cols,Scalar **vals) 281289bc588SBarry Smith { 28278b31e54SBarry Smith if (cols) { PETSCFREE(*cols); } 28378b31e54SBarry Smith if (vals) { PETSCFREE(*vals); } 284289bc588SBarry Smith return 0; 285289bc588SBarry Smith } 286289bc588SBarry Smith /* ----------------------------------------------------------------*/ 287ec8511deSBarry Smith static int MatInsert_SeqDense(Mat matin,int m,int *indexm,int n, 288e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 289289bc588SBarry Smith { 290ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 291289bc588SBarry Smith int i,j; 292d6dfbf8fSBarry Smith 293289bc588SBarry Smith if (!mat->roworiented) { 294dbb450caSBarry Smith if (addv == INSERT_VALUES) { 295289bc588SBarry Smith for ( j=0; j<n; j++ ) { 296289bc588SBarry Smith for ( i=0; i<m; i++ ) { 297289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 298289bc588SBarry Smith } 299289bc588SBarry Smith } 300289bc588SBarry Smith } 301289bc588SBarry Smith else { 302289bc588SBarry Smith for ( j=0; j<n; j++ ) { 303289bc588SBarry Smith for ( i=0; i<m; i++ ) { 304289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 305289bc588SBarry Smith } 306289bc588SBarry Smith } 307289bc588SBarry Smith } 308e8d4e0b9SBarry Smith } 309e8d4e0b9SBarry Smith else { 310dbb450caSBarry Smith if (addv == INSERT_VALUES) { 311e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 312e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 313e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 314e8d4e0b9SBarry Smith } 315e8d4e0b9SBarry Smith } 316e8d4e0b9SBarry Smith } 317289bc588SBarry Smith else { 318289bc588SBarry Smith for ( i=0; i<m; i++ ) { 319289bc588SBarry Smith for ( j=0; j<n; j++ ) { 320289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 321289bc588SBarry Smith } 322289bc588SBarry Smith } 323289bc588SBarry Smith } 324e8d4e0b9SBarry Smith } 325289bc588SBarry Smith return 0; 326289bc588SBarry Smith } 327e8d4e0b9SBarry Smith 328289bc588SBarry Smith /* -----------------------------------------------------------------*/ 329ec8511deSBarry Smith static int MatCopyPrivate_SeqDense(Mat matin,Mat *newmat) 330289bc588SBarry Smith { 331*48d91487SBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data,*l; 332289bc588SBarry Smith int ierr; 333289bc588SBarry Smith Mat newi; 334*48d91487SBarry Smith 335*48d91487SBarry Smith ierr = MatCreateSeqDense(matin->comm,mat->m,mat->n,&newi);CHKERRQ(ierr); 336ec8511deSBarry Smith l = (Mat_SeqDense *) newi->data; 337416022c9SBarry Smith PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 338289bc588SBarry Smith *newmat = newi; 339289bc588SBarry Smith return 0; 340289bc588SBarry Smith } 341da3a660dSBarry Smith #include "viewer.h" 342289bc588SBarry Smith 343ec8511deSBarry Smith int MatView_SeqDense(PetscObject obj,Viewer ptr) 344289bc588SBarry Smith { 345289bc588SBarry Smith Mat matin = (Mat) obj; 346ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 347289bc588SBarry Smith Scalar *v; 348d636dbe3SBarry Smith int i,j,ierr; 34923242f5aSBarry Smith PetscObject vobj = (PetscObject) ptr; 350da3a660dSBarry Smith 35123242f5aSBarry Smith if (ptr == 0) { 35221b0d8fbSLois Curfman McInnes ptr = STDOUT_VIEWER_SELF; vobj = (PetscObject) ptr; 35323242f5aSBarry Smith } 35423242f5aSBarry Smith if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) { 3556f75cfa4SBarry Smith return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v); 356da3a660dSBarry Smith } 357da3a660dSBarry Smith else { 358d636dbe3SBarry Smith FILE *fd; 359d636dbe3SBarry Smith int format; 360d636dbe3SBarry Smith ierr = ViewerFileGetPointer_Private(ptr,&fd); CHKERRQ(ierr); 361d636dbe3SBarry Smith ierr = ViewerFileGetFormat_Private(ptr,&format); CHKERRQ(ierr); 362f72585f2SLois Curfman McInnes if (format == FILE_FORMAT_INFO) { 363f72585f2SLois Curfman McInnes /* do nothing for now */ 364f72585f2SLois Curfman McInnes } 365f72585f2SLois Curfman McInnes else { 366289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 367289bc588SBarry Smith v = mat->v + i; 368289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 369289bc588SBarry Smith #if defined(PETSC_COMPLEX) 3708e37dc09SBarry Smith fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 371289bc588SBarry Smith #else 3728e37dc09SBarry Smith fprintf(fd,"%6.4e ",*v); v += mat->m; 373289bc588SBarry Smith #endif 374289bc588SBarry Smith } 3758e37dc09SBarry Smith fprintf(fd,"\n"); 376289bc588SBarry Smith } 377da3a660dSBarry Smith } 378f72585f2SLois Curfman McInnes } 379289bc588SBarry Smith return 0; 380289bc588SBarry Smith } 381289bc588SBarry Smith 382289bc588SBarry Smith 383ec8511deSBarry Smith static int MatDestroy_SeqDense(PetscObject obj) 384289bc588SBarry Smith { 385289bc588SBarry Smith Mat mat = (Mat) obj; 386ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 387a5a9c739SBarry Smith #if defined(PETSC_LOG) 388a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 389a5a9c739SBarry Smith #endif 39078b31e54SBarry Smith if (l->pivots) PETSCFREE(l->pivots); 39178b31e54SBarry Smith PETSCFREE(l); 392a5a9c739SBarry Smith PLogObjectDestroy(mat); 3939e25ed09SBarry Smith PETSCHEADERDESTROY(mat); 394289bc588SBarry Smith return 0; 395289bc588SBarry Smith } 396289bc588SBarry Smith 397ec8511deSBarry Smith static int MatTranspose_SeqDense(Mat matin,Mat *matout) 398289bc588SBarry Smith { 399ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 400d3e5ee88SLois Curfman McInnes int k, j, m, n; 401d3e5ee88SLois Curfman McInnes Scalar *v, tmp; 40248b35521SBarry Smith 403d3e5ee88SLois Curfman McInnes v = mat->v; m = mat->m; n = mat->n; 40448b35521SBarry Smith if (!matout) { /* in place transpose */ 405*48d91487SBarry Smith if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Not for rectangular matrix in place"); 406d3e5ee88SLois Curfman McInnes for ( j=0; j<m; j++ ) { 407289bc588SBarry Smith for ( k=0; k<j; k++ ) { 408d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 409d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 410d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 411289bc588SBarry Smith } 412289bc588SBarry Smith } 41348b35521SBarry Smith } 414d3e5ee88SLois Curfman McInnes else { /* out-of-place transpose */ 415d3e5ee88SLois Curfman McInnes int ierr; 416d3e5ee88SLois Curfman McInnes Mat tmat; 417ec8511deSBarry Smith Mat_SeqDense *tmatd; 418d3e5ee88SLois Curfman McInnes Scalar *v2; 419fafbff53SBarry Smith ierr = MatCreateSeqDense(matin->comm,mat->n,mat->m,&tmat); CHKERRQ(ierr); 420ec8511deSBarry Smith tmatd = (Mat_SeqDense *) tmat->data; 4210de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 422d3e5ee88SLois Curfman McInnes for ( j=0; j<n; j++ ) { 4230de55854SLois Curfman McInnes for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 424d3e5ee88SLois Curfman McInnes } 425d3e5ee88SLois Curfman McInnes ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 426d3e5ee88SLois Curfman McInnes ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 427d3e5ee88SLois Curfman McInnes *matout = tmat; 42848b35521SBarry Smith } 429289bc588SBarry Smith return 0; 430289bc588SBarry Smith } 431289bc588SBarry Smith 432ec8511deSBarry Smith static int MatEqual_SeqDense(Mat matin1,Mat matin2) 433289bc588SBarry Smith { 434ec8511deSBarry Smith Mat_SeqDense *mat1 = (Mat_SeqDense *) matin1->data; 435ec8511deSBarry Smith Mat_SeqDense *mat2 = (Mat_SeqDense *) matin2->data; 436289bc588SBarry Smith int i; 437289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 438289bc588SBarry Smith if (mat1->m != mat2->m) return 0; 439289bc588SBarry Smith if (mat1->n != mat2->n) return 0; 440289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 441289bc588SBarry Smith if (*v1 != *v2) return 0; 442289bc588SBarry Smith v1++; v2++; 443289bc588SBarry Smith } 444289bc588SBarry Smith return 1; 445289bc588SBarry Smith } 446289bc588SBarry Smith 447ec8511deSBarry Smith static int MatGetDiagonal_SeqDense(Mat matin,Vec v) 448289bc588SBarry Smith { 449ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 4506abc6512SBarry Smith int i, n; 4516abc6512SBarry Smith Scalar *x; 452289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 453ec8511deSBarry Smith if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 454289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 455289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 456289bc588SBarry Smith } 457289bc588SBarry Smith return 0; 458289bc588SBarry Smith } 459289bc588SBarry Smith 460ec8511deSBarry Smith static int MatScale_SeqDense(Mat matin,Vec ll,Vec rr) 461289bc588SBarry Smith { 462ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 463da3a660dSBarry Smith Scalar *l,*r,x,*v; 464da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 46528988994SBarry Smith if (ll) { 466da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 467ec8511deSBarry Smith if (m != mat->m) SETERRQ(1,"MatScale_SeqDense:Left scaling vec wrong size"); 468da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 469da3a660dSBarry Smith x = l[i]; 470da3a660dSBarry Smith v = mat->v + i; 471da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 472da3a660dSBarry Smith } 473da3a660dSBarry Smith } 47428988994SBarry Smith if (rr) { 475da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 476ec8511deSBarry Smith if (n != mat->n) SETERRQ(1,"MatScale_SeqDense:Right scaling vec wrong size"); 477da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 478da3a660dSBarry Smith x = r[i]; 479da3a660dSBarry Smith v = mat->v + i*m; 480da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 481da3a660dSBarry Smith } 482da3a660dSBarry Smith } 483289bc588SBarry Smith return 0; 484289bc588SBarry Smith } 485289bc588SBarry Smith 486ec8511deSBarry Smith static int MatNorm_SeqDense(Mat matin,MatNormType type,double *norm) 487289bc588SBarry Smith { 488ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 489289bc588SBarry Smith Scalar *v = mat->v; 490289bc588SBarry Smith double sum = 0.0; 491289bc588SBarry Smith int i, j; 492289bc588SBarry Smith if (type == NORM_FROBENIUS) { 493289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 494289bc588SBarry Smith #if defined(PETSC_COMPLEX) 495289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 496289bc588SBarry Smith #else 497289bc588SBarry Smith sum += (*v)*(*v); v++; 498289bc588SBarry Smith #endif 499289bc588SBarry Smith } 500289bc588SBarry Smith *norm = sqrt(sum); 501289bc588SBarry Smith } 502289bc588SBarry Smith else if (type == NORM_1) { 503289bc588SBarry Smith *norm = 0.0; 504289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 505289bc588SBarry Smith sum = 0.0; 506289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 507289bc588SBarry Smith #if defined(PETSC_COMPLEX) 508289bc588SBarry Smith sum += abs(*v++); 509289bc588SBarry Smith #else 510289bc588SBarry Smith sum += fabs(*v++); 511289bc588SBarry Smith #endif 512289bc588SBarry Smith } 513289bc588SBarry Smith if (sum > *norm) *norm = sum; 514289bc588SBarry Smith } 515289bc588SBarry Smith } 516289bc588SBarry Smith else if (type == NORM_INFINITY) { 517289bc588SBarry Smith *norm = 0.0; 518289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 519289bc588SBarry Smith v = mat->v + j; 520289bc588SBarry Smith sum = 0.0; 521289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 522289bc588SBarry Smith #if defined(PETSC_COMPLEX) 523289bc588SBarry Smith sum += abs(*v); v += mat->m; 524289bc588SBarry Smith #else 525289bc588SBarry Smith sum += fabs(*v); v += mat->m; 526289bc588SBarry Smith #endif 527289bc588SBarry Smith } 528289bc588SBarry Smith if (sum > *norm) *norm = sum; 529289bc588SBarry Smith } 530289bc588SBarry Smith } 531289bc588SBarry Smith else { 532*48d91487SBarry Smith SETERRQ(1,"MatNorm_SeqDense:No two norm"); 533289bc588SBarry Smith } 534289bc588SBarry Smith return 0; 535289bc588SBarry Smith } 536289bc588SBarry Smith 537ec8511deSBarry Smith static int MatSetOption_SeqDense(Mat aijin,MatOption op) 538289bc588SBarry Smith { 539ec8511deSBarry Smith Mat_SeqDense *aij = (Mat_SeqDense *) aijin->data; 540289bc588SBarry Smith if (op == ROW_ORIENTED) aij->roworiented = 1; 541289bc588SBarry Smith else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 542289bc588SBarry Smith /* doesn't care about sorted rows or columns */ 543289bc588SBarry Smith return 0; 544289bc588SBarry Smith } 545289bc588SBarry Smith 546ec8511deSBarry Smith static int MatZeroEntries_SeqDense(Mat A) 5476f0a148fSBarry Smith { 548ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 549416022c9SBarry Smith PetscZero(l->v,l->m*l->n*sizeof(Scalar)); 5506f0a148fSBarry Smith return 0; 5516f0a148fSBarry Smith } 5526f0a148fSBarry Smith 553ec8511deSBarry Smith static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 5546f0a148fSBarry Smith { 555ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 5566abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 5576f0a148fSBarry Smith Scalar *slot; 55878b31e54SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 55978b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 5606f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5616f0a148fSBarry Smith slot = l->v + rows[i]; 5626f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 5636f0a148fSBarry Smith } 5646f0a148fSBarry Smith if (diag) { 5656f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5666f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 567260925b8SBarry Smith *slot = *diag; 5686f0a148fSBarry Smith } 5696f0a148fSBarry Smith } 570260925b8SBarry Smith ISRestoreIndices(is,&rows); 5716f0a148fSBarry Smith return 0; 5726f0a148fSBarry Smith } 573557bce09SLois Curfman McInnes 574ec8511deSBarry Smith static int MatGetSize_SeqDense(Mat matin,int *m,int *n) 575557bce09SLois Curfman McInnes { 576ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 577557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 578557bce09SLois Curfman McInnes return 0; 579557bce09SLois Curfman McInnes } 580557bce09SLois Curfman McInnes 581ec8511deSBarry Smith static int MatGetOwnershipRange_SeqDense(Mat matin,int *m,int *n) 582d3e5ee88SLois Curfman McInnes { 583ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 584d3e5ee88SLois Curfman McInnes *m = 0; *n = mat->m; 585d3e5ee88SLois Curfman McInnes return 0; 586d3e5ee88SLois Curfman McInnes } 587d3e5ee88SLois Curfman McInnes 588ec8511deSBarry Smith static int MatGetArray_SeqDense(Mat matin,Scalar **array) 58964e87e97SBarry Smith { 590ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 59164e87e97SBarry Smith *array = mat->v; 59264e87e97SBarry Smith return 0; 59364e87e97SBarry Smith } 5940754003eSLois Curfman McInnes 5950754003eSLois Curfman McInnes 596ec8511deSBarry Smith static int MatGetSubMatrixInPlace_SeqDense(Mat matin,IS isrow,IS iscol) 5970754003eSLois Curfman McInnes { 598ec8511deSBarry Smith SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done"); 5990754003eSLois Curfman McInnes } 6000754003eSLois Curfman McInnes 601ec8511deSBarry Smith static int MatGetSubMatrix_SeqDense(Mat matin,IS isrow,IS iscol,Mat *submat) 6020754003eSLois Curfman McInnes { 603ec8511deSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) matin->data; 6040754003eSLois Curfman McInnes int nznew, *smap, i, j, ierr, oldcols = mat->n; 605160018dcSBarry Smith int *irow, *icol, nrows, ncols, *cwork; 6060754003eSLois Curfman McInnes Scalar *vwork, *val; 6070754003eSLois Curfman McInnes Mat newmat; 6080754003eSLois Curfman McInnes 60978b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 61078b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 61178b31e54SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 61278b31e54SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 6130754003eSLois Curfman McInnes 61478b31e54SBarry Smith smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap); 61578b31e54SBarry Smith cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork); 61678b31e54SBarry Smith vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 617416022c9SBarry Smith PetscZero((char*)smap,oldcols*sizeof(int)); 6180754003eSLois Curfman McInnes for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 6190754003eSLois Curfman McInnes 6200754003eSLois Curfman McInnes /* Create and fill new matrix */ 621*48d91487SBarry Smith ierr = MatCreateSeqDense(matin->comm,nrows,ncols,&newmat);CHKERRQ(ierr); 6220754003eSLois Curfman McInnes for (i=0; i<nrows; i++) { 6230754003eSLois Curfman McInnes nznew = 0; 6240754003eSLois Curfman McInnes val = mat->v + irow[i]; 6250754003eSLois Curfman McInnes for (j=0; j<oldcols; j++) { 6260754003eSLois Curfman McInnes if (smap[j]) { 6270754003eSLois Curfman McInnes cwork[nznew] = smap[j] - 1; 6280754003eSLois Curfman McInnes vwork[nznew++] = val[j * mat->m]; 6290754003eSLois Curfman McInnes } 6300754003eSLois Curfman McInnes } 631dbb450caSBarry Smith ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES); 63278b31e54SBarry Smith CHKERRQ(ierr); 6330754003eSLois Curfman McInnes } 63478b31e54SBarry Smith ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 63578b31e54SBarry Smith ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 6360754003eSLois Curfman McInnes 6370754003eSLois Curfman McInnes /* Free work space */ 63878b31e54SBarry Smith PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork); 63978b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 64078b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 6410754003eSLois Curfman McInnes *submat = newmat; 6420754003eSLois Curfman McInnes return 0; 6430754003eSLois Curfman McInnes } 6440754003eSLois Curfman McInnes 645289bc588SBarry Smith /* -------------------------------------------------------------------*/ 646ec8511deSBarry Smith static struct _MatOps MatOps = {MatInsert_SeqDense, 647ec8511deSBarry Smith MatGetRow_SeqDense, MatRestoreRow_SeqDense, 648ec8511deSBarry Smith MatMult_SeqDense, MatMultAdd_SeqDense, 649ec8511deSBarry Smith MatMultTrans_SeqDense, MatMultTransAdd_SeqDense, 650ec8511deSBarry Smith MatSolve_SeqDense,MatSolveAdd_SeqDense, 651ec8511deSBarry Smith MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense, 652ec8511deSBarry Smith MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense, 653ec8511deSBarry Smith MatRelax_SeqDense, 654ec8511deSBarry Smith MatTranspose_SeqDense, 655ec8511deSBarry Smith MatGetInfo_SeqDense,MatEqual_SeqDense, 656ec8511deSBarry Smith MatGetDiagonal_SeqDense,MatScale_SeqDense,MatNorm_SeqDense, 657289bc588SBarry Smith 0,0, 658ec8511deSBarry Smith 0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0, 659ec8511deSBarry Smith MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense, 660ec8511deSBarry Smith MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense, 661ec8511deSBarry Smith MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense, 662ec8511deSBarry Smith 0,0,MatGetArray_SeqDense,0,0, 663ec8511deSBarry Smith MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense, 664ec8511deSBarry Smith MatCopyPrivate_SeqDense}; 6650754003eSLois Curfman McInnes 6664b828684SBarry Smith /*@C 667fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 668d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 669d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 670289bc588SBarry Smith 67120563c6bSBarry Smith Input Parameters: 6720c775827SLois Curfman McInnes . comm - MPI communicator, set to MPI_COMM_SELF 6730c775827SLois Curfman McInnes . m - number of rows 6740c775827SLois Curfman McInnes . n - number of column 67520563c6bSBarry Smith 67620563c6bSBarry Smith Output Parameter: 6770c775827SLois Curfman McInnes . newmat - the matrix 67820563c6bSBarry Smith 679dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 680d65003e9SLois Curfman McInnes 681dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues() 68220563c6bSBarry Smith @*/ 683fafbff53SBarry Smith int MatCreateSeqDense(MPI_Comm comm,int m,int n,Mat *newmat) 684289bc588SBarry Smith { 685ec8511deSBarry Smith int size = sizeof(Mat_SeqDense) + m*n*sizeof(Scalar); 686289bc588SBarry Smith Mat mat; 687ec8511deSBarry Smith Mat_SeqDense *l; 688289bc588SBarry Smith *newmat = 0; 689ec8511deSBarry Smith PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm); 690a5a9c739SBarry Smith PLogObjectCreate(mat); 691ec8511deSBarry Smith l = (Mat_SeqDense *) PETSCMALLOC(size); CHKPTRQ(l); 692416022c9SBarry Smith PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 693ec8511deSBarry Smith mat->destroy = MatDestroy_SeqDense; 694ec8511deSBarry Smith mat->view = MatView_SeqDense; 695289bc588SBarry Smith mat->data = (void *) l; 696289bc588SBarry Smith mat->factor = 0; 697752f5567SLois Curfman McInnes PLogObjectMemory(mat,sizeof(struct _Mat) + size); 698d6dfbf8fSBarry Smith 699289bc588SBarry Smith l->m = m; 700289bc588SBarry Smith l->n = n; 701289bc588SBarry Smith l->v = (Scalar *) (l + 1); 702289bc588SBarry Smith l->pivots = 0; 703289bc588SBarry Smith l->roworiented = 1; 704d6dfbf8fSBarry Smith 705416022c9SBarry Smith PetscZero(l->v,m*n*sizeof(Scalar)); 706289bc588SBarry Smith *newmat = mat; 707289bc588SBarry Smith return 0; 708289bc588SBarry Smith } 709289bc588SBarry Smith 710ec8511deSBarry Smith int MatCreate_SeqDense(Mat matin,Mat *newmat) 711289bc588SBarry Smith { 712ec8511deSBarry Smith Mat_SeqDense *m = (Mat_SeqDense *) matin->data; 713fafbff53SBarry Smith return MatCreateSeqDense(matin->comm,m->m,m->n,newmat); 714289bc588SBarry Smith } 715