1cb512458SBarry Smith #ifndef lint 2*d3e5ee88SLois Curfman McInnes static char vcid[] = "$Id: dense.c,v 1.54 1995/09/06 03:05:24 bsmith Exp curfman $"; 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 */ 2044a69424SLois Curfman McInnes } Mat_Dense; 21289bc588SBarry Smith 2220e1a0d4SLois Curfman McInnes static int MatGetInfo_Dense(Mat matin,MatInfoType flag,int *nz, 2320e1a0d4SLois Curfman McInnes int *nzalloc,int *mem) 24289bc588SBarry Smith { 2544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 26289bc588SBarry Smith int i,N = mat->m*mat->n,count = 0; 27289bc588SBarry Smith Scalar *v = mat->v; 28289bc588SBarry Smith for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 29752f5567SLois Curfman McInnes *nz = count; *nzalloc = N; *mem = (int)matin->mem; 30fa9ec3f1SLois Curfman McInnes return 0; 31289bc588SBarry Smith } 32289bc588SBarry Smith 33289bc588SBarry Smith /* ---------------------------------------------------------------*/ 34289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots, 35289bc588SBarry Smith rather than put it in the Mat->col slot.*/ 3649d8b64dSBarry Smith static int MatLUFactor_Dense(Mat matin,IS row,IS col,double f) 37289bc588SBarry Smith { 3844a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 396abc6512SBarry Smith int info; 40289bc588SBarry Smith if (!mat->pivots) { 4178b31e54SBarry Smith mat->pivots = (int *) PETSCMALLOC(mat->m*sizeof(int)); 4278b31e54SBarry Smith CHKPTRQ(mat->pivots); 43752f5567SLois Curfman McInnes PLogObjectMemory(matin,mat->m*sizeof(int)); 44289bc588SBarry Smith } 45289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 46bbb6d6a8SBarry Smith if (info) SETERRQ(1,"MatLUFactor_Dense: Bad LU factorization"); 47289bc588SBarry Smith matin->factor = FACTOR_LU; 48289bc588SBarry Smith return 0; 49289bc588SBarry Smith } 5049d8b64dSBarry Smith static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,double f, 5149d8b64dSBarry Smith Mat *fact) 52289bc588SBarry Smith { 53289bc588SBarry Smith int ierr; 54bbb6d6a8SBarry Smith ierr = MatConvert(matin,MATSAME,fact); CHKERRQ(ierr); 55289bc588SBarry Smith return 0; 56289bc588SBarry Smith } 5744a69424SLois Curfman McInnes static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact) 58289bc588SBarry Smith { 5949d8b64dSBarry Smith return MatLUFactor(*fact,0,0,1.0); 60289bc588SBarry Smith } 6149d8b64dSBarry Smith static int MatCholeskyFactorSymbolic_Dense(Mat matin,IS row,double f,Mat *fact) 62289bc588SBarry Smith { 63289bc588SBarry Smith int ierr; 64bbb6d6a8SBarry Smith ierr = MatConvert(matin,MATSAME,fact); CHKERRQ(ierr); 65289bc588SBarry Smith return 0; 66289bc588SBarry Smith } 6746fc4a8fSLois Curfman McInnes static int MatCholeskyFactorNumeric_Dense(Mat matin,Mat *fact) 68289bc588SBarry Smith { 6949d8b64dSBarry Smith return MatCholeskyFactor(*fact,0,1.0); 70289bc588SBarry Smith } 7149d8b64dSBarry Smith static int MatCholeskyFactor_Dense(Mat matin,IS perm,double f) 72289bc588SBarry Smith { 7344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 746abc6512SBarry Smith int info; 75752f5567SLois Curfman McInnes if (mat->pivots) { 76752f5567SLois Curfman McInnes PETSCFREE(mat->pivots); 77752f5567SLois Curfman McInnes PLogObjectMemory(matin,-mat->m*sizeof(int)); 78752f5567SLois Curfman McInnes mat->pivots = 0; 79752f5567SLois Curfman McInnes } 80289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 81bbb6d6a8SBarry Smith if (info) SETERRQ(1,"MatCholeskyFactor_Dense: Bad factorization"); 82289bc588SBarry Smith matin->factor = FACTOR_CHOLESKY; 83289bc588SBarry Smith return 0; 84289bc588SBarry Smith } 85289bc588SBarry Smith 8644a69424SLois Curfman McInnes static int MatSolve_Dense(Mat matin,Vec xx,Vec yy) 87289bc588SBarry Smith { 8844a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 896abc6512SBarry Smith int one = 1, info; 906abc6512SBarry Smith Scalar *x, *y; 91289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 9278b31e54SBarry Smith PETSCMEMCPY(y,x,mat->m*sizeof(Scalar)); 939e25ed09SBarry Smith if (matin->factor == FACTOR_LU) { 94289bc588SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 95289bc588SBarry Smith y, &mat->m, &info ); 96289bc588SBarry Smith } 979e25ed09SBarry Smith else if (matin->factor == FACTOR_CHOLESKY){ 98289bc588SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 99289bc588SBarry Smith y, &mat->m, &info ); 100289bc588SBarry Smith } 101bbb6d6a8SBarry Smith else SETERRQ(1,"MatSolve_Dense: Matrix must be factored to solve"); 102bbb6d6a8SBarry Smith if (info) SETERRQ(1,"MatSolve_Dense: Bad solve"); 103289bc588SBarry Smith return 0; 104289bc588SBarry Smith } 10544a69424SLois Curfman McInnes static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy) 106da3a660dSBarry Smith { 10744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1086abc6512SBarry Smith int one = 1, info; 1096abc6512SBarry Smith Scalar *x, *y; 110da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 11178b31e54SBarry Smith PETSCMEMCPY(y,x,mat->m*sizeof(Scalar)); 112752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 113da3a660dSBarry Smith if (mat->pivots) { 114da3a660dSBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 115da3a660dSBarry Smith y, &mat->m, &info ); 116da3a660dSBarry Smith } 117da3a660dSBarry Smith else { 118da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 119da3a660dSBarry Smith y, &mat->m, &info ); 120da3a660dSBarry Smith } 121bbb6d6a8SBarry Smith if (info) SETERRQ(1,"MatSolveTrans_Dense: Bad solve"); 122da3a660dSBarry Smith return 0; 123da3a660dSBarry Smith } 12444a69424SLois Curfman McInnes static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 125da3a660dSBarry Smith { 12644a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1276abc6512SBarry Smith int one = 1, info,ierr; 1286abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 129da3a660dSBarry Smith Vec tmp = 0; 130da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 131da3a660dSBarry Smith if (yy == zz) { 13278b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 133464493b3SBarry Smith PLogObjectParent(matin,tmp); 13478b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 135da3a660dSBarry Smith } 13678b31e54SBarry Smith PETSCMEMCPY(y,x,mat->m*sizeof(Scalar)); 137752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 138da3a660dSBarry Smith if (mat->pivots) { 139da3a660dSBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 140da3a660dSBarry Smith y, &mat->m, &info ); 141da3a660dSBarry Smith } 142da3a660dSBarry Smith else { 143da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 144da3a660dSBarry Smith y, &mat->m, &info ); 145da3a660dSBarry Smith } 146bbb6d6a8SBarry Smith if (info) SETERRQ(1,"MatSolveAdd_Dense: Bad solve"); 147da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 148da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 149da3a660dSBarry Smith return 0; 150da3a660dSBarry Smith } 15144a69424SLois Curfman McInnes static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy) 152da3a660dSBarry Smith { 15344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1546abc6512SBarry Smith int one = 1, info,ierr; 1556abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 156da3a660dSBarry Smith Vec tmp; 157da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 158da3a660dSBarry Smith if (yy == zz) { 15978b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 160464493b3SBarry Smith PLogObjectParent(matin,tmp); 16178b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 162da3a660dSBarry Smith } 16378b31e54SBarry Smith PETSCMEMCPY(y,x,mat->m*sizeof(Scalar)); 164752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 165da3a660dSBarry Smith if (mat->pivots) { 166da3a660dSBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 167da3a660dSBarry Smith y, &mat->m, &info ); 168da3a660dSBarry Smith } 169da3a660dSBarry Smith else { 170da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 171da3a660dSBarry Smith y, &mat->m, &info ); 172da3a660dSBarry Smith } 173bbb6d6a8SBarry Smith if (info) SETERRQ(1,"MatSolveTransAdd_Dense: Bad solve"); 174da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 175da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 176da3a660dSBarry Smith return 0; 177da3a660dSBarry Smith } 178289bc588SBarry Smith /* ------------------------------------------------------------------*/ 17920e1a0d4SLois Curfman McInnes static int MatRelax_Dense(Mat matin,Vec bb,double omega,MatSORType flag, 18020e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 181289bc588SBarry Smith { 18244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 183289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 1846abc6512SBarry Smith int o = 1,ierr, m = mat->m, i; 185289bc588SBarry Smith 186289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 187289bc588SBarry Smith /* this is a hack fix, should have another version without 188289bc588SBarry Smith the second BLdot */ 189bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx); CHKERRQ(ierr); 190289bc588SBarry Smith } 191289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 192289bc588SBarry Smith while (its--) { 193289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 194289bc588SBarry Smith for ( i=0; i<m; i++ ) { 195f1747703SBarry Smith #if defined(PETSC_COMPLEX) 196f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 197f1747703SBarry Smith not happy about returning a double complex */ 198f1747703SBarry Smith int _i; 199f1747703SBarry Smith Scalar sum = b[i]; 200f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 201f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 202f1747703SBarry Smith } 203f1747703SBarry Smith xt = sum; 204f1747703SBarry Smith #else 205289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 206f1747703SBarry Smith #endif 207d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 208289bc588SBarry Smith } 209289bc588SBarry Smith } 210289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 211289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 212f1747703SBarry Smith #if defined(PETSC_COMPLEX) 213f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 214f1747703SBarry Smith not happy about returning a double complex */ 215f1747703SBarry Smith int _i; 216f1747703SBarry Smith Scalar sum = b[i]; 217f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 218f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 219f1747703SBarry Smith } 220f1747703SBarry Smith xt = sum; 221f1747703SBarry Smith #else 222289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 223f1747703SBarry Smith #endif 224d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 225289bc588SBarry Smith } 226289bc588SBarry Smith } 227289bc588SBarry Smith } 228289bc588SBarry Smith return 0; 229289bc588SBarry Smith } 230289bc588SBarry Smith 231289bc588SBarry Smith /* -----------------------------------------------------------------*/ 23244a69424SLois Curfman McInnes static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy) 233289bc588SBarry Smith { 23444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 235289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 236289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 237289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 238289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 239289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 240289bc588SBarry Smith return 0; 241289bc588SBarry Smith } 24244a69424SLois Curfman McInnes static int MatMult_Dense(Mat matin,Vec xx,Vec yy) 243289bc588SBarry Smith { 24444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 245289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 246289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 247289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 248289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 249289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 250289bc588SBarry Smith return 0; 251289bc588SBarry Smith } 25244a69424SLois Curfman McInnes static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 253289bc588SBarry Smith { 25444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 255289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 2566abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 257289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 25878b31e54SBarry Smith if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar)); 259289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 260289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 261289bc588SBarry Smith return 0; 262289bc588SBarry Smith } 26344a69424SLois Curfman McInnes static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 264289bc588SBarry Smith { 26544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 266289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 267289bc588SBarry Smith int _One=1; 2686abc6512SBarry Smith Scalar _DOne=1.0; 269289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 270289bc588SBarry Smith VecGetArray(zz,&z); 27178b31e54SBarry Smith if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar)); 272289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 273289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 274289bc588SBarry Smith return 0; 275289bc588SBarry Smith } 276289bc588SBarry Smith 277289bc588SBarry Smith /* -----------------------------------------------------------------*/ 27844a69424SLois Curfman McInnes static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols, 279289bc588SBarry Smith Scalar **vals) 280289bc588SBarry Smith { 28144a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 282289bc588SBarry Smith Scalar *v; 283289bc588SBarry Smith int i; 284289bc588SBarry Smith *ncols = mat->n; 285289bc588SBarry Smith if (cols) { 28678b31e54SBarry Smith *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols); 287289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) *cols[i] = i; 288289bc588SBarry Smith } 289289bc588SBarry Smith if (vals) { 29078b31e54SBarry Smith *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 291289bc588SBarry Smith v = mat->v + row; 292289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 293289bc588SBarry Smith } 294289bc588SBarry Smith return 0; 295289bc588SBarry Smith } 29644a69424SLois Curfman McInnes static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols, 297289bc588SBarry Smith Scalar **vals) 298289bc588SBarry Smith { 29978b31e54SBarry Smith if (cols) { PETSCFREE(*cols); } 30078b31e54SBarry Smith if (vals) { PETSCFREE(*vals); } 301289bc588SBarry Smith return 0; 302289bc588SBarry Smith } 303289bc588SBarry Smith /* ----------------------------------------------------------------*/ 30444a69424SLois Curfman McInnes static int MatInsert_Dense(Mat matin,int m,int *indexm,int n, 305e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 306289bc588SBarry Smith { 30744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 308289bc588SBarry Smith int i,j; 309d6dfbf8fSBarry Smith 310289bc588SBarry Smith if (!mat->roworiented) { 311edae2e7dSBarry Smith if (addv == INSERTVALUES) { 312289bc588SBarry Smith for ( j=0; j<n; j++ ) { 313289bc588SBarry Smith for ( i=0; i<m; i++ ) { 314289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 315289bc588SBarry Smith } 316289bc588SBarry Smith } 317289bc588SBarry Smith } 318289bc588SBarry Smith else { 319289bc588SBarry Smith for ( j=0; j<n; j++ ) { 320289bc588SBarry Smith for ( i=0; i<m; i++ ) { 321289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 322289bc588SBarry Smith } 323289bc588SBarry Smith } 324289bc588SBarry Smith } 325e8d4e0b9SBarry Smith } 326e8d4e0b9SBarry Smith else { 327edae2e7dSBarry Smith if (addv == INSERTVALUES) { 328e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 329e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 330e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 331e8d4e0b9SBarry Smith } 332e8d4e0b9SBarry Smith } 333e8d4e0b9SBarry Smith } 334289bc588SBarry Smith else { 335289bc588SBarry Smith for ( i=0; i<m; i++ ) { 336289bc588SBarry Smith for ( j=0; j<n; j++ ) { 337289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 338289bc588SBarry Smith } 339289bc588SBarry Smith } 340289bc588SBarry Smith } 341e8d4e0b9SBarry Smith } 342289bc588SBarry Smith return 0; 343289bc588SBarry Smith } 344e8d4e0b9SBarry Smith 345289bc588SBarry Smith /* -----------------------------------------------------------------*/ 346ff7509f2SBarry Smith static int MatCopyPrivate_Dense(Mat matin,Mat *newmat) 347289bc588SBarry Smith { 34844a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 349289bc588SBarry Smith int ierr; 350289bc588SBarry Smith Mat newi; 35144a69424SLois Curfman McInnes Mat_Dense *l; 352bbb6d6a8SBarry Smith ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi); 353bbb6d6a8SBarry Smith CHKERRQ(ierr); 35444a69424SLois Curfman McInnes l = (Mat_Dense *) newi->data; 35578b31e54SBarry Smith PETSCMEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 356289bc588SBarry Smith *newmat = newi; 357289bc588SBarry Smith return 0; 358289bc588SBarry Smith } 359da3a660dSBarry Smith #include "viewer.h" 360289bc588SBarry Smith 36144a69424SLois Curfman McInnes int MatView_Dense(PetscObject obj,Viewer ptr) 362289bc588SBarry Smith { 363289bc588SBarry Smith Mat matin = (Mat) obj; 36444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 365289bc588SBarry Smith Scalar *v; 366d636dbe3SBarry Smith int i,j,ierr; 36723242f5aSBarry Smith PetscObject vobj = (PetscObject) ptr; 368da3a660dSBarry Smith 36923242f5aSBarry Smith if (ptr == 0) { 37021b0d8fbSLois Curfman McInnes ptr = STDOUT_VIEWER_SELF; vobj = (PetscObject) ptr; 37123242f5aSBarry Smith } 37223242f5aSBarry Smith if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) { 3736f75cfa4SBarry Smith return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v); 374da3a660dSBarry Smith } 375da3a660dSBarry Smith else { 376d636dbe3SBarry Smith FILE *fd; 377d636dbe3SBarry Smith int format; 378d636dbe3SBarry Smith ierr = ViewerFileGetPointer_Private(ptr,&fd); CHKERRQ(ierr); 379d636dbe3SBarry Smith ierr = ViewerFileGetFormat_Private(ptr,&format); CHKERRQ(ierr); 380f72585f2SLois Curfman McInnes if (format == FILE_FORMAT_INFO) { 381f72585f2SLois Curfman McInnes /* do nothing for now */ 382f72585f2SLois Curfman McInnes } 383f72585f2SLois Curfman McInnes else { 384289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 385289bc588SBarry Smith v = mat->v + i; 386289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 387289bc588SBarry Smith #if defined(PETSC_COMPLEX) 3888e37dc09SBarry Smith fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 389289bc588SBarry Smith #else 3908e37dc09SBarry Smith fprintf(fd,"%6.4e ",*v); v += mat->m; 391289bc588SBarry Smith #endif 392289bc588SBarry Smith } 3938e37dc09SBarry Smith fprintf(fd,"\n"); 394289bc588SBarry Smith } 395da3a660dSBarry Smith } 396f72585f2SLois Curfman McInnes } 397289bc588SBarry Smith return 0; 398289bc588SBarry Smith } 399289bc588SBarry Smith 400289bc588SBarry Smith 40144a69424SLois Curfman McInnes static int MatDestroy_Dense(PetscObject obj) 402289bc588SBarry Smith { 403289bc588SBarry Smith Mat mat = (Mat) obj; 40444a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) mat->data; 405a5a9c739SBarry Smith #if defined(PETSC_LOG) 406a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 407a5a9c739SBarry Smith #endif 40878b31e54SBarry Smith if (l->pivots) PETSCFREE(l->pivots); 40978b31e54SBarry Smith PETSCFREE(l); 410a5a9c739SBarry Smith PLogObjectDestroy(mat); 4119e25ed09SBarry Smith PETSCHEADERDESTROY(mat); 412289bc588SBarry Smith return 0; 413289bc588SBarry Smith } 414289bc588SBarry Smith 41548b35521SBarry Smith static int MatTranspose_Dense(Mat matin,Mat *matout) 416289bc588SBarry Smith { 41744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 418*d3e5ee88SLois Curfman McInnes int k, j, m, n; 419*d3e5ee88SLois Curfman McInnes Scalar *v, tmp; 42048b35521SBarry Smith 421*d3e5ee88SLois Curfman McInnes v = mat->v; m = mat->m; n = mat->n; 42248b35521SBarry Smith if (!matout) { /* in place transpose */ 423*d3e5ee88SLois Curfman McInnes if (m != n) SETERRQ(1, 424*d3e5ee88SLois Curfman McInnes "MatTranspose_Dense: Cannot transpose rectangular matrix in place"); 425*d3e5ee88SLois Curfman McInnes for ( j=0; j<m; j++ ) { 426289bc588SBarry Smith for ( k=0; k<j; k++ ) { 427*d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 428*d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 429*d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 430289bc588SBarry Smith } 431289bc588SBarry Smith } 43248b35521SBarry Smith } 433*d3e5ee88SLois Curfman McInnes else { /* out-of-place transpose */ 434*d3e5ee88SLois Curfman McInnes int ierr; 435*d3e5ee88SLois Curfman McInnes Mat tmat; 436*d3e5ee88SLois Curfman McInnes Mat_Dense *tmatd; 437*d3e5ee88SLois Curfman McInnes Scalar *v2; 438*d3e5ee88SLois Curfman McInnes ierr = MatCreateSequentialDense(matin->comm,mat->n,mat->m,&tmat); CHKERRQ(ierr); 439*d3e5ee88SLois Curfman McInnes tmatd = (Mat_Dense *) tmat->data; 440*d3e5ee88SLois Curfman McInnes v = mat->v; v2 = tmatd->v; m = tmatd->m; n = tmatd->n; 441*d3e5ee88SLois Curfman McInnes for ( j=0; j<n; j++ ) { 442*d3e5ee88SLois Curfman McInnes for ( k=0; k<m; k++ ) { 443*d3e5ee88SLois Curfman McInnes v2[j + k*n] = v[k + j*m]; 444*d3e5ee88SLois Curfman McInnes } 445*d3e5ee88SLois Curfman McInnes } 446*d3e5ee88SLois Curfman McInnes ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 447*d3e5ee88SLois Curfman McInnes ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 448*d3e5ee88SLois Curfman McInnes *matout = tmat; 44948b35521SBarry Smith } 450289bc588SBarry Smith return 0; 451289bc588SBarry Smith } 452289bc588SBarry Smith 45344a69424SLois Curfman McInnes static int MatEqual_Dense(Mat matin1,Mat matin2) 454289bc588SBarry Smith { 45544a69424SLois Curfman McInnes Mat_Dense *mat1 = (Mat_Dense *) matin1->data; 45644a69424SLois Curfman McInnes Mat_Dense *mat2 = (Mat_Dense *) matin2->data; 457289bc588SBarry Smith int i; 458289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 459289bc588SBarry Smith if (mat1->m != mat2->m) return 0; 460289bc588SBarry Smith if (mat1->n != mat2->n) return 0; 461289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 462289bc588SBarry Smith if (*v1 != *v2) return 0; 463289bc588SBarry Smith v1++; v2++; 464289bc588SBarry Smith } 465289bc588SBarry Smith return 1; 466289bc588SBarry Smith } 467289bc588SBarry Smith 46846fc4a8fSLois Curfman McInnes static int MatGetDiagonal_Dense(Mat matin,Vec v) 469289bc588SBarry Smith { 47044a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 4716abc6512SBarry Smith int i, n; 4726abc6512SBarry Smith Scalar *x; 473289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 474bbb6d6a8SBarry Smith if (n != mat->m) SETERRQ(1,"MatGetDiagonal_Dense:Nonconforming mat and vec"); 475289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 476289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 477289bc588SBarry Smith } 478289bc588SBarry Smith return 0; 479289bc588SBarry Smith } 480289bc588SBarry Smith 48144a69424SLois Curfman McInnes static int MatScale_Dense(Mat matin,Vec ll,Vec rr) 482289bc588SBarry Smith { 48344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 484da3a660dSBarry Smith Scalar *l,*r,x,*v; 485da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 48628988994SBarry Smith if (ll) { 487da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 488bbb6d6a8SBarry Smith if (m != mat->m) SETERRQ(1,"MatScale_Dense:Left scaling vec wrong size"); 489da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 490da3a660dSBarry Smith x = l[i]; 491da3a660dSBarry Smith v = mat->v + i; 492da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 493da3a660dSBarry Smith } 494da3a660dSBarry Smith } 49528988994SBarry Smith if (rr) { 496da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 497bbb6d6a8SBarry Smith if (n != mat->n) SETERRQ(1,"MatScale_Dense:Right scaling vec wrong size"); 498da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 499da3a660dSBarry Smith x = r[i]; 500da3a660dSBarry Smith v = mat->v + i*m; 501da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 502da3a660dSBarry Smith } 503da3a660dSBarry Smith } 504289bc588SBarry Smith return 0; 505289bc588SBarry Smith } 506289bc588SBarry Smith 5077c17b2cfSLois Curfman McInnes static int MatNorm_Dense(Mat matin,MatNormType type,double *norm) 508289bc588SBarry Smith { 50944a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 510289bc588SBarry Smith Scalar *v = mat->v; 511289bc588SBarry Smith double sum = 0.0; 512289bc588SBarry Smith int i, j; 513289bc588SBarry Smith if (type == NORM_FROBENIUS) { 514289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 515289bc588SBarry Smith #if defined(PETSC_COMPLEX) 516289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 517289bc588SBarry Smith #else 518289bc588SBarry Smith sum += (*v)*(*v); v++; 519289bc588SBarry Smith #endif 520289bc588SBarry Smith } 521289bc588SBarry Smith *norm = sqrt(sum); 522289bc588SBarry Smith } 523289bc588SBarry Smith else if (type == NORM_1) { 524289bc588SBarry Smith *norm = 0.0; 525289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 526289bc588SBarry Smith sum = 0.0; 527289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 528289bc588SBarry Smith #if defined(PETSC_COMPLEX) 529289bc588SBarry Smith sum += abs(*v++); 530289bc588SBarry Smith #else 531289bc588SBarry Smith sum += fabs(*v++); 532289bc588SBarry Smith #endif 533289bc588SBarry Smith } 534289bc588SBarry Smith if (sum > *norm) *norm = sum; 535289bc588SBarry Smith } 536289bc588SBarry Smith } 537289bc588SBarry Smith else if (type == NORM_INFINITY) { 538289bc588SBarry Smith *norm = 0.0; 539289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 540289bc588SBarry Smith v = mat->v + j; 541289bc588SBarry Smith sum = 0.0; 542289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 543289bc588SBarry Smith #if defined(PETSC_COMPLEX) 544289bc588SBarry Smith sum += abs(*v); v += mat->m; 545289bc588SBarry Smith #else 546289bc588SBarry Smith sum += fabs(*v); v += mat->m; 547289bc588SBarry Smith #endif 548289bc588SBarry Smith } 549289bc588SBarry Smith if (sum > *norm) *norm = sum; 550289bc588SBarry Smith } 551289bc588SBarry Smith } 552289bc588SBarry Smith else { 553bbb6d6a8SBarry Smith SETERRQ(1,"MatNorm_Dense:No two norm yet"); 554289bc588SBarry Smith } 555289bc588SBarry Smith return 0; 556289bc588SBarry Smith } 557289bc588SBarry Smith 55820e1a0d4SLois Curfman McInnes static int MatSetOption_Dense(Mat aijin,MatOption op) 559289bc588SBarry Smith { 56044a69424SLois Curfman McInnes Mat_Dense *aij = (Mat_Dense *) aijin->data; 561289bc588SBarry Smith if (op == ROW_ORIENTED) aij->roworiented = 1; 562289bc588SBarry Smith else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 563289bc588SBarry Smith /* doesn't care about sorted rows or columns */ 564289bc588SBarry Smith return 0; 565289bc588SBarry Smith } 566289bc588SBarry Smith 567*d3e5ee88SLois Curfman McInnes static int MatZeroEntries_Dense(Mat A) 5686f0a148fSBarry Smith { 56944a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) A->data; 57078b31e54SBarry Smith PETSCMEMSET(l->v,0,l->m*l->n*sizeof(Scalar)); 5716f0a148fSBarry Smith return 0; 5726f0a148fSBarry Smith } 5736f0a148fSBarry Smith 57444a69424SLois Curfman McInnes static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag) 5756f0a148fSBarry Smith { 57644a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) A->data; 5776abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 5786f0a148fSBarry Smith Scalar *slot; 57978b31e54SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 58078b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 5816f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5826f0a148fSBarry Smith slot = l->v + rows[i]; 5836f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 5846f0a148fSBarry Smith } 5856f0a148fSBarry Smith if (diag) { 5866f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5876f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 588260925b8SBarry Smith *slot = *diag; 5896f0a148fSBarry Smith } 5906f0a148fSBarry Smith } 591260925b8SBarry Smith ISRestoreIndices(is,&rows); 5926f0a148fSBarry Smith return 0; 5936f0a148fSBarry Smith } 594557bce09SLois Curfman McInnes 595fa9ec3f1SLois Curfman McInnes static int MatGetSize_Dense(Mat matin,int *m,int *n) 596557bce09SLois Curfman McInnes { 59744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 598557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 599557bce09SLois Curfman McInnes return 0; 600557bce09SLois Curfman McInnes } 601557bce09SLois Curfman McInnes 602*d3e5ee88SLois Curfman McInnes static int MatGetOwnershipRange_Dense(Mat matin,int *m,int *n) 603*d3e5ee88SLois Curfman McInnes { 604*d3e5ee88SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 605*d3e5ee88SLois Curfman McInnes *m = 0; *n = mat->m; 606*d3e5ee88SLois Curfman McInnes return 0; 607*d3e5ee88SLois Curfman McInnes } 608*d3e5ee88SLois Curfman McInnes 60944a69424SLois Curfman McInnes static int MatGetArray_Dense(Mat matin,Scalar **array) 61064e87e97SBarry Smith { 61144a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 61264e87e97SBarry Smith *array = mat->v; 61364e87e97SBarry Smith return 0; 61464e87e97SBarry Smith } 6150754003eSLois Curfman McInnes 6160754003eSLois Curfman McInnes 6170754003eSLois Curfman McInnes static int MatGetSubMatrixInPlace_Dense(Mat matin,IS isrow,IS iscol) 6180754003eSLois Curfman McInnes { 619bbb6d6a8SBarry Smith SETERRQ(1,"MatGetSubMatrixInPlace_Dense: not done"); 6200754003eSLois Curfman McInnes } 6210754003eSLois Curfman McInnes 6220754003eSLois Curfman McInnes static int MatGetSubMatrix_Dense(Mat matin,IS isrow,IS iscol,Mat *submat) 6230754003eSLois Curfman McInnes { 6240754003eSLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 6250754003eSLois Curfman McInnes int nznew, *smap, i, j, ierr, oldcols = mat->n; 626160018dcSBarry Smith int *irow, *icol, nrows, ncols, *cwork; 6270754003eSLois Curfman McInnes Scalar *vwork, *val; 6280754003eSLois Curfman McInnes Mat newmat; 6290754003eSLois Curfman McInnes 63078b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 63178b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 63278b31e54SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 63378b31e54SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 6340754003eSLois Curfman McInnes 63578b31e54SBarry Smith smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap); 63678b31e54SBarry Smith cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork); 63778b31e54SBarry Smith vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 63878b31e54SBarry Smith PETSCMEMSET((char*)smap,0,oldcols*sizeof(int)); 6390754003eSLois Curfman McInnes for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 6400754003eSLois Curfman McInnes 6410754003eSLois Curfman McInnes /* Create and fill new matrix */ 6420754003eSLois Curfman McInnes ierr = MatCreateSequentialDense(matin->comm,nrows,ncols,&newmat); 64378b31e54SBarry Smith CHKERRQ(ierr); 6440754003eSLois Curfman McInnes for (i=0; i<nrows; i++) { 6450754003eSLois Curfman McInnes nznew = 0; 6460754003eSLois Curfman McInnes val = mat->v + irow[i]; 6470754003eSLois Curfman McInnes for (j=0; j<oldcols; j++) { 6480754003eSLois Curfman McInnes if (smap[j]) { 6490754003eSLois Curfman McInnes cwork[nznew] = smap[j] - 1; 6500754003eSLois Curfman McInnes vwork[nznew++] = val[j * mat->m]; 6510754003eSLois Curfman McInnes } 6520754003eSLois Curfman McInnes } 653edae2e7dSBarry Smith ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES); 65478b31e54SBarry Smith CHKERRQ(ierr); 6550754003eSLois Curfman McInnes } 65678b31e54SBarry Smith ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 65778b31e54SBarry Smith ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 6580754003eSLois Curfman McInnes 6590754003eSLois Curfman McInnes /* Free work space */ 66078b31e54SBarry Smith PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork); 66178b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 66278b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 6630754003eSLois Curfman McInnes *submat = newmat; 6640754003eSLois Curfman McInnes return 0; 6650754003eSLois Curfman McInnes } 6660754003eSLois Curfman McInnes 667289bc588SBarry Smith /* -------------------------------------------------------------------*/ 66844a69424SLois Curfman McInnes static struct _MatOps MatOps = {MatInsert_Dense, 66944a69424SLois Curfman McInnes MatGetRow_Dense, MatRestoreRow_Dense, 67044a69424SLois Curfman McInnes MatMult_Dense, MatMultAdd_Dense, 67144a69424SLois Curfman McInnes MatMultTrans_Dense, MatMultTransAdd_Dense, 67244a69424SLois Curfman McInnes MatSolve_Dense,MatSolveAdd_Dense, 67344a69424SLois Curfman McInnes MatSolveTrans_Dense,MatSolveTransAdd_Dense, 67446fc4a8fSLois Curfman McInnes MatLUFactor_Dense,MatCholeskyFactor_Dense, 67544a69424SLois Curfman McInnes MatRelax_Dense, 67648b35521SBarry Smith MatTranspose_Dense, 677fa9ec3f1SLois Curfman McInnes MatGetInfo_Dense,MatEqual_Dense, 67846fc4a8fSLois Curfman McInnes MatGetDiagonal_Dense,MatScale_Dense,MatNorm_Dense, 679289bc588SBarry Smith 0,0, 680*d3e5ee88SLois Curfman McInnes 0, MatSetOption_Dense,MatZeroEntries_Dense,MatZeroRows_Dense,0, 68144a69424SLois Curfman McInnes MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense, 68246fc4a8fSLois Curfman McInnes MatCholeskyFactorSymbolic_Dense,MatCholeskyFactorNumeric_Dense, 683*d3e5ee88SLois Curfman McInnes MatGetSize_Dense,MatGetSize_Dense,MatGetOwnershipRange_Dense, 68407a0e7adSLois Curfman McInnes 0,0,MatGetArray_Dense,0,0, 68507a0e7adSLois Curfman McInnes MatGetSubMatrix_Dense,MatGetSubMatrixInPlace_Dense, 686ff7509f2SBarry Smith MatCopyPrivate_Dense}; 6870754003eSLois Curfman McInnes 6884b828684SBarry Smith /*@C 68920563c6bSBarry Smith MatCreateSequentialDense - Creates a sequential dense matrix that 690d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 691d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 692289bc588SBarry Smith 69320563c6bSBarry Smith Input Parameters: 6940c775827SLois Curfman McInnes . comm - MPI communicator, set to MPI_COMM_SELF 6950c775827SLois Curfman McInnes . m - number of rows 6960c775827SLois Curfman McInnes . n - number of column 69720563c6bSBarry Smith 69820563c6bSBarry Smith Output Parameter: 6990c775827SLois Curfman McInnes . newmat - the matrix 70020563c6bSBarry Smith 701dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 702d65003e9SLois Curfman McInnes 703dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues() 70420563c6bSBarry Smith @*/ 7056b5873e3SBarry Smith int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat) 706289bc588SBarry Smith { 70744a69424SLois Curfman McInnes int size = sizeof(Mat_Dense) + m*n*sizeof(Scalar); 708289bc588SBarry Smith Mat mat; 70944a69424SLois Curfman McInnes Mat_Dense *l; 710289bc588SBarry Smith *newmat = 0; 7116b5873e3SBarry Smith PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm); 712a5a9c739SBarry Smith PLogObjectCreate(mat); 71378b31e54SBarry Smith l = (Mat_Dense *) PETSCMALLOC(size); CHKPTRQ(l); 714289bc588SBarry Smith mat->ops = &MatOps; 71544a69424SLois Curfman McInnes mat->destroy = MatDestroy_Dense; 71644a69424SLois Curfman McInnes mat->view = MatView_Dense; 717289bc588SBarry Smith mat->data = (void *) l; 718289bc588SBarry Smith mat->factor = 0; 719752f5567SLois Curfman McInnes PLogObjectMemory(mat,sizeof(struct _Mat) + size); 720d6dfbf8fSBarry Smith 721289bc588SBarry Smith l->m = m; 722289bc588SBarry Smith l->n = n; 723289bc588SBarry Smith l->v = (Scalar *) (l + 1); 724289bc588SBarry Smith l->pivots = 0; 725289bc588SBarry Smith l->roworiented = 1; 726d6dfbf8fSBarry Smith 72778b31e54SBarry Smith PETSCMEMSET(l->v,0,m*n*sizeof(Scalar)); 728289bc588SBarry Smith *newmat = mat; 729289bc588SBarry Smith return 0; 730289bc588SBarry Smith } 731289bc588SBarry Smith 73244a69424SLois Curfman McInnes int MatCreate_Dense(Mat matin,Mat *newmat) 733289bc588SBarry Smith { 73444a69424SLois Curfman McInnes Mat_Dense *m = (Mat_Dense *) matin->data; 7356b5873e3SBarry Smith return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat); 736289bc588SBarry Smith } 737