1cb512458SBarry Smith #ifndef lint 2*6b5873e3SBarry Smith static char vcid[] = "$Id: dense.c,v 1.20 1995/03/29 20:56:54 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" 13289bc588SBarry Smith 14289bc588SBarry Smith typedef struct { 15289bc588SBarry Smith Scalar *v; 16289bc588SBarry Smith int roworiented; 17289bc588SBarry Smith int m,n,pad; 18289bc588SBarry Smith int *pivots; /* pivots in LU factorization */ 1944a69424SLois Curfman McInnes } Mat_Dense; 20289bc588SBarry Smith 21289bc588SBarry Smith 22fa9ec3f1SLois Curfman McInnes static int MatGetInfo_Dense(Mat matin,int flag,int *nz,int *nzalloc,int *mem) 23289bc588SBarry Smith { 2444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 25289bc588SBarry Smith int i,N = mat->m*mat->n,count = 0; 26289bc588SBarry Smith Scalar *v = mat->v; 27289bc588SBarry Smith for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 28fa9ec3f1SLois Curfman McInnes *nz = count; *nzalloc = N; *mem = N*sizeof(Scalar); 29fa9ec3f1SLois Curfman McInnes return 0; 30289bc588SBarry Smith } 31289bc588SBarry Smith 32289bc588SBarry Smith /* ---------------------------------------------------------------*/ 33289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots, 34289bc588SBarry Smith rather than put it in the Mat->col slot.*/ 3544a69424SLois Curfman McInnes static int MatLUFactor_Dense(Mat matin,IS row,IS col) 36289bc588SBarry Smith { 3744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 386abc6512SBarry Smith int info; 39289bc588SBarry Smith if (!mat->pivots) { 40289bc588SBarry Smith mat->pivots = (int *) MALLOC( mat->m*sizeof(int) ); 41289bc588SBarry Smith CHKPTR(mat->pivots); 42289bc588SBarry Smith } 43289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 44289bc588SBarry Smith if (info) SETERR(1,"Bad LU factorization"); 45289bc588SBarry Smith matin->factor = FACTOR_LU; 46289bc588SBarry Smith return 0; 47289bc588SBarry Smith } 4844a69424SLois Curfman McInnes static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,Mat *fact) 49289bc588SBarry Smith { 50289bc588SBarry Smith int ierr; 516abc6512SBarry Smith if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0); 52289bc588SBarry Smith return 0; 53289bc588SBarry Smith } 5444a69424SLois Curfman McInnes static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact) 55289bc588SBarry Smith { 5620563c6bSBarry Smith return MatLUFactor(*fact,0,0); 57289bc588SBarry Smith } 5844a69424SLois Curfman McInnes static int MatChFactorSymbolic_Dense(Mat matin,IS row,Mat *fact) 59289bc588SBarry Smith { 60289bc588SBarry Smith int ierr; 616abc6512SBarry Smith if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0); 62289bc588SBarry Smith return 0; 63289bc588SBarry Smith } 6444a69424SLois Curfman McInnes static int MatChFactorNumeric_Dense(Mat matin,Mat *fact) 65289bc588SBarry Smith { 6620563c6bSBarry Smith return MatCholeskyFactor(*fact,0); 67289bc588SBarry Smith } 6844a69424SLois Curfman McInnes static int MatChFactor_Dense(Mat matin,IS perm) 69289bc588SBarry Smith { 7044a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 716abc6512SBarry Smith int info; 72289bc588SBarry Smith if (mat->pivots) {FREE(mat->pivots); mat->pivots = 0;} 73289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 74289bc588SBarry Smith if (info) SETERR(1,"Bad Cholesky factorization"); 75289bc588SBarry Smith matin->factor = FACTOR_CHOLESKY; 76289bc588SBarry Smith return 0; 77289bc588SBarry Smith } 78289bc588SBarry Smith 7944a69424SLois Curfman McInnes static int MatSolve_Dense(Mat matin,Vec xx,Vec yy) 80289bc588SBarry Smith { 8144a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 826abc6512SBarry Smith int one = 1, info; 836abc6512SBarry Smith Scalar *x, *y; 84289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 85289bc588SBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 869e25ed09SBarry Smith if (matin->factor == FACTOR_LU) { 87289bc588SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 88289bc588SBarry Smith y, &mat->m, &info ); 89289bc588SBarry Smith } 909e25ed09SBarry Smith else if (matin->factor == FACTOR_CHOLESKY){ 91289bc588SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 92289bc588SBarry Smith y, &mat->m, &info ); 93289bc588SBarry Smith } 949e25ed09SBarry Smith else SETERR(1,"Matrix must be factored to solve"); 95289bc588SBarry Smith if (info) SETERR(1,"Bad solve"); 96289bc588SBarry Smith return 0; 97289bc588SBarry Smith } 9844a69424SLois Curfman McInnes static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy) 99da3a660dSBarry Smith { 10044a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1016abc6512SBarry Smith int one = 1, info; 1026abc6512SBarry Smith Scalar *x, *y; 103da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 104da3a660dSBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 105da3a660dSBarry Smith /* assume if pivots exist then LU else Cholesky */ 106da3a660dSBarry Smith if (mat->pivots) { 107da3a660dSBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 108da3a660dSBarry Smith y, &mat->m, &info ); 109da3a660dSBarry Smith } 110da3a660dSBarry Smith else { 111da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 112da3a660dSBarry Smith y, &mat->m, &info ); 113da3a660dSBarry Smith } 114da3a660dSBarry Smith if (info) SETERR(1,"Bad solve"); 115da3a660dSBarry Smith return 0; 116da3a660dSBarry Smith } 11744a69424SLois Curfman McInnes static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 118da3a660dSBarry Smith { 11944a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1206abc6512SBarry Smith int one = 1, info,ierr; 1216abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 122da3a660dSBarry Smith Vec tmp = 0; 123da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 124da3a660dSBarry Smith if (yy == zz) { 125da3a660dSBarry Smith ierr = VecCreate(yy,&tmp); CHKERR(ierr); 126da3a660dSBarry Smith ierr = VecCopy(yy,tmp); CHKERR(ierr); 127da3a660dSBarry Smith } 128da3a660dSBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 129da3a660dSBarry Smith /* assume if pivots exist then LU else Cholesky */ 130da3a660dSBarry Smith if (mat->pivots) { 131da3a660dSBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 132da3a660dSBarry Smith y, &mat->m, &info ); 133da3a660dSBarry Smith } 134da3a660dSBarry Smith else { 135da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 136da3a660dSBarry Smith y, &mat->m, &info ); 137da3a660dSBarry Smith } 138da3a660dSBarry Smith if (info) SETERR(1,"Bad solve"); 139da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 140da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 141da3a660dSBarry Smith return 0; 142da3a660dSBarry Smith } 14344a69424SLois Curfman McInnes static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy) 144da3a660dSBarry Smith { 14544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1466abc6512SBarry Smith int one = 1, info,ierr; 1476abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 148da3a660dSBarry Smith Vec tmp; 149da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 150da3a660dSBarry Smith if (yy == zz) { 151da3a660dSBarry Smith ierr = VecCreate(yy,&tmp); CHKERR(ierr); 152da3a660dSBarry Smith ierr = VecCopy(yy,tmp); CHKERR(ierr); 153da3a660dSBarry Smith } 154da3a660dSBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 155da3a660dSBarry Smith /* assume if pivots exist then LU else Cholesky */ 156da3a660dSBarry Smith if (mat->pivots) { 157da3a660dSBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 158da3a660dSBarry Smith y, &mat->m, &info ); 159da3a660dSBarry Smith } 160da3a660dSBarry Smith else { 161da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 162da3a660dSBarry Smith y, &mat->m, &info ); 163da3a660dSBarry Smith } 164da3a660dSBarry Smith if (info) SETERR(1,"Bad solve"); 165da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 166da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 167da3a660dSBarry Smith return 0; 168da3a660dSBarry Smith } 169289bc588SBarry Smith /* ------------------------------------------------------------------*/ 17044a69424SLois Curfman McInnes static int MatRelax_Dense(Mat matin,Vec bb,double omega,int flag,double shift, 171289bc588SBarry Smith int its,Vec xx) 172289bc588SBarry Smith { 17344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 174289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 1756abc6512SBarry Smith int o = 1,ierr, m = mat->m, i; 176289bc588SBarry Smith 177289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 178289bc588SBarry Smith /* this is a hack fix, should have another version without 179289bc588SBarry Smith the second BLdot */ 1806abc6512SBarry Smith if ((ierr = VecSet(&zero,xx))) SETERR(ierr,0); 181289bc588SBarry Smith } 182289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 183289bc588SBarry Smith while (its--) { 184289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 185289bc588SBarry Smith for ( i=0; i<m; i++ ) { 186289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 187d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 188289bc588SBarry Smith } 189289bc588SBarry Smith } 190289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 191289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 192289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 193d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 194289bc588SBarry Smith } 195289bc588SBarry Smith } 196289bc588SBarry Smith } 197289bc588SBarry Smith return 0; 198289bc588SBarry Smith } 199289bc588SBarry Smith 200289bc588SBarry Smith /* -----------------------------------------------------------------*/ 20144a69424SLois Curfman McInnes static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy) 202289bc588SBarry Smith { 20344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 204289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 205289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 206289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 207289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 208289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 209289bc588SBarry Smith return 0; 210289bc588SBarry Smith } 21144a69424SLois Curfman McInnes static int MatMult_Dense(Mat matin,Vec xx,Vec yy) 212289bc588SBarry Smith { 21344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 214289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 215289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 216289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 217289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 218289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 219289bc588SBarry Smith return 0; 220289bc588SBarry Smith } 22144a69424SLois Curfman McInnes static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 222289bc588SBarry Smith { 22344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 224289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 2256abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 226289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 227289bc588SBarry Smith if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 228289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 229289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 230289bc588SBarry Smith return 0; 231289bc588SBarry Smith } 23244a69424SLois Curfman McInnes static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 233289bc588SBarry Smith { 23444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 235289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 236289bc588SBarry Smith int _One=1; 2376abc6512SBarry Smith Scalar _DOne=1.0; 238289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 239289bc588SBarry Smith VecGetArray(zz,&z); 240289bc588SBarry Smith if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 241289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 242289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 243289bc588SBarry Smith return 0; 244289bc588SBarry Smith } 245289bc588SBarry Smith 246289bc588SBarry Smith /* -----------------------------------------------------------------*/ 24744a69424SLois Curfman McInnes static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols, 248289bc588SBarry Smith Scalar **vals) 249289bc588SBarry Smith { 25044a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 251289bc588SBarry Smith Scalar *v; 252289bc588SBarry Smith int i; 253289bc588SBarry Smith *ncols = mat->n; 254289bc588SBarry Smith if (cols) { 255289bc588SBarry Smith *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols); 256289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) *cols[i] = i; 257289bc588SBarry Smith } 258289bc588SBarry Smith if (vals) { 259289bc588SBarry Smith *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals); 260289bc588SBarry Smith v = mat->v + row; 261289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 262289bc588SBarry Smith } 263289bc588SBarry Smith return 0; 264289bc588SBarry Smith } 26544a69424SLois Curfman McInnes static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols, 266289bc588SBarry Smith Scalar **vals) 267289bc588SBarry Smith { 268289bc588SBarry Smith if (cols) { FREE(*cols); } 269289bc588SBarry Smith if (vals) { FREE(*vals); } 270289bc588SBarry Smith return 0; 271289bc588SBarry Smith } 272289bc588SBarry Smith /* ----------------------------------------------------------------*/ 27344a69424SLois Curfman McInnes static int MatInsert_Dense(Mat matin,int m,int *indexm,int n, 274e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 275289bc588SBarry Smith { 27644a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 277289bc588SBarry Smith int i,j; 278d6dfbf8fSBarry Smith 279289bc588SBarry Smith if (!mat->roworiented) { 280e8d4e0b9SBarry Smith if (addv == InsertValues) { 281289bc588SBarry Smith for ( j=0; j<n; j++ ) { 282d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v += m; continue;} 283289bc588SBarry Smith for ( i=0; i<m; i++ ) { 284d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v++; continue;} 285289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 286289bc588SBarry Smith } 287289bc588SBarry Smith } 288289bc588SBarry Smith } 289289bc588SBarry Smith else { 290289bc588SBarry Smith for ( j=0; j<n; j++ ) { 291d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v += m; continue;} 292289bc588SBarry Smith for ( i=0; i<m; i++ ) { 293d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v++; continue;} 294289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 295289bc588SBarry Smith } 296289bc588SBarry Smith } 297289bc588SBarry Smith } 298e8d4e0b9SBarry Smith } 299e8d4e0b9SBarry Smith else { 300e8d4e0b9SBarry Smith if (addv == InsertValues) { 301e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 302d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v += n; continue;} 303e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 304d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v++; continue;} 305e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 306e8d4e0b9SBarry Smith } 307e8d4e0b9SBarry Smith } 308e8d4e0b9SBarry Smith } 309289bc588SBarry Smith else { 310289bc588SBarry Smith for ( i=0; i<m; i++ ) { 311d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v += n; continue;} 312289bc588SBarry Smith for ( j=0; j<n; j++ ) { 313d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v++; continue;} 314289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 315289bc588SBarry Smith } 316289bc588SBarry Smith } 317289bc588SBarry Smith } 318e8d4e0b9SBarry Smith } 319289bc588SBarry Smith return 0; 320289bc588SBarry Smith } 321e8d4e0b9SBarry Smith 322289bc588SBarry Smith /* -----------------------------------------------------------------*/ 32344a69424SLois Curfman McInnes static int MatCopy_Dense(Mat matin,Mat *newmat) 324289bc588SBarry Smith { 32544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 326289bc588SBarry Smith int ierr; 327289bc588SBarry Smith Mat newi; 32844a69424SLois Curfman McInnes Mat_Dense *l; 329*6b5873e3SBarry Smith if ((ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi))) 330*6b5873e3SBarry Smith SETERR(ierr,0); 33144a69424SLois Curfman McInnes l = (Mat_Dense *) newi->data; 332289bc588SBarry Smith MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 333289bc588SBarry Smith *newmat = newi; 334289bc588SBarry Smith return 0; 335289bc588SBarry Smith } 336da3a660dSBarry Smith #include "viewer.h" 337289bc588SBarry Smith 33844a69424SLois Curfman McInnes int MatView_Dense(PetscObject obj,Viewer ptr) 339289bc588SBarry Smith { 340289bc588SBarry Smith Mat matin = (Mat) obj; 34144a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 342289bc588SBarry Smith Scalar *v; 343289bc588SBarry Smith int i,j; 344da3a660dSBarry Smith PetscObject ojb = (PetscObject) ptr; 345da3a660dSBarry Smith 3466abc6512SBarry Smith if (ojb && ojb->cookie == VIEWER_COOKIE && ojb->type == MATLAB_VIEWER) { 347da3a660dSBarry Smith return ViewerMatlabPutArray(ptr,mat->m,mat->n,mat->v); 348da3a660dSBarry Smith } 349da3a660dSBarry Smith else { 3508e37dc09SBarry Smith FILE *fd = ViewerFileGetPointer(ptr); 351289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 352289bc588SBarry Smith v = mat->v + i; 353289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 354289bc588SBarry Smith #if defined(PETSC_COMPLEX) 3558e37dc09SBarry Smith fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 356289bc588SBarry Smith #else 3578e37dc09SBarry Smith fprintf(fd,"%6.4e ",*v); v += mat->m; 358289bc588SBarry Smith #endif 359289bc588SBarry Smith } 3608e37dc09SBarry Smith fprintf(fd,"\n"); 361289bc588SBarry Smith } 362da3a660dSBarry Smith } 363289bc588SBarry Smith return 0; 364289bc588SBarry Smith } 365289bc588SBarry Smith 366289bc588SBarry Smith 36744a69424SLois Curfman McInnes static int MatDestroy_Dense(PetscObject obj) 368289bc588SBarry Smith { 369289bc588SBarry Smith Mat mat = (Mat) obj; 37044a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) mat->data; 371a5a9c739SBarry Smith #if defined(PETSC_LOG) 372a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 373a5a9c739SBarry Smith #endif 374289bc588SBarry Smith if (l->pivots) FREE(l->pivots); 375289bc588SBarry Smith FREE(l); 376a5a9c739SBarry Smith PLogObjectDestroy(mat); 3779e25ed09SBarry Smith PETSCHEADERDESTROY(mat); 378289bc588SBarry Smith return 0; 379289bc588SBarry Smith } 380289bc588SBarry Smith 38144a69424SLois Curfman McInnes static int MatTrans_Dense(Mat matin) 382289bc588SBarry Smith { 38344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 384289bc588SBarry Smith int k,j; 385289bc588SBarry Smith Scalar *v = mat->v, tmp; 386289bc588SBarry Smith if (mat->m != mat->n) { 387289bc588SBarry Smith SETERR(1,"Cannot transpose rectangular dense matrix"); 388289bc588SBarry Smith } 389289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 390289bc588SBarry Smith for ( k=0; k<j; k++ ) { 391289bc588SBarry Smith tmp = v[j + k*mat->n]; 392289bc588SBarry Smith v[j + k*mat->n] = v[k + j*mat->n]; 393289bc588SBarry Smith v[k + j*mat->n] = tmp; 394289bc588SBarry Smith } 395289bc588SBarry Smith } 396289bc588SBarry Smith return 0; 397289bc588SBarry Smith } 398289bc588SBarry Smith 39944a69424SLois Curfman McInnes static int MatEqual_Dense(Mat matin1,Mat matin2) 400289bc588SBarry Smith { 40144a69424SLois Curfman McInnes Mat_Dense *mat1 = (Mat_Dense *) matin1->data; 40244a69424SLois Curfman McInnes Mat_Dense *mat2 = (Mat_Dense *) matin2->data; 403289bc588SBarry Smith int i; 404289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 405289bc588SBarry Smith if (mat1->m != mat2->m) return 0; 406289bc588SBarry Smith if (mat1->n != mat2->n) return 0; 407289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 408289bc588SBarry Smith if (*v1 != *v2) return 0; 409289bc588SBarry Smith v1++; v2++; 410289bc588SBarry Smith } 411289bc588SBarry Smith return 1; 412289bc588SBarry Smith } 413289bc588SBarry Smith 41444a69424SLois Curfman McInnes static int MatGetDiag_Dense(Mat matin,Vec v) 415289bc588SBarry Smith { 41644a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 4176abc6512SBarry Smith int i, n; 4186abc6512SBarry Smith Scalar *x; 419289bc588SBarry Smith CHKTYPE(v,SEQVECTOR); 420289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 421289bc588SBarry Smith if (n != mat->m) SETERR(1,"Nonconforming matrix and vector"); 422289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 423289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 424289bc588SBarry Smith } 425289bc588SBarry Smith return 0; 426289bc588SBarry Smith } 427289bc588SBarry Smith 42844a69424SLois Curfman McInnes static int MatScale_Dense(Mat matin,Vec ll,Vec rr) 429289bc588SBarry Smith { 43044a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 431da3a660dSBarry Smith Scalar *l,*r,x,*v; 432da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 43328988994SBarry Smith if (ll) { 434da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 435da3a660dSBarry Smith if (m != mat->m) SETERR(1,"Left scaling vector wrong length"); 436da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 437da3a660dSBarry Smith x = l[i]; 438da3a660dSBarry Smith v = mat->v + i; 439da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 440da3a660dSBarry Smith } 441da3a660dSBarry Smith } 44228988994SBarry Smith if (rr) { 443da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 444da3a660dSBarry Smith if (n != mat->n) SETERR(1,"Right scaling vector wrong length"); 445da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 446da3a660dSBarry Smith x = r[i]; 447da3a660dSBarry Smith v = mat->v + i*m; 448da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 449da3a660dSBarry Smith } 450da3a660dSBarry Smith } 451289bc588SBarry Smith return 0; 452289bc588SBarry Smith } 453289bc588SBarry Smith 454da3a660dSBarry Smith 45544a69424SLois Curfman McInnes static int MatNorm_Dense(Mat matin,int type,double *norm) 456289bc588SBarry Smith { 45744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 458289bc588SBarry Smith Scalar *v = mat->v; 459289bc588SBarry Smith double sum = 0.0; 460289bc588SBarry Smith int i, j; 461289bc588SBarry Smith if (type == NORM_FROBENIUS) { 462289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 463289bc588SBarry Smith #if defined(PETSC_COMPLEX) 464289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 465289bc588SBarry Smith #else 466289bc588SBarry Smith sum += (*v)*(*v); v++; 467289bc588SBarry Smith #endif 468289bc588SBarry Smith } 469289bc588SBarry Smith *norm = sqrt(sum); 470289bc588SBarry Smith } 471289bc588SBarry Smith else if (type == NORM_1) { 472289bc588SBarry Smith *norm = 0.0; 473289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 474289bc588SBarry Smith sum = 0.0; 475289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 476289bc588SBarry Smith #if defined(PETSC_COMPLEX) 477289bc588SBarry Smith sum += abs(*v++); 478289bc588SBarry Smith #else 479289bc588SBarry Smith sum += fabs(*v++); 480289bc588SBarry Smith #endif 481289bc588SBarry Smith } 482289bc588SBarry Smith if (sum > *norm) *norm = sum; 483289bc588SBarry Smith } 484289bc588SBarry Smith } 485289bc588SBarry Smith else if (type == NORM_INFINITY) { 486289bc588SBarry Smith *norm = 0.0; 487289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 488289bc588SBarry Smith v = mat->v + j; 489289bc588SBarry Smith sum = 0.0; 490289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 491289bc588SBarry Smith #if defined(PETSC_COMPLEX) 492289bc588SBarry Smith sum += abs(*v); v += mat->m; 493289bc588SBarry Smith #else 494289bc588SBarry Smith sum += fabs(*v); v += mat->m; 495289bc588SBarry Smith #endif 496289bc588SBarry Smith } 497289bc588SBarry Smith if (sum > *norm) *norm = sum; 498289bc588SBarry Smith } 499289bc588SBarry Smith } 500289bc588SBarry Smith else { 501289bc588SBarry Smith SETERR(1,"No support for the two norm yet"); 502289bc588SBarry Smith } 503289bc588SBarry Smith return 0; 504289bc588SBarry Smith } 505289bc588SBarry Smith 50644a69424SLois Curfman McInnes static int MatInsOpt_Dense(Mat aijin,int op) 507289bc588SBarry Smith { 50844a69424SLois Curfman McInnes Mat_Dense *aij = (Mat_Dense *) aijin->data; 509289bc588SBarry Smith if (op == ROW_ORIENTED) aij->roworiented = 1; 510289bc588SBarry Smith else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 511289bc588SBarry Smith /* doesn't care about sorted rows or columns */ 512289bc588SBarry Smith return 0; 513289bc588SBarry Smith } 514289bc588SBarry Smith 51544a69424SLois Curfman McInnes static int MatZero_Dense(Mat A) 5166f0a148fSBarry Smith { 51744a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) A->data; 5186f0a148fSBarry Smith MEMSET(l->v,0,l->m*l->n*sizeof(Scalar)); 5196f0a148fSBarry Smith return 0; 5206f0a148fSBarry Smith } 5216f0a148fSBarry Smith 52244a69424SLois Curfman McInnes static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag) 5236f0a148fSBarry Smith { 52444a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) A->data; 5256abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 5266f0a148fSBarry Smith Scalar *slot; 527260925b8SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERR(ierr); 528260925b8SBarry Smith ierr = ISGetIndices(is,&rows); CHKERR(ierr); 5296f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5306f0a148fSBarry Smith slot = l->v + rows[i]; 5316f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 5326f0a148fSBarry Smith } 5336f0a148fSBarry Smith if (diag) { 5346f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5356f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 536260925b8SBarry Smith *slot = *diag; 5376f0a148fSBarry Smith } 5386f0a148fSBarry Smith } 539260925b8SBarry Smith ISRestoreIndices(is,&rows); 5406f0a148fSBarry Smith return 0; 5416f0a148fSBarry Smith } 542557bce09SLois Curfman McInnes 543fa9ec3f1SLois Curfman McInnes static int MatGetSize_Dense(Mat matin,int *m,int *n) 544557bce09SLois Curfman McInnes { 54544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 546557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 547557bce09SLois Curfman McInnes return 0; 548557bce09SLois Curfman McInnes } 549557bce09SLois Curfman McInnes 55044a69424SLois Curfman McInnes static int MatGetArray_Dense(Mat matin,Scalar **array) 55164e87e97SBarry Smith { 55244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 55364e87e97SBarry Smith *array = mat->v; 55464e87e97SBarry Smith return 0; 55564e87e97SBarry Smith } 556289bc588SBarry Smith /* -------------------------------------------------------------------*/ 55744a69424SLois Curfman McInnes static struct _MatOps MatOps = {MatInsert_Dense, 55844a69424SLois Curfman McInnes MatGetRow_Dense, MatRestoreRow_Dense, 55944a69424SLois Curfman McInnes MatMult_Dense, MatMultAdd_Dense, 56044a69424SLois Curfman McInnes MatMultTrans_Dense, MatMultTransAdd_Dense, 56144a69424SLois Curfman McInnes MatSolve_Dense,MatSolveAdd_Dense, 56244a69424SLois Curfman McInnes MatSolveTrans_Dense,MatSolveTransAdd_Dense, 56344a69424SLois Curfman McInnes MatLUFactor_Dense,MatChFactor_Dense, 56444a69424SLois Curfman McInnes MatRelax_Dense, 56544a69424SLois Curfman McInnes MatTrans_Dense, 566fa9ec3f1SLois Curfman McInnes MatGetInfo_Dense,MatEqual_Dense, 56744a69424SLois Curfman McInnes MatCopy_Dense, 56844a69424SLois Curfman McInnes MatGetDiag_Dense,MatScale_Dense,MatNorm_Dense, 569289bc588SBarry Smith 0,0, 57044a69424SLois Curfman McInnes 0, MatInsOpt_Dense,MatZero_Dense,MatZeroRows_Dense,0, 57144a69424SLois Curfman McInnes MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense, 57244a69424SLois Curfman McInnes MatChFactorSymbolic_Dense,MatChFactorNumeric_Dense, 573fa9ec3f1SLois Curfman McInnes MatGetSize_Dense,MatGetSize_Dense,0, 57444a69424SLois Curfman McInnes 0,0,MatGetArray_Dense 575289bc588SBarry Smith }; 57620563c6bSBarry Smith /*@ 57720563c6bSBarry Smith MatCreateSequentialDense - Creates a sequential dense matrix that 57820563c6bSBarry Smith is stored in the usual Fortran 77 manner. Many of the matrix 57920563c6bSBarry Smith operations use the BLAS and LAPACK routines. 580289bc588SBarry Smith 58120563c6bSBarry Smith Input Parameters: 58220563c6bSBarry Smith . m, n - the number of rows and columns in the matrix. 583*6b5873e3SBarry Smith . comm - should be MPI_COMM_SELF 58420563c6bSBarry Smith 58520563c6bSBarry Smith Output Parameter: 58620563c6bSBarry Smith . newmat - the matrix created. 58720563c6bSBarry Smith 58820563c6bSBarry Smith Keywords: dense matrix, lapack, blas 58920563c6bSBarry Smith @*/ 590*6b5873e3SBarry Smith int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat) 591289bc588SBarry Smith { 59244a69424SLois Curfman McInnes int size = sizeof(Mat_Dense) + m*n*sizeof(Scalar); 593289bc588SBarry Smith Mat mat; 59444a69424SLois Curfman McInnes Mat_Dense *l; 595289bc588SBarry Smith *newmat = 0; 596*6b5873e3SBarry Smith PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm); 597a5a9c739SBarry Smith PLogObjectCreate(mat); 59844a69424SLois Curfman McInnes l = (Mat_Dense *) MALLOC(size); CHKPTR(l); 599289bc588SBarry Smith mat->ops = &MatOps; 60044a69424SLois Curfman McInnes mat->destroy = MatDestroy_Dense; 60144a69424SLois Curfman McInnes mat->view = MatView_Dense; 602289bc588SBarry Smith mat->data = (void *) l; 603289bc588SBarry Smith mat->factor = 0; 604289bc588SBarry Smith mat->col = 0; 605289bc588SBarry Smith mat->row = 0; 606d6dfbf8fSBarry Smith 607289bc588SBarry Smith l->m = m; 608289bc588SBarry Smith l->n = n; 609289bc588SBarry Smith l->v = (Scalar *) (l + 1); 610289bc588SBarry Smith l->pivots = 0; 611289bc588SBarry Smith l->roworiented = 1; 612d6dfbf8fSBarry Smith 613289bc588SBarry Smith MEMSET(l->v,0,m*n*sizeof(Scalar)); 614289bc588SBarry Smith *newmat = mat; 615289bc588SBarry Smith return 0; 616289bc588SBarry Smith } 617289bc588SBarry Smith 61844a69424SLois Curfman McInnes int MatCreate_Dense(Mat matin,Mat *newmat) 619289bc588SBarry Smith { 62044a69424SLois Curfman McInnes Mat_Dense *m = (Mat_Dense *) matin->data; 621*6b5873e3SBarry Smith return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat); 622289bc588SBarry Smith } 623