15392566eSBarry Smith 25392566eSBarry Smith 3cb512458SBarry Smith #ifndef lint 4*49d8b64dSBarry Smith static char vcid[] = "$Id: dense.c,v 1.42 1995/07/06 17:19:31 bsmith Exp bsmith $"; 5cb512458SBarry Smith #endif 6289bc588SBarry Smith 7289bc588SBarry Smith /* 8289bc588SBarry Smith Standard Fortran style matrices 9289bc588SBarry Smith */ 10289bc588SBarry Smith #include "ptscimpl.h" 11289bc588SBarry Smith #include "plapack.h" 12289bc588SBarry Smith #include "matimpl.h" 13289bc588SBarry Smith #include "math.h" 14289bc588SBarry Smith #include "vec/vecimpl.h" 1520e1a0d4SLois Curfman McInnes #include "pviewer.h" 16289bc588SBarry Smith 17289bc588SBarry Smith typedef struct { 18289bc588SBarry Smith Scalar *v; 19289bc588SBarry Smith int roworiented; 20289bc588SBarry Smith int m,n,pad; 21289bc588SBarry Smith int *pivots; /* pivots in LU factorization */ 2244a69424SLois Curfman McInnes } Mat_Dense; 23289bc588SBarry Smith 2420e1a0d4SLois Curfman McInnes static int MatGetInfo_Dense(Mat matin,MatInfoType flag,int *nz, 2520e1a0d4SLois Curfman McInnes int *nzalloc,int *mem) 26289bc588SBarry Smith { 2744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 28289bc588SBarry Smith int i,N = mat->m*mat->n,count = 0; 29289bc588SBarry Smith Scalar *v = mat->v; 30289bc588SBarry Smith for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 31fa9ec3f1SLois Curfman McInnes *nz = count; *nzalloc = N; *mem = N*sizeof(Scalar); 32fa9ec3f1SLois Curfman McInnes return 0; 33289bc588SBarry Smith } 34289bc588SBarry Smith 35289bc588SBarry Smith /* ---------------------------------------------------------------*/ 36289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots, 37289bc588SBarry Smith rather than put it in the Mat->col slot.*/ 38*49d8b64dSBarry Smith static int MatLUFactor_Dense(Mat matin,IS row,IS col,double f) 39289bc588SBarry Smith { 4044a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 416abc6512SBarry Smith int info; 42289bc588SBarry Smith if (!mat->pivots) { 4378b31e54SBarry Smith mat->pivots = (int *) PETSCMALLOC( mat->m*sizeof(int) ); 4478b31e54SBarry Smith CHKPTRQ(mat->pivots); 45289bc588SBarry Smith } 46289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 4778b31e54SBarry Smith if (info) SETERRQ(1,"Bad LU factorization"); 48289bc588SBarry Smith matin->factor = FACTOR_LU; 49289bc588SBarry Smith return 0; 50289bc588SBarry Smith } 51*49d8b64dSBarry Smith static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,double f, 52*49d8b64dSBarry Smith Mat *fact) 53289bc588SBarry Smith { 54289bc588SBarry Smith int ierr; 5578b31e54SBarry Smith if ((ierr = MatConvert(matin,MATSAME,fact))) SETERRQ(ierr,0); 56289bc588SBarry Smith return 0; 57289bc588SBarry Smith } 5844a69424SLois Curfman McInnes static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact) 59289bc588SBarry Smith { 60*49d8b64dSBarry Smith return MatLUFactor(*fact,0,0,1.0); 61289bc588SBarry Smith } 62*49d8b64dSBarry Smith static int MatCholeskyFactorSymbolic_Dense(Mat matin,IS row,double f,Mat *fact) 63289bc588SBarry Smith { 64289bc588SBarry Smith int ierr; 6578b31e54SBarry Smith if ((ierr = MatConvert(matin,MATSAME,fact))) SETERRQ(ierr,0); 66289bc588SBarry Smith return 0; 67289bc588SBarry Smith } 6846fc4a8fSLois Curfman McInnes static int MatCholeskyFactorNumeric_Dense(Mat matin,Mat *fact) 69289bc588SBarry Smith { 70*49d8b64dSBarry Smith return MatCholeskyFactor(*fact,0,1.0); 71289bc588SBarry Smith } 72*49d8b64dSBarry Smith static int MatCholeskyFactor_Dense(Mat matin,IS perm,double f) 73289bc588SBarry Smith { 7444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 756abc6512SBarry Smith int info; 7678b31e54SBarry Smith if (mat->pivots) {PETSCFREE(mat->pivots); mat->pivots = 0;} 77289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 7878b31e54SBarry Smith if (info) SETERRQ(1,"Bad Cholesky factorization"); 79289bc588SBarry Smith matin->factor = FACTOR_CHOLESKY; 80289bc588SBarry Smith return 0; 81289bc588SBarry Smith } 82289bc588SBarry Smith 8344a69424SLois Curfman McInnes static int MatSolve_Dense(Mat matin,Vec xx,Vec yy) 84289bc588SBarry Smith { 8544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 866abc6512SBarry Smith int one = 1, info; 876abc6512SBarry Smith Scalar *x, *y; 88289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 8978b31e54SBarry Smith PETSCMEMCPY(y,x,mat->m*sizeof(Scalar)); 909e25ed09SBarry Smith if (matin->factor == FACTOR_LU) { 91289bc588SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 92289bc588SBarry Smith y, &mat->m, &info ); 93289bc588SBarry Smith } 949e25ed09SBarry Smith else if (matin->factor == FACTOR_CHOLESKY){ 95289bc588SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 96289bc588SBarry Smith y, &mat->m, &info ); 97289bc588SBarry Smith } 9878b31e54SBarry Smith else SETERRQ(1,"Matrix must be factored to solve"); 9978b31e54SBarry Smith if (info) SETERRQ(1,"Bad solve"); 100289bc588SBarry Smith return 0; 101289bc588SBarry Smith } 10244a69424SLois Curfman McInnes static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy) 103da3a660dSBarry Smith { 10444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1056abc6512SBarry Smith int one = 1, info; 1066abc6512SBarry Smith Scalar *x, *y; 107da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 10878b31e54SBarry Smith PETSCMEMCPY(y,x,mat->m*sizeof(Scalar)); 109da3a660dSBarry Smith /* assume if pivots exist then LU else Cholesky */ 110da3a660dSBarry Smith if (mat->pivots) { 111da3a660dSBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 112da3a660dSBarry Smith y, &mat->m, &info ); 113da3a660dSBarry Smith } 114da3a660dSBarry Smith else { 115da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 116da3a660dSBarry Smith y, &mat->m, &info ); 117da3a660dSBarry Smith } 11878b31e54SBarry Smith if (info) SETERRQ(1,"Bad solve"); 119da3a660dSBarry Smith return 0; 120da3a660dSBarry Smith } 12144a69424SLois Curfman McInnes static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 122da3a660dSBarry Smith { 12344a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1246abc6512SBarry Smith int one = 1, info,ierr; 1256abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 126da3a660dSBarry Smith Vec tmp = 0; 127da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 128da3a660dSBarry Smith if (yy == zz) { 12978b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 13078b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 131da3a660dSBarry Smith } 13278b31e54SBarry Smith PETSCMEMCPY(y,x,mat->m*sizeof(Scalar)); 133da3a660dSBarry Smith /* assume if pivots exist then LU else Cholesky */ 134da3a660dSBarry Smith if (mat->pivots) { 135da3a660dSBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 136da3a660dSBarry Smith y, &mat->m, &info ); 137da3a660dSBarry Smith } 138da3a660dSBarry Smith else { 139da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 140da3a660dSBarry Smith y, &mat->m, &info ); 141da3a660dSBarry Smith } 14278b31e54SBarry Smith if (info) SETERRQ(1,"Bad solve"); 143da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 144da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 145da3a660dSBarry Smith return 0; 146da3a660dSBarry Smith } 14744a69424SLois Curfman McInnes static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy) 148da3a660dSBarry Smith { 14944a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1506abc6512SBarry Smith int one = 1, info,ierr; 1516abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 152da3a660dSBarry Smith Vec tmp; 153da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 154da3a660dSBarry Smith if (yy == zz) { 15578b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 15678b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 157da3a660dSBarry Smith } 15878b31e54SBarry Smith PETSCMEMCPY(y,x,mat->m*sizeof(Scalar)); 159da3a660dSBarry Smith /* assume if pivots exist then LU else Cholesky */ 160da3a660dSBarry Smith if (mat->pivots) { 161da3a660dSBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 162da3a660dSBarry Smith y, &mat->m, &info ); 163da3a660dSBarry Smith } 164da3a660dSBarry Smith else { 165da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 166da3a660dSBarry Smith y, &mat->m, &info ); 167da3a660dSBarry Smith } 16878b31e54SBarry Smith if (info) SETERRQ(1,"Bad solve"); 169da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 170da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 171da3a660dSBarry Smith return 0; 172da3a660dSBarry Smith } 173289bc588SBarry Smith /* ------------------------------------------------------------------*/ 17420e1a0d4SLois Curfman McInnes static int MatRelax_Dense(Mat matin,Vec bb,double omega,MatSORType flag, 17520e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 176289bc588SBarry Smith { 17744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 178289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 1796abc6512SBarry Smith int o = 1,ierr, m = mat->m, i; 180289bc588SBarry Smith 181289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 182289bc588SBarry Smith /* this is a hack fix, should have another version without 183289bc588SBarry Smith the second BLdot */ 18478b31e54SBarry Smith if ((ierr = VecSet(&zero,xx))) SETERRQ(ierr,0); 185289bc588SBarry Smith } 186289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 187289bc588SBarry Smith while (its--) { 188289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 189289bc588SBarry Smith for ( i=0; i<m; i++ ) { 190289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 191d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 192289bc588SBarry Smith } 193289bc588SBarry Smith } 194289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 195289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 196289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 197d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 198289bc588SBarry Smith } 199289bc588SBarry Smith } 200289bc588SBarry Smith } 201289bc588SBarry Smith return 0; 202289bc588SBarry Smith } 203289bc588SBarry Smith 204289bc588SBarry Smith /* -----------------------------------------------------------------*/ 20544a69424SLois Curfman McInnes static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy) 206289bc588SBarry Smith { 20744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 208289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 209289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 210289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 211289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 212289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 213289bc588SBarry Smith return 0; 214289bc588SBarry Smith } 21544a69424SLois Curfman McInnes static int MatMult_Dense(Mat matin,Vec xx,Vec yy) 216289bc588SBarry Smith { 21744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 218289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 219289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 220289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 221289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 222289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 223289bc588SBarry Smith return 0; 224289bc588SBarry Smith } 22544a69424SLois Curfman McInnes static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 226289bc588SBarry Smith { 22744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 228289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 2296abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 230289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 23178b31e54SBarry Smith if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar)); 232289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 233289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 234289bc588SBarry Smith return 0; 235289bc588SBarry Smith } 23644a69424SLois Curfman McInnes static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 237289bc588SBarry Smith { 23844a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 239289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 240289bc588SBarry Smith int _One=1; 2416abc6512SBarry Smith Scalar _DOne=1.0; 242289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 243289bc588SBarry Smith VecGetArray(zz,&z); 24478b31e54SBarry Smith if (zz != yy) PETSCMEMCPY(y,z,mat->m*sizeof(Scalar)); 245289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 246289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 247289bc588SBarry Smith return 0; 248289bc588SBarry Smith } 249289bc588SBarry Smith 250289bc588SBarry Smith /* -----------------------------------------------------------------*/ 25144a69424SLois Curfman McInnes static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols, 252289bc588SBarry Smith Scalar **vals) 253289bc588SBarry Smith { 25444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 255289bc588SBarry Smith Scalar *v; 256289bc588SBarry Smith int i; 257289bc588SBarry Smith *ncols = mat->n; 258289bc588SBarry Smith if (cols) { 25978b31e54SBarry Smith *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols); 260289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) *cols[i] = i; 261289bc588SBarry Smith } 262289bc588SBarry Smith if (vals) { 26378b31e54SBarry Smith *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 264289bc588SBarry Smith v = mat->v + row; 265289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 266289bc588SBarry Smith } 267289bc588SBarry Smith return 0; 268289bc588SBarry Smith } 26944a69424SLois Curfman McInnes static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols, 270289bc588SBarry Smith Scalar **vals) 271289bc588SBarry Smith { 27278b31e54SBarry Smith if (cols) { PETSCFREE(*cols); } 27378b31e54SBarry Smith if (vals) { PETSCFREE(*vals); } 274289bc588SBarry Smith return 0; 275289bc588SBarry Smith } 276289bc588SBarry Smith /* ----------------------------------------------------------------*/ 27744a69424SLois Curfman McInnes static int MatInsert_Dense(Mat matin,int m,int *indexm,int n, 278e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 279289bc588SBarry Smith { 28044a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 281289bc588SBarry Smith int i,j; 282d6dfbf8fSBarry Smith 283289bc588SBarry Smith if (!mat->roworiented) { 284edae2e7dSBarry Smith if (addv == INSERTVALUES) { 285289bc588SBarry Smith for ( j=0; j<n; j++ ) { 286289bc588SBarry Smith for ( i=0; i<m; i++ ) { 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++ ) { 293289bc588SBarry Smith for ( i=0; i<m; i++ ) { 294289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 295289bc588SBarry Smith } 296289bc588SBarry Smith } 297289bc588SBarry Smith } 298e8d4e0b9SBarry Smith } 299e8d4e0b9SBarry Smith else { 300edae2e7dSBarry Smith if (addv == INSERTVALUES) { 301e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 302e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 303e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 304e8d4e0b9SBarry Smith } 305e8d4e0b9SBarry Smith } 306e8d4e0b9SBarry Smith } 307289bc588SBarry Smith else { 308289bc588SBarry Smith for ( i=0; i<m; i++ ) { 309289bc588SBarry Smith for ( j=0; j<n; j++ ) { 310289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 311289bc588SBarry Smith } 312289bc588SBarry Smith } 313289bc588SBarry Smith } 314e8d4e0b9SBarry Smith } 315289bc588SBarry Smith return 0; 316289bc588SBarry Smith } 317e8d4e0b9SBarry Smith 318289bc588SBarry Smith /* -----------------------------------------------------------------*/ 319ff7509f2SBarry Smith static int MatCopyPrivate_Dense(Mat matin,Mat *newmat) 320289bc588SBarry Smith { 32144a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 322289bc588SBarry Smith int ierr; 323289bc588SBarry Smith Mat newi; 32444a69424SLois Curfman McInnes Mat_Dense *l; 3256b5873e3SBarry Smith if ((ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi))) 32678b31e54SBarry Smith SETERRQ(ierr,0); 32744a69424SLois Curfman McInnes l = (Mat_Dense *) newi->data; 32878b31e54SBarry Smith PETSCMEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 329289bc588SBarry Smith *newmat = newi; 330289bc588SBarry Smith return 0; 331289bc588SBarry Smith } 332da3a660dSBarry Smith #include "viewer.h" 333289bc588SBarry Smith 33444a69424SLois Curfman McInnes int MatView_Dense(PetscObject obj,Viewer ptr) 335289bc588SBarry Smith { 336289bc588SBarry Smith Mat matin = (Mat) obj; 33744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 338289bc588SBarry Smith Scalar *v; 339289bc588SBarry Smith int i,j; 34023242f5aSBarry Smith PetscObject vobj = (PetscObject) ptr; 341da3a660dSBarry Smith 34223242f5aSBarry Smith if (ptr == 0) { 34323242f5aSBarry Smith ptr = STDOUT_VIEWER; vobj = (PetscObject) ptr; 34423242f5aSBarry Smith } 34523242f5aSBarry Smith if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) { 3466f75cfa4SBarry Smith return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v); 347da3a660dSBarry Smith } 348da3a660dSBarry Smith else { 3494a0ce102SLois Curfman McInnes FILE *fd = ViewerFileGetPointer_Private(ptr); 350289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 351289bc588SBarry Smith v = mat->v + i; 352289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 353289bc588SBarry Smith #if defined(PETSC_COMPLEX) 3548e37dc09SBarry Smith fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 355289bc588SBarry Smith #else 3568e37dc09SBarry Smith fprintf(fd,"%6.4e ",*v); v += mat->m; 357289bc588SBarry Smith #endif 358289bc588SBarry Smith } 3598e37dc09SBarry Smith fprintf(fd,"\n"); 360289bc588SBarry Smith } 361da3a660dSBarry Smith } 362289bc588SBarry Smith return 0; 363289bc588SBarry Smith } 364289bc588SBarry Smith 365289bc588SBarry Smith 36644a69424SLois Curfman McInnes static int MatDestroy_Dense(PetscObject obj) 367289bc588SBarry Smith { 368289bc588SBarry Smith Mat mat = (Mat) obj; 36944a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) mat->data; 370a5a9c739SBarry Smith #if defined(PETSC_LOG) 371a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 372a5a9c739SBarry Smith #endif 37378b31e54SBarry Smith if (l->pivots) PETSCFREE(l->pivots); 37478b31e54SBarry Smith PETSCFREE(l); 375a5a9c739SBarry Smith PLogObjectDestroy(mat); 3769e25ed09SBarry Smith PETSCHEADERDESTROY(mat); 377289bc588SBarry Smith return 0; 378289bc588SBarry Smith } 379289bc588SBarry Smith 38048b35521SBarry Smith static int MatTranspose_Dense(Mat matin,Mat *matout) 381289bc588SBarry Smith { 38244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 383289bc588SBarry Smith int k,j; 384289bc588SBarry Smith Scalar *v = mat->v, tmp; 38548b35521SBarry Smith 38648b35521SBarry Smith if (!matout) { /* in place transpose */ 387289bc588SBarry Smith if (mat->m != mat->n) { 38848b35521SBarry Smith SETERRQ(1,"MatTranspose_Dense:Cannot transpose rectangular matrix"); 389289bc588SBarry Smith } 390289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 391289bc588SBarry Smith for ( k=0; k<j; k++ ) { 392289bc588SBarry Smith tmp = v[j + k*mat->n]; 393289bc588SBarry Smith v[j + k*mat->n] = v[k + j*mat->n]; 394289bc588SBarry Smith v[k + j*mat->n] = tmp; 395289bc588SBarry Smith } 396289bc588SBarry Smith } 39748b35521SBarry Smith } 39848b35521SBarry Smith else { 39948b35521SBarry Smith SETERRQ(1,"MatTranspose_Dense:not out of place transpose yet"); 40048b35521SBarry 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 VecGetArray(v,&x); VecGetSize(v,&n); 42578b31e54SBarry Smith if (n != mat->m) SETERRQ(1,"Nonconforming matrix and vector"); 426289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 427289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 428289bc588SBarry Smith } 429289bc588SBarry Smith return 0; 430289bc588SBarry Smith } 431289bc588SBarry Smith 43244a69424SLois Curfman McInnes static int MatScale_Dense(Mat matin,Vec ll,Vec rr) 433289bc588SBarry Smith { 43444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 435da3a660dSBarry Smith Scalar *l,*r,x,*v; 436da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 43728988994SBarry Smith if (ll) { 438da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 43978b31e54SBarry Smith if (m != mat->m) SETERRQ(1,"Left scaling vector wrong length"); 440da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 441da3a660dSBarry Smith x = l[i]; 442da3a660dSBarry Smith v = mat->v + i; 443da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 444da3a660dSBarry Smith } 445da3a660dSBarry Smith } 44628988994SBarry Smith if (rr) { 447da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 44878b31e54SBarry Smith if (n != mat->n) SETERRQ(1,"Right scaling vector wrong length"); 449da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 450da3a660dSBarry Smith x = r[i]; 451da3a660dSBarry Smith v = mat->v + i*m; 452da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 453da3a660dSBarry Smith } 454da3a660dSBarry Smith } 455289bc588SBarry Smith return 0; 456289bc588SBarry Smith } 457289bc588SBarry Smith 458da3a660dSBarry Smith 45944a69424SLois Curfman McInnes static int MatNorm_Dense(Mat matin,int type,double *norm) 460289bc588SBarry Smith { 46144a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 462289bc588SBarry Smith Scalar *v = mat->v; 463289bc588SBarry Smith double sum = 0.0; 464289bc588SBarry Smith int i, j; 465289bc588SBarry Smith if (type == NORM_FROBENIUS) { 466289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 467289bc588SBarry Smith #if defined(PETSC_COMPLEX) 468289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 469289bc588SBarry Smith #else 470289bc588SBarry Smith sum += (*v)*(*v); v++; 471289bc588SBarry Smith #endif 472289bc588SBarry Smith } 473289bc588SBarry Smith *norm = sqrt(sum); 474289bc588SBarry Smith } 475289bc588SBarry Smith else if (type == NORM_1) { 476289bc588SBarry Smith *norm = 0.0; 477289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 478289bc588SBarry Smith sum = 0.0; 479289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 480289bc588SBarry Smith #if defined(PETSC_COMPLEX) 481289bc588SBarry Smith sum += abs(*v++); 482289bc588SBarry Smith #else 483289bc588SBarry Smith sum += fabs(*v++); 484289bc588SBarry Smith #endif 485289bc588SBarry Smith } 486289bc588SBarry Smith if (sum > *norm) *norm = sum; 487289bc588SBarry Smith } 488289bc588SBarry Smith } 489289bc588SBarry Smith else if (type == NORM_INFINITY) { 490289bc588SBarry Smith *norm = 0.0; 491289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 492289bc588SBarry Smith v = mat->v + j; 493289bc588SBarry Smith sum = 0.0; 494289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 495289bc588SBarry Smith #if defined(PETSC_COMPLEX) 496289bc588SBarry Smith sum += abs(*v); v += mat->m; 497289bc588SBarry Smith #else 498289bc588SBarry Smith sum += fabs(*v); v += mat->m; 499289bc588SBarry Smith #endif 500289bc588SBarry Smith } 501289bc588SBarry Smith if (sum > *norm) *norm = sum; 502289bc588SBarry Smith } 503289bc588SBarry Smith } 504289bc588SBarry Smith else { 50578b31e54SBarry Smith SETERRQ(1,"No support for the two norm yet"); 506289bc588SBarry Smith } 507289bc588SBarry Smith return 0; 508289bc588SBarry Smith } 509289bc588SBarry Smith 51020e1a0d4SLois Curfman McInnes static int MatSetOption_Dense(Mat aijin,MatOption op) 511289bc588SBarry Smith { 51244a69424SLois Curfman McInnes Mat_Dense *aij = (Mat_Dense *) aijin->data; 513289bc588SBarry Smith if (op == ROW_ORIENTED) aij->roworiented = 1; 514289bc588SBarry Smith else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 515289bc588SBarry Smith /* doesn't care about sorted rows or columns */ 516289bc588SBarry Smith return 0; 517289bc588SBarry Smith } 518289bc588SBarry Smith 51944a69424SLois Curfman McInnes static int MatZero_Dense(Mat A) 5206f0a148fSBarry Smith { 52144a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) A->data; 52278b31e54SBarry Smith PETSCMEMSET(l->v,0,l->m*l->n*sizeof(Scalar)); 5236f0a148fSBarry Smith return 0; 5246f0a148fSBarry Smith } 5256f0a148fSBarry Smith 52644a69424SLois Curfman McInnes static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag) 5276f0a148fSBarry Smith { 52844a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) A->data; 5296abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 5306f0a148fSBarry Smith Scalar *slot; 53178b31e54SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 53278b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 5336f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5346f0a148fSBarry Smith slot = l->v + rows[i]; 5356f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 5366f0a148fSBarry Smith } 5376f0a148fSBarry Smith if (diag) { 5386f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5396f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 540260925b8SBarry Smith *slot = *diag; 5416f0a148fSBarry Smith } 5426f0a148fSBarry Smith } 543260925b8SBarry Smith ISRestoreIndices(is,&rows); 5446f0a148fSBarry Smith return 0; 5456f0a148fSBarry Smith } 546557bce09SLois Curfman McInnes 547fa9ec3f1SLois Curfman McInnes static int MatGetSize_Dense(Mat matin,int *m,int *n) 548557bce09SLois Curfman McInnes { 54944a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 550557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 551557bce09SLois Curfman McInnes return 0; 552557bce09SLois Curfman McInnes } 553557bce09SLois Curfman McInnes 55444a69424SLois Curfman McInnes static int MatGetArray_Dense(Mat matin,Scalar **array) 55564e87e97SBarry Smith { 55644a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 55764e87e97SBarry Smith *array = mat->v; 55864e87e97SBarry Smith return 0; 55964e87e97SBarry Smith } 5600754003eSLois Curfman McInnes 5610754003eSLois Curfman McInnes 5620754003eSLois Curfman McInnes static int MatGetSubMatrixInPlace_Dense(Mat matin,IS isrow,IS iscol) 5630754003eSLois Curfman McInnes { 56478b31e54SBarry Smith SETERRQ(1,"MatGetSubMatrixInPlace_Dense not yet done."); 5650754003eSLois Curfman McInnes } 5660754003eSLois Curfman McInnes 5670754003eSLois Curfman McInnes static int MatGetSubMatrix_Dense(Mat matin,IS isrow,IS iscol,Mat *submat) 5680754003eSLois Curfman McInnes { 5690754003eSLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 5700754003eSLois Curfman McInnes int nznew, *smap, i, j, ierr, oldcols = mat->n; 571160018dcSBarry Smith int *irow, *icol, nrows, ncols, *cwork; 5720754003eSLois Curfman McInnes Scalar *vwork, *val; 5730754003eSLois Curfman McInnes Mat newmat; 5740754003eSLois Curfman McInnes 57578b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 57678b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 57778b31e54SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 57878b31e54SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 5790754003eSLois Curfman McInnes 58078b31e54SBarry Smith smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap); 58178b31e54SBarry Smith cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork); 58278b31e54SBarry Smith vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 58378b31e54SBarry Smith PETSCMEMSET((char*)smap,0,oldcols*sizeof(int)); 5840754003eSLois Curfman McInnes for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 5850754003eSLois Curfman McInnes 5860754003eSLois Curfman McInnes /* Create and fill new matrix */ 5870754003eSLois Curfman McInnes ierr = MatCreateSequentialDense(matin->comm,nrows,ncols,&newmat); 58878b31e54SBarry Smith CHKERRQ(ierr); 5890754003eSLois Curfman McInnes for (i=0; i<nrows; i++) { 5900754003eSLois Curfman McInnes nznew = 0; 5910754003eSLois Curfman McInnes val = mat->v + irow[i]; 5920754003eSLois Curfman McInnes for (j=0; j<oldcols; j++) { 5930754003eSLois Curfman McInnes if (smap[j]) { 5940754003eSLois Curfman McInnes cwork[nznew] = smap[j] - 1; 5950754003eSLois Curfman McInnes vwork[nznew++] = val[j * mat->m]; 5960754003eSLois Curfman McInnes } 5970754003eSLois Curfman McInnes } 598edae2e7dSBarry Smith ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERTVALUES); 59978b31e54SBarry Smith CHKERRQ(ierr); 6000754003eSLois Curfman McInnes } 60178b31e54SBarry Smith ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 60278b31e54SBarry Smith ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 6030754003eSLois Curfman McInnes 6040754003eSLois Curfman McInnes /* Free work space */ 60578b31e54SBarry Smith PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork); 60678b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 60778b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 6080754003eSLois Curfman McInnes *submat = newmat; 6090754003eSLois Curfman McInnes return 0; 6100754003eSLois Curfman McInnes } 6110754003eSLois Curfman McInnes 612289bc588SBarry Smith /* -------------------------------------------------------------------*/ 61344a69424SLois Curfman McInnes static struct _MatOps MatOps = {MatInsert_Dense, 61444a69424SLois Curfman McInnes MatGetRow_Dense, MatRestoreRow_Dense, 61544a69424SLois Curfman McInnes MatMult_Dense, MatMultAdd_Dense, 61644a69424SLois Curfman McInnes MatMultTrans_Dense, MatMultTransAdd_Dense, 61744a69424SLois Curfman McInnes MatSolve_Dense,MatSolveAdd_Dense, 61844a69424SLois Curfman McInnes MatSolveTrans_Dense,MatSolveTransAdd_Dense, 61946fc4a8fSLois Curfman McInnes MatLUFactor_Dense,MatCholeskyFactor_Dense, 62044a69424SLois Curfman McInnes MatRelax_Dense, 62148b35521SBarry Smith MatTranspose_Dense, 622fa9ec3f1SLois Curfman McInnes MatGetInfo_Dense,MatEqual_Dense, 62346fc4a8fSLois Curfman McInnes MatGetDiagonal_Dense,MatScale_Dense,MatNorm_Dense, 624289bc588SBarry Smith 0,0, 625681d1451SLois Curfman McInnes 0, MatSetOption_Dense,MatZero_Dense,MatZeroRows_Dense,0, 62644a69424SLois Curfman McInnes MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense, 62746fc4a8fSLois Curfman McInnes MatCholeskyFactorSymbolic_Dense,MatCholeskyFactorNumeric_Dense, 628fa9ec3f1SLois Curfman McInnes MatGetSize_Dense,MatGetSize_Dense,0, 62907a0e7adSLois Curfman McInnes 0,0,MatGetArray_Dense,0,0, 63007a0e7adSLois Curfman McInnes MatGetSubMatrix_Dense,MatGetSubMatrixInPlace_Dense, 631ff7509f2SBarry Smith MatCopyPrivate_Dense}; 6320754003eSLois Curfman McInnes 63320563c6bSBarry Smith /*@ 63420563c6bSBarry Smith MatCreateSequentialDense - Creates a sequential dense matrix that 635d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 636d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 637289bc588SBarry Smith 63820563c6bSBarry Smith Input Parameters: 6390c775827SLois Curfman McInnes . comm - MPI communicator, set to MPI_COMM_SELF 6400c775827SLois Curfman McInnes . m - number of rows 6410c775827SLois Curfman McInnes . n - number of column 64220563c6bSBarry Smith 64320563c6bSBarry Smith Output Parameter: 6440c775827SLois Curfman McInnes . newmat - the matrix 64520563c6bSBarry Smith 646dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 647d65003e9SLois Curfman McInnes 648dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues() 64920563c6bSBarry Smith @*/ 6506b5873e3SBarry Smith int MatCreateSequentialDense(MPI_Comm comm,int m,int n,Mat *newmat) 651289bc588SBarry Smith { 65244a69424SLois Curfman McInnes int size = sizeof(Mat_Dense) + m*n*sizeof(Scalar); 653289bc588SBarry Smith Mat mat; 65444a69424SLois Curfman McInnes Mat_Dense *l; 655289bc588SBarry Smith *newmat = 0; 6566b5873e3SBarry Smith PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,comm); 657a5a9c739SBarry Smith PLogObjectCreate(mat); 65878b31e54SBarry Smith l = (Mat_Dense *) PETSCMALLOC(size); CHKPTRQ(l); 659289bc588SBarry Smith mat->ops = &MatOps; 66044a69424SLois Curfman McInnes mat->destroy = MatDestroy_Dense; 66144a69424SLois Curfman McInnes mat->view = MatView_Dense; 662289bc588SBarry Smith mat->data = (void *) l; 663289bc588SBarry Smith mat->factor = 0; 664d6dfbf8fSBarry Smith 665289bc588SBarry Smith l->m = m; 666289bc588SBarry Smith l->n = n; 667289bc588SBarry Smith l->v = (Scalar *) (l + 1); 668289bc588SBarry Smith l->pivots = 0; 669289bc588SBarry Smith l->roworiented = 1; 670d6dfbf8fSBarry Smith 67178b31e54SBarry Smith PETSCMEMSET(l->v,0,m*n*sizeof(Scalar)); 672289bc588SBarry Smith *newmat = mat; 673289bc588SBarry Smith return 0; 674289bc588SBarry Smith } 675289bc588SBarry Smith 67644a69424SLois Curfman McInnes int MatCreate_Dense(Mat matin,Mat *newmat) 677289bc588SBarry Smith { 67844a69424SLois Curfman McInnes Mat_Dense *m = (Mat_Dense *) matin->data; 6796b5873e3SBarry Smith return MatCreateSequentialDense(matin->comm,m->m,m->n,newmat); 680289bc588SBarry Smith } 681