1cb512458SBarry Smith #ifndef lint 2*752f5567SLois Curfman McInnes static char vcid[] = "$Id: dense.c,v 1.49 1995/08/17 13:29:40 curfman Exp $"; 3cb512458SBarry Smith #endif 4289bc588SBarry Smith 5289bc588SBarry Smith /* 6289bc588SBarry Smith Standard Fortran style matrices 7289bc588SBarry Smith */ 819b02663SBarry Smith #include "petsc.h" 9289bc588SBarry Smith #include "plapack.h" 10289bc588SBarry Smith #include "matimpl.h" 11289bc588SBarry Smith #include "math.h" 12289bc588SBarry Smith #include "vec/vecimpl.h" 1320e1a0d4SLois Curfman McInnes #include "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++;} 29*752f5567SLois 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); 43*752f5567SLois 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; 75*752f5567SLois Curfman McInnes if (mat->pivots) { 76*752f5567SLois Curfman McInnes PETSCFREE(mat->pivots); 77*752f5567SLois Curfman McInnes PLogObjectMemory(matin,-mat->m*sizeof(int)); 78*752f5567SLois Curfman McInnes mat->pivots = 0; 79*752f5567SLois 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)); 112*752f5567SLois 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)); 137*752f5567SLois 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)); 164*752f5567SLois 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++ ) { 195289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 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-- ) { 201289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 202d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 203289bc588SBarry Smith } 204289bc588SBarry Smith } 205289bc588SBarry Smith } 206289bc588SBarry Smith return 0; 207289bc588SBarry Smith } 208289bc588SBarry Smith 209289bc588SBarry Smith /* -----------------------------------------------------------------*/ 21044a69424SLois Curfman McInnes static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy) 211289bc588SBarry Smith { 21244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 213289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 214289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 215289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 216289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 217289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 218289bc588SBarry Smith return 0; 219289bc588SBarry Smith } 22044a69424SLois Curfman McInnes static int MatMult_Dense(Mat matin,Vec xx,Vec yy) 221289bc588SBarry Smith { 22244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 223289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 224289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 225289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 226289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 227289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 228289bc588SBarry Smith return 0; 229289bc588SBarry Smith } 23044a69424SLois Curfman McInnes static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 231289bc588SBarry Smith { 23244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 233289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 2346abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 235289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 23678b31e54SBarry Smith if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar)); 237289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 238289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 239289bc588SBarry Smith return 0; 240289bc588SBarry Smith } 24144a69424SLois Curfman McInnes static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 242289bc588SBarry Smith { 24344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 244289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 245289bc588SBarry Smith int _One=1; 2466abc6512SBarry Smith Scalar _DOne=1.0; 247289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 248289bc588SBarry Smith VecGetArray(zz,&z); 24978b31e54SBarry Smith if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar)); 250289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 251289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 252289bc588SBarry Smith return 0; 253289bc588SBarry Smith } 254289bc588SBarry Smith 255289bc588SBarry Smith /* -----------------------------------------------------------------*/ 25644a69424SLois Curfman McInnes static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols, 257289bc588SBarry Smith Scalar **vals) 258289bc588SBarry Smith { 25944a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 260289bc588SBarry Smith Scalar *v; 261289bc588SBarry Smith int i; 262289bc588SBarry Smith *ncols = mat->n; 263289bc588SBarry Smith if (cols) { 26478b31e54SBarry Smith *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols); 265289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) *cols[i] = i; 266289bc588SBarry Smith } 267289bc588SBarry Smith if (vals) { 26878b31e54SBarry Smith *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 269289bc588SBarry Smith v = mat->v + row; 270289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 271289bc588SBarry Smith } 272289bc588SBarry Smith return 0; 273289bc588SBarry Smith } 27444a69424SLois Curfman McInnes static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols, 275289bc588SBarry Smith Scalar **vals) 276289bc588SBarry Smith { 27778b31e54SBarry Smith if (cols) { PETSCFREE(*cols); } 27878b31e54SBarry Smith if (vals) { PETSCFREE(*vals); } 279289bc588SBarry Smith return 0; 280289bc588SBarry Smith } 281289bc588SBarry Smith /* ----------------------------------------------------------------*/ 28244a69424SLois Curfman McInnes static int MatInsert_Dense(Mat matin,int m,int *indexm,int n, 283e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 284289bc588SBarry Smith { 28544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 286289bc588SBarry Smith int i,j; 287d6dfbf8fSBarry Smith 288289bc588SBarry Smith if (!mat->roworiented) { 289edae2e7dSBarry Smith if (addv == INSERTVALUES) { 290289bc588SBarry Smith for ( j=0; j<n; j++ ) { 291289bc588SBarry Smith for ( i=0; i<m; i++ ) { 292289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 293289bc588SBarry Smith } 294289bc588SBarry Smith } 295289bc588SBarry Smith } 296289bc588SBarry Smith else { 297289bc588SBarry Smith for ( j=0; j<n; j++ ) { 298289bc588SBarry Smith for ( i=0; i<m; i++ ) { 299289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 300289bc588SBarry Smith } 301289bc588SBarry Smith } 302289bc588SBarry Smith } 303e8d4e0b9SBarry Smith } 304e8d4e0b9SBarry Smith else { 305edae2e7dSBarry Smith if (addv == INSERTVALUES) { 306e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 307e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 308e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 309e8d4e0b9SBarry Smith } 310e8d4e0b9SBarry Smith } 311e8d4e0b9SBarry Smith } 312289bc588SBarry Smith else { 313289bc588SBarry Smith for ( i=0; i<m; i++ ) { 314289bc588SBarry Smith for ( j=0; j<n; j++ ) { 315289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 316289bc588SBarry Smith } 317289bc588SBarry Smith } 318289bc588SBarry Smith } 319e8d4e0b9SBarry Smith } 320289bc588SBarry Smith return 0; 321289bc588SBarry Smith } 322e8d4e0b9SBarry Smith 323289bc588SBarry Smith /* -----------------------------------------------------------------*/ 324ff7509f2SBarry Smith static int MatCopyPrivate_Dense(Mat matin,Mat *newmat) 325289bc588SBarry Smith { 32644a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 327289bc588SBarry Smith int ierr; 328289bc588SBarry Smith Mat newi; 32944a69424SLois Curfman McInnes Mat_Dense *l; 330bbb6d6a8SBarry Smith ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi); 331bbb6d6a8SBarry Smith CHKERRQ(ierr); 33244a69424SLois Curfman McInnes l = (Mat_Dense *) newi->data; 33378b31e54SBarry Smith PETSCMEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 334289bc588SBarry Smith *newmat = newi; 335289bc588SBarry Smith return 0; 336289bc588SBarry Smith } 337da3a660dSBarry Smith #include "viewer.h" 338289bc588SBarry Smith 33944a69424SLois Curfman McInnes int MatView_Dense(PetscObject obj,Viewer ptr) 340289bc588SBarry Smith { 341289bc588SBarry Smith Mat matin = (Mat) obj; 34244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 343289bc588SBarry Smith Scalar *v; 344289bc588SBarry Smith int i,j; 34523242f5aSBarry Smith PetscObject vobj = (PetscObject) ptr; 346da3a660dSBarry Smith 34723242f5aSBarry Smith if (ptr == 0) { 34823242f5aSBarry Smith ptr = STDOUT_VIEWER; vobj = (PetscObject) ptr; 34923242f5aSBarry Smith } 35023242f5aSBarry Smith if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) { 3516f75cfa4SBarry Smith return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v); 352da3a660dSBarry Smith } 353da3a660dSBarry Smith else { 3544a0ce102SLois Curfman McInnes FILE *fd = ViewerFileGetPointer_Private(ptr); 355289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 356289bc588SBarry Smith v = mat->v + i; 357289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 358289bc588SBarry Smith #if defined(PETSC_COMPLEX) 3598e37dc09SBarry Smith fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 360289bc588SBarry Smith #else 3618e37dc09SBarry Smith fprintf(fd,"%6.4e ",*v); v += mat->m; 362289bc588SBarry Smith #endif 363289bc588SBarry Smith } 3648e37dc09SBarry Smith fprintf(fd,"\n"); 365289bc588SBarry Smith } 366da3a660dSBarry Smith } 367289bc588SBarry Smith return 0; 368289bc588SBarry Smith } 369289bc588SBarry Smith 370289bc588SBarry Smith 37144a69424SLois Curfman McInnes static int MatDestroy_Dense(PetscObject obj) 372289bc588SBarry Smith { 373289bc588SBarry Smith Mat mat = (Mat) obj; 37444a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) mat->data; 375a5a9c739SBarry Smith #if defined(PETSC_LOG) 376a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 377a5a9c739SBarry Smith #endif 37878b31e54SBarry Smith if (l->pivots) PETSCFREE(l->pivots); 37978b31e54SBarry Smith PETSCFREE(l); 380a5a9c739SBarry Smith PLogObjectDestroy(mat); 3819e25ed09SBarry Smith PETSCHEADERDESTROY(mat); 382289bc588SBarry Smith return 0; 383289bc588SBarry Smith } 384289bc588SBarry Smith 38548b35521SBarry Smith static int MatTranspose_Dense(Mat matin,Mat *matout) 386289bc588SBarry Smith { 38744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 388289bc588SBarry Smith int k,j; 389289bc588SBarry Smith Scalar *v = mat->v, tmp; 39048b35521SBarry Smith 39148b35521SBarry Smith if (!matout) { /* in place transpose */ 392289bc588SBarry Smith if (mat->m != mat->n) { 39348b35521SBarry Smith SETERRQ(1,"MatTranspose_Dense:Cannot transpose rectangular matrix"); 394289bc588SBarry Smith } 395289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 396289bc588SBarry Smith for ( k=0; k<j; k++ ) { 397289bc588SBarry Smith tmp = v[j + k*mat->n]; 398289bc588SBarry Smith v[j + k*mat->n] = v[k + j*mat->n]; 399289bc588SBarry Smith v[k + j*mat->n] = tmp; 400289bc588SBarry Smith } 401289bc588SBarry Smith } 40248b35521SBarry Smith } 40348b35521SBarry Smith else { 40448b35521SBarry Smith SETERRQ(1,"MatTranspose_Dense:not out of place transpose yet"); 40548b35521SBarry Smith } 406289bc588SBarry Smith return 0; 407289bc588SBarry Smith } 408289bc588SBarry Smith 40944a69424SLois Curfman McInnes static int MatEqual_Dense(Mat matin1,Mat matin2) 410289bc588SBarry Smith { 41144a69424SLois Curfman McInnes Mat_Dense *mat1 = (Mat_Dense *) matin1->data; 41244a69424SLois Curfman McInnes Mat_Dense *mat2 = (Mat_Dense *) matin2->data; 413289bc588SBarry Smith int i; 414289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 415289bc588SBarry Smith if (mat1->m != mat2->m) return 0; 416289bc588SBarry Smith if (mat1->n != mat2->n) return 0; 417289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 418289bc588SBarry Smith if (*v1 != *v2) return 0; 419289bc588SBarry Smith v1++; v2++; 420289bc588SBarry Smith } 421289bc588SBarry Smith return 1; 422289bc588SBarry Smith } 423289bc588SBarry Smith 42446fc4a8fSLois Curfman McInnes static int MatGetDiagonal_Dense(Mat matin,Vec v) 425289bc588SBarry Smith { 42644a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 4276abc6512SBarry Smith int i, n; 4286abc6512SBarry Smith Scalar *x; 429289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 430bbb6d6a8SBarry Smith if (n != mat->m) SETERRQ(1,"MatGetDiagonal_Dense:Nonconforming mat and vec"); 431289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 432289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 433289bc588SBarry Smith } 434289bc588SBarry Smith return 0; 435289bc588SBarry Smith } 436289bc588SBarry Smith 43744a69424SLois Curfman McInnes static int MatScale_Dense(Mat matin,Vec ll,Vec rr) 438289bc588SBarry Smith { 43944a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 440da3a660dSBarry Smith Scalar *l,*r,x,*v; 441da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 44228988994SBarry Smith if (ll) { 443da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 444bbb6d6a8SBarry Smith if (m != mat->m) SETERRQ(1,"MatScale_Dense:Left scaling vec wrong size"); 445da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 446da3a660dSBarry Smith x = l[i]; 447da3a660dSBarry Smith v = mat->v + i; 448da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 449da3a660dSBarry Smith } 450da3a660dSBarry Smith } 45128988994SBarry Smith if (rr) { 452da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 453bbb6d6a8SBarry Smith if (n != mat->n) SETERRQ(1,"MatScale_Dense:Right scaling vec wrong size"); 454da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 455da3a660dSBarry Smith x = r[i]; 456da3a660dSBarry Smith v = mat->v + i*m; 457da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 458da3a660dSBarry Smith } 459da3a660dSBarry Smith } 460289bc588SBarry Smith return 0; 461289bc588SBarry Smith } 462289bc588SBarry Smith 4637c17b2cfSLois Curfman McInnes static int MatNorm_Dense(Mat matin,MatNormType type,double *norm) 464289bc588SBarry Smith { 46544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 466289bc588SBarry Smith Scalar *v = mat->v; 467289bc588SBarry Smith double sum = 0.0; 468289bc588SBarry Smith int i, j; 469289bc588SBarry Smith if (type == NORM_FROBENIUS) { 470289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 471289bc588SBarry Smith #if defined(PETSC_COMPLEX) 472289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 473289bc588SBarry Smith #else 474289bc588SBarry Smith sum += (*v)*(*v); v++; 475289bc588SBarry Smith #endif 476289bc588SBarry Smith } 477289bc588SBarry Smith *norm = sqrt(sum); 478289bc588SBarry Smith } 479289bc588SBarry Smith else if (type == NORM_1) { 480289bc588SBarry Smith *norm = 0.0; 481289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 482289bc588SBarry Smith sum = 0.0; 483289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 484289bc588SBarry Smith #if defined(PETSC_COMPLEX) 485289bc588SBarry Smith sum += abs(*v++); 486289bc588SBarry Smith #else 487289bc588SBarry Smith sum += fabs(*v++); 488289bc588SBarry Smith #endif 489289bc588SBarry Smith } 490289bc588SBarry Smith if (sum > *norm) *norm = sum; 491289bc588SBarry Smith } 492289bc588SBarry Smith } 493289bc588SBarry Smith else if (type == NORM_INFINITY) { 494289bc588SBarry Smith *norm = 0.0; 495289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 496289bc588SBarry Smith v = mat->v + j; 497289bc588SBarry Smith sum = 0.0; 498289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 499289bc588SBarry Smith #if defined(PETSC_COMPLEX) 500289bc588SBarry Smith sum += abs(*v); v += mat->m; 501289bc588SBarry Smith #else 502289bc588SBarry Smith sum += fabs(*v); v += mat->m; 503289bc588SBarry Smith #endif 504289bc588SBarry Smith } 505289bc588SBarry Smith if (sum > *norm) *norm = sum; 506289bc588SBarry Smith } 507289bc588SBarry Smith } 508289bc588SBarry Smith else { 509bbb6d6a8SBarry Smith SETERRQ(1,"MatNorm_Dense:No two norm yet"); 510289bc588SBarry Smith } 511289bc588SBarry Smith return 0; 512289bc588SBarry Smith } 513289bc588SBarry Smith 51420e1a0d4SLois Curfman McInnes static int MatSetOption_Dense(Mat aijin,MatOption op) 515289bc588SBarry Smith { 51644a69424SLois Curfman McInnes Mat_Dense *aij = (Mat_Dense *) aijin->data; 517289bc588SBarry Smith if (op == ROW_ORIENTED) aij->roworiented = 1; 518289bc588SBarry Smith else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 519289bc588SBarry Smith /* doesn't care about sorted rows or columns */ 520289bc588SBarry Smith return 0; 521289bc588SBarry Smith } 522289bc588SBarry Smith 52344a69424SLois Curfman McInnes static int MatZero_Dense(Mat A) 5246f0a148fSBarry Smith { 52544a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) A->data; 52678b31e54SBarry Smith PETSCMEMSET(l->v,0,l->m*l->n*sizeof(Scalar)); 5276f0a148fSBarry Smith return 0; 5286f0a148fSBarry Smith } 5296f0a148fSBarry Smith 53044a69424SLois Curfman McInnes static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag) 5316f0a148fSBarry Smith { 53244a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) A->data; 5336abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 5346f0a148fSBarry Smith Scalar *slot; 53578b31e54SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 53678b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 5376f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5386f0a148fSBarry Smith slot = l->v + rows[i]; 5396f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 5406f0a148fSBarry Smith } 5416f0a148fSBarry Smith if (diag) { 5426f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5436f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 544260925b8SBarry Smith *slot = *diag; 5456f0a148fSBarry Smith } 5466f0a148fSBarry Smith } 547260925b8SBarry Smith ISRestoreIndices(is,&rows); 5486f0a148fSBarry Smith return 0; 5496f0a148fSBarry Smith } 550557bce09SLois Curfman McInnes 551fa9ec3f1SLois Curfman McInnes static int MatGetSize_Dense(Mat matin,int *m,int *n) 552557bce09SLois Curfman McInnes { 55344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 554557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 555557bce09SLois Curfman McInnes return 0; 556557bce09SLois Curfman McInnes } 557557bce09SLois Curfman McInnes 55844a69424SLois Curfman McInnes static int MatGetArray_Dense(Mat matin,Scalar **array) 55964e87e97SBarry Smith { 56044a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 56164e87e97SBarry Smith *array = mat->v; 56264e87e97SBarry Smith return 0; 56364e87e97SBarry Smith } 5640754003eSLois Curfman McInnes 5650754003eSLois Curfman McInnes 5660754003eSLois Curfman McInnes static int MatGetSubMatrixInPlace_Dense(Mat matin,IS isrow,IS iscol) 5670754003eSLois Curfman McInnes { 568bbb6d6a8SBarry Smith SETERRQ(1,"MatGetSubMatrixInPlace_Dense: not done"); 5690754003eSLois Curfman McInnes } 5700754003eSLois Curfman McInnes 5710754003eSLois Curfman McInnes static int MatGetSubMatrix_Dense(Mat matin,IS isrow,IS iscol,Mat *submat) 5720754003eSLois Curfman McInnes { 5730754003eSLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 5740754003eSLois Curfman McInnes int nznew, *smap, i, j, ierr, oldcols = mat->n; 575160018dcSBarry Smith int *irow, *icol, nrows, ncols, *cwork; 5760754003eSLois Curfman McInnes Scalar *vwork, *val; 5770754003eSLois Curfman McInnes Mat newmat; 5780754003eSLois Curfman McInnes 57978b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 58078b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 58178b31e54SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 58278b31e54SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 5830754003eSLois Curfman McInnes 58478b31e54SBarry Smith smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap); 58578b31e54SBarry Smith cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork); 58678b31e54SBarry Smith vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 58778b31e54SBarry Smith PETSCMEMSET((char*)smap,0,oldcols*sizeof(int)); 5880754003eSLois Curfman McInnes for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 5890754003eSLois Curfman McInnes 5900754003eSLois Curfman McInnes /* Create and fill new matrix */ 5910754003eSLois Curfman McInnes ierr = MatCreateSequentialDense(matin->comm,nrows,ncols,&newmat); 59278b31e54SBarry Smith CHKERRQ(ierr); 5930754003eSLois Curfman McInnes for (i=0; i<nrows; i++) { 5940754003eSLois Curfman McInnes nznew = 0; 5950754003eSLois Curfman McInnes val = mat->v + irow[i]; 5960754003eSLois Curfman McInnes for (j=0; j<oldcols; j++) { 5970754003eSLois Curfman McInnes if (smap[j]) { 5980754003eSLois Curfman McInnes cwork[nznew] = smap[j] - 1; 5990754003eSLois Curfman McInnes vwork[nznew++] = val[j * mat->m]; 6000754003eSLois Curfman McInnes } 6010754003eSLois Curfman McInnes } 602edae2e7dSBarry Smith ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES); 60378b31e54SBarry Smith CHKERRQ(ierr); 6040754003eSLois Curfman McInnes } 60578b31e54SBarry Smith ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 60678b31e54SBarry Smith ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 6070754003eSLois Curfman McInnes 6080754003eSLois Curfman McInnes /* Free work space */ 60978b31e54SBarry Smith PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork); 61078b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 61178b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 6120754003eSLois Curfman McInnes *submat = newmat; 6130754003eSLois Curfman McInnes return 0; 6140754003eSLois Curfman McInnes } 6150754003eSLois Curfman McInnes 616289bc588SBarry Smith /* -------------------------------------------------------------------*/ 61744a69424SLois Curfman McInnes static struct _MatOps MatOps = {MatInsert_Dense, 61844a69424SLois Curfman McInnes MatGetRow_Dense, MatRestoreRow_Dense, 61944a69424SLois Curfman McInnes MatMult_Dense, MatMultAdd_Dense, 62044a69424SLois Curfman McInnes MatMultTrans_Dense, MatMultTransAdd_Dense, 62144a69424SLois Curfman McInnes MatSolve_Dense,MatSolveAdd_Dense, 62244a69424SLois Curfman McInnes MatSolveTrans_Dense,MatSolveTransAdd_Dense, 62346fc4a8fSLois Curfman McInnes MatLUFactor_Dense,MatCholeskyFactor_Dense, 62444a69424SLois Curfman McInnes MatRelax_Dense, 62548b35521SBarry Smith MatTranspose_Dense, 626fa9ec3f1SLois Curfman McInnes MatGetInfo_Dense,MatEqual_Dense, 62746fc4a8fSLois Curfman McInnes MatGetDiagonal_Dense,MatScale_Dense,MatNorm_Dense, 628289bc588SBarry Smith 0,0, 629681d1451SLois Curfman McInnes 0, MatSetOption_Dense,MatZero_Dense,MatZeroRows_Dense,0, 63044a69424SLois Curfman McInnes MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense, 63146fc4a8fSLois Curfman McInnes MatCholeskyFactorSymbolic_Dense,MatCholeskyFactorNumeric_Dense, 632fa9ec3f1SLois Curfman McInnes MatGetSize_Dense,MatGetSize_Dense,0, 63307a0e7adSLois Curfman McInnes 0,0,MatGetArray_Dense,0,0, 63407a0e7adSLois Curfman McInnes MatGetSubMatrix_Dense,MatGetSubMatrixInPlace_Dense, 635ff7509f2SBarry Smith MatCopyPrivate_Dense}; 6360754003eSLois Curfman McInnes 63720563c6bSBarry Smith /*@ 63820563c6bSBarry Smith MatCreateSequentialDense - Creates a sequential dense matrix that 639d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 640d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 641289bc588SBarry Smith 64220563c6bSBarry Smith Input Parameters: 6430c775827SLois Curfman McInnes . comm - MPI communicator, set to MPI_COMM_SELF 6440c775827SLois Curfman McInnes . m - number of rows 6450c775827SLois Curfman McInnes . n - number of column 64620563c6bSBarry Smith 64720563c6bSBarry Smith Output Parameter: 6480c775827SLois Curfman McInnes . newmat - the matrix 64920563c6bSBarry Smith 650dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 651d65003e9SLois Curfman McInnes 652dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues() 65320563c6bSBarry Smith @*/ 6546b5873e3SBarry Smith int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat) 655289bc588SBarry Smith { 65644a69424SLois Curfman McInnes int size = sizeof(Mat_Dense) + m*n*sizeof(Scalar); 657289bc588SBarry Smith Mat mat; 65844a69424SLois Curfman McInnes Mat_Dense *l; 659289bc588SBarry Smith *newmat = 0; 6606b5873e3SBarry Smith PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm); 661a5a9c739SBarry Smith PLogObjectCreate(mat); 66278b31e54SBarry Smith l = (Mat_Dense *) PETSCMALLOC(size); CHKPTRQ(l); 663289bc588SBarry Smith mat->ops = &MatOps; 66444a69424SLois Curfman McInnes mat->destroy = MatDestroy_Dense; 66544a69424SLois Curfman McInnes mat->view = MatView_Dense; 666289bc588SBarry Smith mat->data = (void *) l; 667289bc588SBarry Smith mat->factor = 0; 668*752f5567SLois Curfman McInnes PLogObjectMemory(mat,sizeof(struct _Mat) + size); 669d6dfbf8fSBarry Smith 670289bc588SBarry Smith l->m = m; 671289bc588SBarry Smith l->n = n; 672289bc588SBarry Smith l->v = (Scalar *) (l + 1); 673289bc588SBarry Smith l->pivots = 0; 674289bc588SBarry Smith l->roworiented = 1; 675d6dfbf8fSBarry Smith 67678b31e54SBarry Smith PETSCMEMSET(l->v,0,m*n*sizeof(Scalar)); 677289bc588SBarry Smith *newmat = mat; 678289bc588SBarry Smith return 0; 679289bc588SBarry Smith } 680289bc588SBarry Smith 68144a69424SLois Curfman McInnes int MatCreate_Dense(Mat matin,Mat *newmat) 682289bc588SBarry Smith { 68344a69424SLois Curfman McInnes Mat_Dense *m = (Mat_Dense *) matin->data; 6846b5873e3SBarry Smith return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat); 685289bc588SBarry Smith } 686