1cb512458SBarry Smith #ifndef lint 2*edae2e7dSBarry Smith static char vcid[] = "$Id: dense.c,v 1.32 1995/05/02 21:11:26 curfman Exp bsmith $"; 3cb512458SBarry Smith #endif 4289bc588SBarry Smith 5289bc588SBarry Smith /* 6289bc588SBarry Smith Standard Fortran style matrices 7289bc588SBarry Smith */ 8289bc588SBarry Smith #include "ptscimpl.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 22289bc588SBarry Smith 2320e1a0d4SLois Curfman McInnes static int MatGetInfo_Dense(Mat matin,MatInfoType flag,int *nz, 2420e1a0d4SLois Curfman McInnes int *nzalloc,int *mem) 25289bc588SBarry Smith { 2644a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 27289bc588SBarry Smith int i,N = mat->m*mat->n,count = 0; 28289bc588SBarry Smith Scalar *v = mat->v; 29289bc588SBarry Smith for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 30fa9ec3f1SLois Curfman McInnes *nz = count; *nzalloc = N; *mem = N*sizeof(Scalar); 31fa9ec3f1SLois Curfman McInnes return 0; 32289bc588SBarry Smith } 33289bc588SBarry Smith 34289bc588SBarry Smith /* ---------------------------------------------------------------*/ 35289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots, 36289bc588SBarry Smith rather than put it in the Mat->col slot.*/ 3744a69424SLois Curfman McInnes static int MatLUFactor_Dense(Mat matin,IS row,IS col) 38289bc588SBarry Smith { 3944a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 406abc6512SBarry Smith int info; 41289bc588SBarry Smith if (!mat->pivots) { 42289bc588SBarry Smith mat->pivots = (int *) MALLOC( mat->m*sizeof(int) ); 43289bc588SBarry Smith CHKPTR(mat->pivots); 44289bc588SBarry Smith } 45289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 46289bc588SBarry Smith if (info) SETERR(1,"Bad LU factorization"); 47289bc588SBarry Smith matin->factor = FACTOR_LU; 48289bc588SBarry Smith return 0; 49289bc588SBarry Smith } 5044a69424SLois Curfman McInnes static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,Mat *fact) 51289bc588SBarry Smith { 52289bc588SBarry Smith int ierr; 536abc6512SBarry Smith if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0); 54289bc588SBarry Smith return 0; 55289bc588SBarry Smith } 5644a69424SLois Curfman McInnes static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact) 57289bc588SBarry Smith { 5820563c6bSBarry Smith return MatLUFactor(*fact,0,0); 59289bc588SBarry Smith } 6046fc4a8fSLois Curfman McInnes static int MatCholeskyFactorSymbolic_Dense(Mat matin,IS row,Mat *fact) 61289bc588SBarry Smith { 62289bc588SBarry Smith int ierr; 636abc6512SBarry Smith if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0); 64289bc588SBarry Smith return 0; 65289bc588SBarry Smith } 6646fc4a8fSLois Curfman McInnes static int MatCholeskyFactorNumeric_Dense(Mat matin,Mat *fact) 67289bc588SBarry Smith { 6820563c6bSBarry Smith return MatCholeskyFactor(*fact,0); 69289bc588SBarry Smith } 7046fc4a8fSLois Curfman McInnes static int MatCholeskyFactor_Dense(Mat matin,IS perm) 71289bc588SBarry Smith { 7244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 736abc6512SBarry Smith int info; 74289bc588SBarry Smith if (mat->pivots) {FREE(mat->pivots); mat->pivots = 0;} 75289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 76289bc588SBarry Smith if (info) SETERR(1,"Bad Cholesky factorization"); 77289bc588SBarry Smith matin->factor = FACTOR_CHOLESKY; 78289bc588SBarry Smith return 0; 79289bc588SBarry Smith } 80289bc588SBarry Smith 8144a69424SLois Curfman McInnes static int MatSolve_Dense(Mat matin,Vec xx,Vec yy) 82289bc588SBarry Smith { 8344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 846abc6512SBarry Smith int one = 1, info; 856abc6512SBarry Smith Scalar *x, *y; 86289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 87289bc588SBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 889e25ed09SBarry Smith if (matin->factor == FACTOR_LU) { 89289bc588SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 90289bc588SBarry Smith y, &mat->m, &info ); 91289bc588SBarry Smith } 929e25ed09SBarry Smith else if (matin->factor == FACTOR_CHOLESKY){ 93289bc588SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 94289bc588SBarry Smith y, &mat->m, &info ); 95289bc588SBarry Smith } 969e25ed09SBarry Smith else SETERR(1,"Matrix must be factored to solve"); 97289bc588SBarry Smith if (info) SETERR(1,"Bad solve"); 98289bc588SBarry Smith return 0; 99289bc588SBarry Smith } 10044a69424SLois Curfman McInnes static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy) 101da3a660dSBarry Smith { 10244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1036abc6512SBarry Smith int one = 1, info; 1046abc6512SBarry Smith Scalar *x, *y; 105da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 106da3a660dSBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 107da3a660dSBarry Smith /* assume if pivots exist then LU else Cholesky */ 108da3a660dSBarry Smith if (mat->pivots) { 109da3a660dSBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 110da3a660dSBarry Smith y, &mat->m, &info ); 111da3a660dSBarry Smith } 112da3a660dSBarry Smith else { 113da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 114da3a660dSBarry Smith y, &mat->m, &info ); 115da3a660dSBarry Smith } 116da3a660dSBarry Smith if (info) SETERR(1,"Bad solve"); 117da3a660dSBarry Smith return 0; 118da3a660dSBarry Smith } 11944a69424SLois Curfman McInnes static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 120da3a660dSBarry Smith { 12144a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1226abc6512SBarry Smith int one = 1, info,ierr; 1236abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 124da3a660dSBarry Smith Vec tmp = 0; 125da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 126da3a660dSBarry Smith if (yy == zz) { 127da3a660dSBarry Smith ierr = VecCreate(yy,&tmp); CHKERR(ierr); 128da3a660dSBarry Smith ierr = VecCopy(yy,tmp); CHKERR(ierr); 129da3a660dSBarry Smith } 130da3a660dSBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 131da3a660dSBarry Smith /* assume if pivots exist then LU else Cholesky */ 132da3a660dSBarry Smith if (mat->pivots) { 133da3a660dSBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 134da3a660dSBarry Smith y, &mat->m, &info ); 135da3a660dSBarry Smith } 136da3a660dSBarry Smith else { 137da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 138da3a660dSBarry Smith y, &mat->m, &info ); 139da3a660dSBarry Smith } 140da3a660dSBarry Smith if (info) SETERR(1,"Bad solve"); 141da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 142da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 143da3a660dSBarry Smith return 0; 144da3a660dSBarry Smith } 14544a69424SLois Curfman McInnes static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy) 146da3a660dSBarry Smith { 14744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1486abc6512SBarry Smith int one = 1, info,ierr; 1496abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 150da3a660dSBarry Smith Vec tmp; 151da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 152da3a660dSBarry Smith if (yy == zz) { 153da3a660dSBarry Smith ierr = VecCreate(yy,&tmp); CHKERR(ierr); 154da3a660dSBarry Smith ierr = VecCopy(yy,tmp); CHKERR(ierr); 155da3a660dSBarry Smith } 156da3a660dSBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 157da3a660dSBarry Smith /* assume if pivots exist then LU else Cholesky */ 158da3a660dSBarry Smith if (mat->pivots) { 159da3a660dSBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 160da3a660dSBarry Smith y, &mat->m, &info ); 161da3a660dSBarry Smith } 162da3a660dSBarry Smith else { 163da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 164da3a660dSBarry Smith y, &mat->m, &info ); 165da3a660dSBarry Smith } 166da3a660dSBarry Smith if (info) SETERR(1,"Bad solve"); 167da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 168da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 169da3a660dSBarry Smith return 0; 170da3a660dSBarry Smith } 171289bc588SBarry Smith /* ------------------------------------------------------------------*/ 17220e1a0d4SLois Curfman McInnes static int MatRelax_Dense(Mat matin,Vec bb,double omega,MatSORType flag, 17320e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 174289bc588SBarry Smith { 17544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 176289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 1776abc6512SBarry Smith int o = 1,ierr, m = mat->m, i; 178289bc588SBarry Smith 179289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 180289bc588SBarry Smith /* this is a hack fix, should have another version without 181289bc588SBarry Smith the second BLdot */ 1826abc6512SBarry Smith if ((ierr = VecSet(&zero,xx))) SETERR(ierr,0); 183289bc588SBarry Smith } 184289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 185289bc588SBarry Smith while (its--) { 186289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 187289bc588SBarry Smith for ( i=0; i<m; i++ ) { 188289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 189d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 190289bc588SBarry Smith } 191289bc588SBarry Smith } 192289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 193289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 194289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 195d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 196289bc588SBarry Smith } 197289bc588SBarry Smith } 198289bc588SBarry Smith } 199289bc588SBarry Smith return 0; 200289bc588SBarry Smith } 201289bc588SBarry Smith 202289bc588SBarry Smith /* -----------------------------------------------------------------*/ 20344a69424SLois Curfman McInnes static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy) 204289bc588SBarry Smith { 20544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 206289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 207289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 208289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 209289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 210289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 211289bc588SBarry Smith return 0; 212289bc588SBarry Smith } 21344a69424SLois Curfman McInnes static int MatMult_Dense(Mat matin,Vec xx,Vec yy) 214289bc588SBarry Smith { 21544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 216289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 217289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 218289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 219289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 220289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 221289bc588SBarry Smith return 0; 222289bc588SBarry Smith } 22344a69424SLois Curfman McInnes static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 224289bc588SBarry Smith { 22544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 226289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 2276abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 228289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 229289bc588SBarry Smith if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 230289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 231289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 232289bc588SBarry Smith return 0; 233289bc588SBarry Smith } 23444a69424SLois Curfman McInnes static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 235289bc588SBarry Smith { 23644a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 237289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 238289bc588SBarry Smith int _One=1; 2396abc6512SBarry Smith Scalar _DOne=1.0; 240289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 241289bc588SBarry Smith VecGetArray(zz,&z); 242289bc588SBarry Smith if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 243289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 244289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 245289bc588SBarry Smith return 0; 246289bc588SBarry Smith } 247289bc588SBarry Smith 248289bc588SBarry Smith /* -----------------------------------------------------------------*/ 24944a69424SLois Curfman McInnes static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols, 250289bc588SBarry Smith Scalar **vals) 251289bc588SBarry Smith { 25244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 253289bc588SBarry Smith Scalar *v; 254289bc588SBarry Smith int i; 255289bc588SBarry Smith *ncols = mat->n; 256289bc588SBarry Smith if (cols) { 257289bc588SBarry Smith *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols); 258289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) *cols[i] = i; 259289bc588SBarry Smith } 260289bc588SBarry Smith if (vals) { 261289bc588SBarry Smith *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals); 262289bc588SBarry Smith v = mat->v + row; 263289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 264289bc588SBarry Smith } 265289bc588SBarry Smith return 0; 266289bc588SBarry Smith } 26744a69424SLois Curfman McInnes static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols, 268289bc588SBarry Smith Scalar **vals) 269289bc588SBarry Smith { 270289bc588SBarry Smith if (cols) { FREE(*cols); } 271289bc588SBarry Smith if (vals) { FREE(*vals); } 272289bc588SBarry Smith return 0; 273289bc588SBarry Smith } 274289bc588SBarry Smith /* ----------------------------------------------------------------*/ 27544a69424SLois Curfman McInnes static int MatInsert_Dense(Mat matin,int m,int *indexm,int n, 276e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 277289bc588SBarry Smith { 27844a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 279289bc588SBarry Smith int i,j; 280d6dfbf8fSBarry Smith 281289bc588SBarry Smith if (!mat->roworiented) { 282*edae2e7dSBarry Smith if (addv == INSERTVALUES) { 283289bc588SBarry Smith for ( j=0; j<n; j++ ) { 284d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v += m; continue;} 285289bc588SBarry Smith for ( i=0; i<m; i++ ) { 286d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v++; continue;} 287289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 288289bc588SBarry Smith } 289289bc588SBarry Smith } 290289bc588SBarry Smith } 291289bc588SBarry Smith else { 292289bc588SBarry Smith for ( j=0; j<n; j++ ) { 293d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v += m; continue;} 294289bc588SBarry Smith for ( i=0; i<m; i++ ) { 295d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v++; continue;} 296289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 297289bc588SBarry Smith } 298289bc588SBarry Smith } 299289bc588SBarry Smith } 300e8d4e0b9SBarry Smith } 301e8d4e0b9SBarry Smith else { 302*edae2e7dSBarry Smith if (addv == INSERTVALUES) { 303e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 304d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v += n; continue;} 305e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 306d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v++; continue;} 307e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 308e8d4e0b9SBarry Smith } 309e8d4e0b9SBarry Smith } 310e8d4e0b9SBarry Smith } 311289bc588SBarry Smith else { 312289bc588SBarry Smith for ( i=0; i<m; i++ ) { 313d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v += n; continue;} 314289bc588SBarry Smith for ( j=0; j<n; j++ ) { 315d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v++; continue;} 316289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 317289bc588SBarry Smith } 318289bc588SBarry Smith } 319289bc588SBarry Smith } 320e8d4e0b9SBarry Smith } 321289bc588SBarry Smith return 0; 322289bc588SBarry Smith } 323e8d4e0b9SBarry Smith 324289bc588SBarry Smith /* -----------------------------------------------------------------*/ 32544a69424SLois Curfman McInnes static int MatCopy_Dense(Mat matin,Mat *newmat) 326289bc588SBarry Smith { 32744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 328289bc588SBarry Smith int ierr; 329289bc588SBarry Smith Mat newi; 33044a69424SLois Curfman McInnes Mat_Dense *l; 3316b5873e3SBarry Smith if ((ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi))) 3326b5873e3SBarry Smith SETERR(ierr,0); 33344a69424SLois Curfman McInnes l = (Mat_Dense *) newi->data; 334289bc588SBarry Smith MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 335289bc588SBarry Smith *newmat = newi; 336289bc588SBarry Smith return 0; 337289bc588SBarry Smith } 338da3a660dSBarry Smith #include "viewer.h" 339289bc588SBarry Smith 34044a69424SLois Curfman McInnes int MatView_Dense(PetscObject obj,Viewer ptr) 341289bc588SBarry Smith { 342289bc588SBarry Smith Mat matin = (Mat) obj; 34344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 344289bc588SBarry Smith Scalar *v; 345289bc588SBarry Smith int i,j; 34623242f5aSBarry Smith PetscObject vobj = (PetscObject) ptr; 347da3a660dSBarry Smith 34823242f5aSBarry Smith if (ptr == 0) { 34923242f5aSBarry Smith ptr = STDOUT_VIEWER; vobj = (PetscObject) ptr; 35023242f5aSBarry Smith } 35123242f5aSBarry Smith if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) { 3526f75cfa4SBarry Smith return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v); 353da3a660dSBarry Smith } 354da3a660dSBarry Smith else { 3554a0ce102SLois Curfman McInnes FILE *fd = ViewerFileGetPointer_Private(ptr); 356289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 357289bc588SBarry Smith v = mat->v + i; 358289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 359289bc588SBarry Smith #if defined(PETSC_COMPLEX) 3608e37dc09SBarry Smith fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 361289bc588SBarry Smith #else 3628e37dc09SBarry Smith fprintf(fd,"%6.4e ",*v); v += mat->m; 363289bc588SBarry Smith #endif 364289bc588SBarry Smith } 3658e37dc09SBarry Smith fprintf(fd,"\n"); 366289bc588SBarry Smith } 367da3a660dSBarry Smith } 368289bc588SBarry Smith return 0; 369289bc588SBarry Smith } 370289bc588SBarry Smith 371289bc588SBarry Smith 37244a69424SLois Curfman McInnes static int MatDestroy_Dense(PetscObject obj) 373289bc588SBarry Smith { 374289bc588SBarry Smith Mat mat = (Mat) obj; 37544a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) mat->data; 376a5a9c739SBarry Smith #if defined(PETSC_LOG) 377a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 378a5a9c739SBarry Smith #endif 379289bc588SBarry Smith if (l->pivots) FREE(l->pivots); 380289bc588SBarry Smith FREE(l); 381a5a9c739SBarry Smith PLogObjectDestroy(mat); 3829e25ed09SBarry Smith PETSCHEADERDESTROY(mat); 383289bc588SBarry Smith return 0; 384289bc588SBarry Smith } 385289bc588SBarry Smith 38644a69424SLois Curfman McInnes static int MatTrans_Dense(Mat matin) 387289bc588SBarry Smith { 38844a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 389289bc588SBarry Smith int k,j; 390289bc588SBarry Smith Scalar *v = mat->v, tmp; 391289bc588SBarry Smith if (mat->m != mat->n) { 392289bc588SBarry Smith SETERR(1,"Cannot transpose rectangular dense matrix"); 393289bc588SBarry Smith } 394289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 395289bc588SBarry Smith for ( k=0; k<j; k++ ) { 396289bc588SBarry Smith tmp = v[j + k*mat->n]; 397289bc588SBarry Smith v[j + k*mat->n] = v[k + j*mat->n]; 398289bc588SBarry Smith v[k + j*mat->n] = tmp; 399289bc588SBarry Smith } 400289bc588SBarry Smith } 401289bc588SBarry Smith return 0; 402289bc588SBarry Smith } 403289bc588SBarry Smith 40444a69424SLois Curfman McInnes static int MatEqual_Dense(Mat matin1,Mat matin2) 405289bc588SBarry Smith { 40644a69424SLois Curfman McInnes Mat_Dense *mat1 = (Mat_Dense *) matin1->data; 40744a69424SLois Curfman McInnes Mat_Dense *mat2 = (Mat_Dense *) matin2->data; 408289bc588SBarry Smith int i; 409289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 410289bc588SBarry Smith if (mat1->m != mat2->m) return 0; 411289bc588SBarry Smith if (mat1->n != mat2->n) return 0; 412289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 413289bc588SBarry Smith if (*v1 != *v2) return 0; 414289bc588SBarry Smith v1++; v2++; 415289bc588SBarry Smith } 416289bc588SBarry Smith return 1; 417289bc588SBarry Smith } 418289bc588SBarry Smith 41946fc4a8fSLois Curfman McInnes static int MatGetDiagonal_Dense(Mat matin,Vec v) 420289bc588SBarry Smith { 42144a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 4226abc6512SBarry Smith int i, n; 4236abc6512SBarry Smith Scalar *x; 424289bc588SBarry Smith CHKTYPE(v,SEQVECTOR); 425289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 426289bc588SBarry Smith if (n != mat->m) SETERR(1,"Nonconforming matrix and vector"); 427289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 428289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 429289bc588SBarry Smith } 430289bc588SBarry Smith return 0; 431289bc588SBarry Smith } 432289bc588SBarry Smith 43344a69424SLois Curfman McInnes static int MatScale_Dense(Mat matin,Vec ll,Vec rr) 434289bc588SBarry Smith { 43544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 436da3a660dSBarry Smith Scalar *l,*r,x,*v; 437da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 43828988994SBarry Smith if (ll) { 439da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 440da3a660dSBarry Smith if (m != mat->m) SETERR(1,"Left scaling vector wrong length"); 441da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 442da3a660dSBarry Smith x = l[i]; 443da3a660dSBarry Smith v = mat->v + i; 444da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 445da3a660dSBarry Smith } 446da3a660dSBarry Smith } 44728988994SBarry Smith if (rr) { 448da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 449da3a660dSBarry Smith if (n != mat->n) SETERR(1,"Right scaling vector wrong length"); 450da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 451da3a660dSBarry Smith x = r[i]; 452da3a660dSBarry Smith v = mat->v + i*m; 453da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 454da3a660dSBarry Smith } 455da3a660dSBarry Smith } 456289bc588SBarry Smith return 0; 457289bc588SBarry Smith } 458289bc588SBarry Smith 459da3a660dSBarry Smith 46044a69424SLois Curfman McInnes static int MatNorm_Dense(Mat matin,int type,double *norm) 461289bc588SBarry Smith { 46244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 463289bc588SBarry Smith Scalar *v = mat->v; 464289bc588SBarry Smith double sum = 0.0; 465289bc588SBarry Smith int i, j; 466289bc588SBarry Smith if (type == NORM_FROBENIUS) { 467289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 468289bc588SBarry Smith #if defined(PETSC_COMPLEX) 469289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 470289bc588SBarry Smith #else 471289bc588SBarry Smith sum += (*v)*(*v); v++; 472289bc588SBarry Smith #endif 473289bc588SBarry Smith } 474289bc588SBarry Smith *norm = sqrt(sum); 475289bc588SBarry Smith } 476289bc588SBarry Smith else if (type == NORM_1) { 477289bc588SBarry Smith *norm = 0.0; 478289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 479289bc588SBarry Smith sum = 0.0; 480289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 481289bc588SBarry Smith #if defined(PETSC_COMPLEX) 482289bc588SBarry Smith sum += abs(*v++); 483289bc588SBarry Smith #else 484289bc588SBarry Smith sum += fabs(*v++); 485289bc588SBarry Smith #endif 486289bc588SBarry Smith } 487289bc588SBarry Smith if (sum > *norm) *norm = sum; 488289bc588SBarry Smith } 489289bc588SBarry Smith } 490289bc588SBarry Smith else if (type == NORM_INFINITY) { 491289bc588SBarry Smith *norm = 0.0; 492289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 493289bc588SBarry Smith v = mat->v + j; 494289bc588SBarry Smith sum = 0.0; 495289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 496289bc588SBarry Smith #if defined(PETSC_COMPLEX) 497289bc588SBarry Smith sum += abs(*v); v += mat->m; 498289bc588SBarry Smith #else 499289bc588SBarry Smith sum += fabs(*v); v += mat->m; 500289bc588SBarry Smith #endif 501289bc588SBarry Smith } 502289bc588SBarry Smith if (sum > *norm) *norm = sum; 503289bc588SBarry Smith } 504289bc588SBarry Smith } 505289bc588SBarry Smith else { 506289bc588SBarry Smith SETERR(1,"No support for the two norm yet"); 507289bc588SBarry Smith } 508289bc588SBarry Smith return 0; 509289bc588SBarry Smith } 510289bc588SBarry Smith 51120e1a0d4SLois Curfman McInnes static int MatSetOption_Dense(Mat aijin,MatOption op) 512289bc588SBarry Smith { 51344a69424SLois Curfman McInnes Mat_Dense *aij = (Mat_Dense *) aijin->data; 514289bc588SBarry Smith if (op == ROW_ORIENTED) aij->roworiented = 1; 515289bc588SBarry Smith else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 516289bc588SBarry Smith /* doesn't care about sorted rows or columns */ 517289bc588SBarry Smith return 0; 518289bc588SBarry Smith } 519289bc588SBarry Smith 52044a69424SLois Curfman McInnes static int MatZero_Dense(Mat A) 5216f0a148fSBarry Smith { 52244a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) A->data; 5236f0a148fSBarry Smith MEMSET(l->v,0,l->m*l->n*sizeof(Scalar)); 5246f0a148fSBarry Smith return 0; 5256f0a148fSBarry Smith } 5266f0a148fSBarry Smith 52744a69424SLois Curfman McInnes static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag) 5286f0a148fSBarry Smith { 52944a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) A->data; 5306abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 5316f0a148fSBarry Smith Scalar *slot; 532260925b8SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERR(ierr); 533260925b8SBarry Smith ierr = ISGetIndices(is,&rows); CHKERR(ierr); 5346f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5356f0a148fSBarry Smith slot = l->v + rows[i]; 5366f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 5376f0a148fSBarry Smith } 5386f0a148fSBarry Smith if (diag) { 5396f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5406f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 541260925b8SBarry Smith *slot = *diag; 5426f0a148fSBarry Smith } 5436f0a148fSBarry Smith } 544260925b8SBarry Smith ISRestoreIndices(is,&rows); 5456f0a148fSBarry Smith return 0; 5466f0a148fSBarry Smith } 547557bce09SLois Curfman McInnes 548fa9ec3f1SLois Curfman McInnes static int MatGetSize_Dense(Mat matin,int *m,int *n) 549557bce09SLois Curfman McInnes { 55044a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 551557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 552557bce09SLois Curfman McInnes return 0; 553557bce09SLois Curfman McInnes } 554557bce09SLois Curfman McInnes 55544a69424SLois Curfman McInnes static int MatGetArray_Dense(Mat matin,Scalar **array) 55664e87e97SBarry Smith { 55744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 55864e87e97SBarry Smith *array = mat->v; 55964e87e97SBarry Smith return 0; 56064e87e97SBarry Smith } 5610754003eSLois Curfman McInnes 5620754003eSLois Curfman McInnes 5630754003eSLois Curfman McInnes static int MatGetSubMatrixInPlace_Dense(Mat matin,IS isrow,IS iscol) 5640754003eSLois Curfman McInnes { 5650754003eSLois Curfman McInnes SETERR(1,"MatGetSubMatrixInPlace_Dense not yet done."); 5660754003eSLois Curfman McInnes } 5670754003eSLois Curfman McInnes 5680754003eSLois Curfman McInnes static int MatGetSubMatrix_Dense(Mat matin,IS isrow,IS iscol,Mat *submat) 5690754003eSLois Curfman McInnes { 5700754003eSLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 5710754003eSLois Curfman McInnes int nznew, *smap, i, j, ierr, oldcols = mat->n; 5720754003eSLois Curfman McInnes int *irow, *icol, nrows, ncols, *cwork, jstart; 5730754003eSLois Curfman McInnes Scalar *vwork, *val; 5740754003eSLois Curfman McInnes Mat newmat; 5750754003eSLois Curfman McInnes 5760754003eSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow); CHKERR(ierr); 5770754003eSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol); CHKERR(ierr); 5780754003eSLois Curfman McInnes ierr = ISGetSize(isrow,&nrows); CHKERR(ierr); 5790754003eSLois Curfman McInnes ierr = ISGetSize(iscol,&ncols); CHKERR(ierr); 5800754003eSLois Curfman McInnes 5810754003eSLois Curfman McInnes smap = (int *) MALLOC(oldcols*sizeof(int)); CHKPTR(smap); 5820754003eSLois Curfman McInnes cwork = (int *) MALLOC(ncols*sizeof(int)); CHKPTR(cwork); 5830754003eSLois Curfman McInnes vwork = (Scalar *) MALLOC(ncols*sizeof(Scalar)); CHKPTR(vwork); 5840754003eSLois Curfman McInnes memset((char*)smap,0,oldcols*sizeof(int)); 5850754003eSLois Curfman McInnes for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 5860754003eSLois Curfman McInnes 5870754003eSLois Curfman McInnes /* Create and fill new matrix */ 5880754003eSLois Curfman McInnes ierr = MatCreateSequentialDense(matin->comm,nrows,ncols,&newmat); 5890754003eSLois Curfman McInnes CHKERR(ierr); 5900754003eSLois Curfman McInnes for (i=0; i<nrows; i++) { 5910754003eSLois Curfman McInnes nznew = 0; 5920754003eSLois Curfman McInnes val = mat->v + irow[i]; 5930754003eSLois Curfman McInnes for (j=0; j<oldcols; j++) { 5940754003eSLois Curfman McInnes if (smap[j]) { 5950754003eSLois Curfman McInnes cwork[nznew] = smap[j] - 1; 5960754003eSLois Curfman McInnes vwork[nznew++] = val[j * mat->m]; 5970754003eSLois Curfman McInnes } 5980754003eSLois Curfman McInnes } 599*edae2e7dSBarry Smith ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES); 6000754003eSLois Curfman McInnes CHKERR(ierr); 6010754003eSLois Curfman McInnes } 602ee50ffe9SBarry Smith ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERR(ierr); 603ee50ffe9SBarry Smith ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERR(ierr); 6040754003eSLois Curfman McInnes 6050754003eSLois Curfman McInnes /* Free work space */ 6060754003eSLois Curfman McInnes FREE(smap); FREE(cwork); FREE(vwork); 6070754003eSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow); CHKERR(ierr); 6080754003eSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol); CHKERR(ierr); 6090754003eSLois Curfman McInnes *submat = newmat; 6100754003eSLois Curfman McInnes return 0; 6110754003eSLois Curfman McInnes } 6120754003eSLois Curfman McInnes 613289bc588SBarry Smith /* -------------------------------------------------------------------*/ 61444a69424SLois Curfman McInnes static struct _MatOps MatOps = {MatInsert_Dense, 61544a69424SLois Curfman McInnes MatGetRow_Dense, MatRestoreRow_Dense, 61644a69424SLois Curfman McInnes MatMult_Dense, MatMultAdd_Dense, 61744a69424SLois Curfman McInnes MatMultTrans_Dense, MatMultTransAdd_Dense, 61844a69424SLois Curfman McInnes MatSolve_Dense,MatSolveAdd_Dense, 61944a69424SLois Curfman McInnes MatSolveTrans_Dense,MatSolveTransAdd_Dense, 62046fc4a8fSLois Curfman McInnes MatLUFactor_Dense,MatCholeskyFactor_Dense, 62144a69424SLois Curfman McInnes MatRelax_Dense, 62244a69424SLois Curfman McInnes MatTrans_Dense, 623fa9ec3f1SLois Curfman McInnes MatGetInfo_Dense,MatEqual_Dense, 62444a69424SLois Curfman McInnes MatCopy_Dense, 62546fc4a8fSLois Curfman McInnes MatGetDiagonal_Dense,MatScale_Dense,MatNorm_Dense, 626289bc588SBarry Smith 0,0, 627681d1451SLois Curfman McInnes 0, MatSetOption_Dense,MatZero_Dense,MatZeroRows_Dense,0, 62844a69424SLois Curfman McInnes MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense, 62946fc4a8fSLois Curfman McInnes MatCholeskyFactorSymbolic_Dense,MatCholeskyFactorNumeric_Dense, 630fa9ec3f1SLois Curfman McInnes MatGetSize_Dense,MatGetSize_Dense,0, 6310754003eSLois Curfman McInnes 0,0,MatGetArray_Dense,0, 6320754003eSLois Curfman McInnes MatGetSubMatrix_Dense,MatGetSubMatrixInPlace_Dense}; 6330754003eSLois Curfman McInnes 63420563c6bSBarry Smith /*@ 63520563c6bSBarry Smith MatCreateSequentialDense - Creates a sequential dense matrix that 636d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 637d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 638289bc588SBarry Smith 63920563c6bSBarry Smith Input Parameters: 6400c775827SLois Curfman McInnes . comm - MPI communicator, set to MPI_COMM_SELF 6410c775827SLois Curfman McInnes . m - number of rows 6420c775827SLois Curfman McInnes . n - number of column 64320563c6bSBarry Smith 64420563c6bSBarry Smith Output Parameter: 6450c775827SLois Curfman McInnes . newmat - the matrix 64620563c6bSBarry Smith 647d65003e9SLois Curfman McInnes .keywords: Mat, dense, matrix, LAPACK, BLAS 648d65003e9SLois Curfman McInnes 649d65003e9SLois Curfman McInnes .seealso: MatCreateSequentialAIJ() 65020563c6bSBarry Smith @*/ 6516b5873e3SBarry Smith int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat) 652289bc588SBarry Smith { 65344a69424SLois Curfman McInnes int size = sizeof(Mat_Dense) + m*n*sizeof(Scalar); 654289bc588SBarry Smith Mat mat; 65544a69424SLois Curfman McInnes Mat_Dense *l; 656289bc588SBarry Smith *newmat = 0; 6576b5873e3SBarry Smith PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm); 658a5a9c739SBarry Smith PLogObjectCreate(mat); 65944a69424SLois Curfman McInnes l = (Mat_Dense *) MALLOC(size); CHKPTR(l); 660289bc588SBarry Smith mat->ops = &MatOps; 66144a69424SLois Curfman McInnes mat->destroy = MatDestroy_Dense; 66244a69424SLois Curfman McInnes mat->view = MatView_Dense; 663289bc588SBarry Smith mat->data = (void *) l; 664289bc588SBarry Smith mat->factor = 0; 665d6dfbf8fSBarry Smith 666289bc588SBarry Smith l->m = m; 667289bc588SBarry Smith l->n = n; 668289bc588SBarry Smith l->v = (Scalar *) (l + 1); 669289bc588SBarry Smith l->pivots = 0; 670289bc588SBarry Smith l->roworiented = 1; 671d6dfbf8fSBarry Smith 672289bc588SBarry Smith MEMSET(l->v,0,m*n*sizeof(Scalar)); 673289bc588SBarry Smith *newmat = mat; 674289bc588SBarry Smith return 0; 675289bc588SBarry Smith } 676289bc588SBarry Smith 67744a69424SLois Curfman McInnes int MatCreate_Dense(Mat matin,Mat *newmat) 678289bc588SBarry Smith { 67944a69424SLois Curfman McInnes Mat_Dense *m = (Mat_Dense *) matin->data; 6806b5873e3SBarry Smith return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat); 681289bc588SBarry Smith } 682