1cb512458SBarry Smith #ifndef lint 2*dbd7a890SLois Curfman McInnes static char vcid[] = "$Id: dense.c,v 1.35 1995/05/03 13:18:03 bsmith Exp curfman $"; 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 2220e1a0d4SLois Curfman McInnes static int MatGetInfo_Dense(Mat matin,MatInfoType flag,int *nz, 2320e1a0d4SLois Curfman McInnes int *nzalloc,int *mem) 24289bc588SBarry Smith { 2544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 26289bc588SBarry Smith int i,N = mat->m*mat->n,count = 0; 27289bc588SBarry Smith Scalar *v = mat->v; 28289bc588SBarry Smith for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 29fa9ec3f1SLois Curfman McInnes *nz = count; *nzalloc = N; *mem = N*sizeof(Scalar); 30fa9ec3f1SLois Curfman McInnes return 0; 31289bc588SBarry Smith } 32289bc588SBarry Smith 33289bc588SBarry Smith /* ---------------------------------------------------------------*/ 34289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots, 35289bc588SBarry Smith rather than put it in the Mat->col slot.*/ 3644a69424SLois Curfman McInnes static int MatLUFactor_Dense(Mat matin,IS row,IS col) 37289bc588SBarry Smith { 3844a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 396abc6512SBarry Smith int info; 40289bc588SBarry Smith if (!mat->pivots) { 41289bc588SBarry Smith mat->pivots = (int *) MALLOC( mat->m*sizeof(int) ); 42289bc588SBarry Smith CHKPTR(mat->pivots); 43289bc588SBarry Smith } 44289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 45289bc588SBarry Smith if (info) SETERR(1,"Bad LU factorization"); 46289bc588SBarry Smith matin->factor = FACTOR_LU; 47289bc588SBarry Smith return 0; 48289bc588SBarry Smith } 4944a69424SLois Curfman McInnes static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,Mat *fact) 50289bc588SBarry Smith { 51289bc588SBarry Smith int ierr; 5207a0e7adSLois Curfman McInnes if ((ierr = MatConvert(matin,MATSAME,fact))) SETERR(ierr,0); 53289bc588SBarry Smith return 0; 54289bc588SBarry Smith } 5544a69424SLois Curfman McInnes static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact) 56289bc588SBarry Smith { 5720563c6bSBarry Smith return MatLUFactor(*fact,0,0); 58289bc588SBarry Smith } 5946fc4a8fSLois Curfman McInnes static int MatCholeskyFactorSymbolic_Dense(Mat matin,IS row,Mat *fact) 60289bc588SBarry Smith { 61289bc588SBarry Smith int ierr; 6207a0e7adSLois Curfman McInnes if ((ierr = MatConvert(matin,MATSAME,fact))) SETERR(ierr,0); 63289bc588SBarry Smith return 0; 64289bc588SBarry Smith } 6546fc4a8fSLois Curfman McInnes static int MatCholeskyFactorNumeric_Dense(Mat matin,Mat *fact) 66289bc588SBarry Smith { 6720563c6bSBarry Smith return MatCholeskyFactor(*fact,0); 68289bc588SBarry Smith } 6946fc4a8fSLois Curfman McInnes static int MatCholeskyFactor_Dense(Mat matin,IS perm) 70289bc588SBarry Smith { 7144a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 726abc6512SBarry Smith int info; 73289bc588SBarry Smith if (mat->pivots) {FREE(mat->pivots); mat->pivots = 0;} 74289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 75289bc588SBarry Smith if (info) SETERR(1,"Bad Cholesky factorization"); 76289bc588SBarry Smith matin->factor = FACTOR_CHOLESKY; 77289bc588SBarry Smith return 0; 78289bc588SBarry Smith } 79289bc588SBarry Smith 8044a69424SLois Curfman McInnes static int MatSolve_Dense(Mat matin,Vec xx,Vec yy) 81289bc588SBarry Smith { 8244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 836abc6512SBarry Smith int one = 1, info; 846abc6512SBarry Smith Scalar *x, *y; 85289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 86289bc588SBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 879e25ed09SBarry Smith if (matin->factor == FACTOR_LU) { 88289bc588SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 89289bc588SBarry Smith y, &mat->m, &info ); 90289bc588SBarry Smith } 919e25ed09SBarry Smith else if (matin->factor == FACTOR_CHOLESKY){ 92289bc588SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 93289bc588SBarry Smith y, &mat->m, &info ); 94289bc588SBarry Smith } 959e25ed09SBarry Smith else SETERR(1,"Matrix must be factored to solve"); 96289bc588SBarry Smith if (info) SETERR(1,"Bad solve"); 97289bc588SBarry Smith return 0; 98289bc588SBarry Smith } 9944a69424SLois Curfman McInnes static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy) 100da3a660dSBarry Smith { 10144a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1026abc6512SBarry Smith int one = 1, info; 1036abc6512SBarry Smith Scalar *x, *y; 104da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 105da3a660dSBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 106da3a660dSBarry Smith /* assume if pivots exist then LU else Cholesky */ 107da3a660dSBarry Smith if (mat->pivots) { 108da3a660dSBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 109da3a660dSBarry Smith y, &mat->m, &info ); 110da3a660dSBarry Smith } 111da3a660dSBarry Smith else { 112da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 113da3a660dSBarry Smith y, &mat->m, &info ); 114da3a660dSBarry Smith } 115da3a660dSBarry Smith if (info) SETERR(1,"Bad solve"); 116da3a660dSBarry Smith return 0; 117da3a660dSBarry Smith } 11844a69424SLois Curfman McInnes static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 119da3a660dSBarry Smith { 12044a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1216abc6512SBarry Smith int one = 1, info,ierr; 1226abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 123da3a660dSBarry Smith Vec tmp = 0; 124da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 125da3a660dSBarry Smith if (yy == zz) { 1266469c4f9SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERR(ierr); 127da3a660dSBarry Smith ierr = VecCopy(yy,tmp); CHKERR(ierr); 128da3a660dSBarry Smith } 129da3a660dSBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 130da3a660dSBarry Smith /* assume if pivots exist then LU else Cholesky */ 131da3a660dSBarry Smith if (mat->pivots) { 132da3a660dSBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 133da3a660dSBarry Smith y, &mat->m, &info ); 134da3a660dSBarry Smith } 135da3a660dSBarry Smith else { 136da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 137da3a660dSBarry Smith y, &mat->m, &info ); 138da3a660dSBarry Smith } 139da3a660dSBarry Smith if (info) SETERR(1,"Bad solve"); 140da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 141da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 142da3a660dSBarry Smith return 0; 143da3a660dSBarry Smith } 14444a69424SLois Curfman McInnes static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy) 145da3a660dSBarry Smith { 14644a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 1476abc6512SBarry Smith int one = 1, info,ierr; 1486abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 149da3a660dSBarry Smith Vec tmp; 150da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 151da3a660dSBarry Smith if (yy == zz) { 1526469c4f9SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERR(ierr); 153da3a660dSBarry Smith ierr = VecCopy(yy,tmp); CHKERR(ierr); 154da3a660dSBarry Smith } 155da3a660dSBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 156da3a660dSBarry Smith /* assume if pivots exist then LU else Cholesky */ 157da3a660dSBarry Smith if (mat->pivots) { 158da3a660dSBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 159da3a660dSBarry Smith y, &mat->m, &info ); 160da3a660dSBarry Smith } 161da3a660dSBarry Smith else { 162da3a660dSBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 163da3a660dSBarry Smith y, &mat->m, &info ); 164da3a660dSBarry Smith } 165da3a660dSBarry Smith if (info) SETERR(1,"Bad solve"); 166da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 167da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 168da3a660dSBarry Smith return 0; 169da3a660dSBarry Smith } 170289bc588SBarry Smith /* ------------------------------------------------------------------*/ 17120e1a0d4SLois Curfman McInnes static int MatRelax_Dense(Mat matin,Vec bb,double omega,MatSORType flag, 17220e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 173289bc588SBarry Smith { 17444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 175289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 1766abc6512SBarry Smith int o = 1,ierr, m = mat->m, i; 177289bc588SBarry Smith 178289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 179289bc588SBarry Smith /* this is a hack fix, should have another version without 180289bc588SBarry Smith the second BLdot */ 1816abc6512SBarry Smith if ((ierr = VecSet(&zero,xx))) SETERR(ierr,0); 182289bc588SBarry Smith } 183289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 184289bc588SBarry Smith while (its--) { 185289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 186289bc588SBarry Smith for ( i=0; i<m; i++ ) { 187289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 188d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 189289bc588SBarry Smith } 190289bc588SBarry Smith } 191289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 192289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 193289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 194d6dfbf8fSBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 195289bc588SBarry Smith } 196289bc588SBarry Smith } 197289bc588SBarry Smith } 198289bc588SBarry Smith return 0; 199289bc588SBarry Smith } 200289bc588SBarry Smith 201289bc588SBarry Smith /* -----------------------------------------------------------------*/ 20244a69424SLois Curfman McInnes static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy) 203289bc588SBarry Smith { 20444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 205289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 206289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 207289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 208289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 209289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 210289bc588SBarry Smith return 0; 211289bc588SBarry Smith } 21244a69424SLois Curfman McInnes static int MatMult_Dense(Mat matin,Vec xx,Vec yy) 213289bc588SBarry Smith { 21444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 215289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 216289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 217289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 218289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 219289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 220289bc588SBarry Smith return 0; 221289bc588SBarry Smith } 22244a69424SLois Curfman McInnes static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 223289bc588SBarry Smith { 22444a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 225289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 2266abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 227289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 228289bc588SBarry Smith if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 229289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 230289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 231289bc588SBarry Smith return 0; 232289bc588SBarry Smith } 23344a69424SLois Curfman McInnes static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 234289bc588SBarry Smith { 23544a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 236289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 237289bc588SBarry Smith int _One=1; 2386abc6512SBarry Smith Scalar _DOne=1.0; 239289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 240289bc588SBarry Smith VecGetArray(zz,&z); 241289bc588SBarry Smith if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 242289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 243289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 244289bc588SBarry Smith return 0; 245289bc588SBarry Smith } 246289bc588SBarry Smith 247289bc588SBarry Smith /* -----------------------------------------------------------------*/ 24844a69424SLois Curfman McInnes static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols, 249289bc588SBarry Smith Scalar **vals) 250289bc588SBarry Smith { 25144a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 252289bc588SBarry Smith Scalar *v; 253289bc588SBarry Smith int i; 254289bc588SBarry Smith *ncols = mat->n; 255289bc588SBarry Smith if (cols) { 256289bc588SBarry Smith *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols); 257289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) *cols[i] = i; 258289bc588SBarry Smith } 259289bc588SBarry Smith if (vals) { 260289bc588SBarry Smith *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals); 261289bc588SBarry Smith v = mat->v + row; 262289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 263289bc588SBarry Smith } 264289bc588SBarry Smith return 0; 265289bc588SBarry Smith } 26644a69424SLois Curfman McInnes static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols, 267289bc588SBarry Smith Scalar **vals) 268289bc588SBarry Smith { 269289bc588SBarry Smith if (cols) { FREE(*cols); } 270289bc588SBarry Smith if (vals) { FREE(*vals); } 271289bc588SBarry Smith return 0; 272289bc588SBarry Smith } 273289bc588SBarry Smith /* ----------------------------------------------------------------*/ 27444a69424SLois Curfman McInnes static int MatInsert_Dense(Mat matin,int m,int *indexm,int n, 275e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 276289bc588SBarry Smith { 27744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 278289bc588SBarry Smith int i,j; 279d6dfbf8fSBarry Smith 280289bc588SBarry Smith if (!mat->roworiented) { 281edae2e7dSBarry Smith if (addv == INSERTVALUES) { 282289bc588SBarry Smith for ( j=0; j<n; j++ ) { 283d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v += m; continue;} 284289bc588SBarry Smith for ( i=0; i<m; i++ ) { 285d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v++; continue;} 286289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 287289bc588SBarry Smith } 288289bc588SBarry Smith } 289289bc588SBarry Smith } 290289bc588SBarry Smith else { 291289bc588SBarry Smith for ( j=0; j<n; j++ ) { 292d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v += m; continue;} 293289bc588SBarry Smith for ( i=0; i<m; i++ ) { 294d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v++; continue;} 295289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 296289bc588SBarry Smith } 297289bc588SBarry Smith } 298289bc588SBarry Smith } 299e8d4e0b9SBarry Smith } 300e8d4e0b9SBarry Smith else { 301edae2e7dSBarry Smith if (addv == INSERTVALUES) { 302e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 303d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v += n; continue;} 304e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 305d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v++; continue;} 306e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 307e8d4e0b9SBarry Smith } 308e8d4e0b9SBarry Smith } 309e8d4e0b9SBarry Smith } 310289bc588SBarry Smith else { 311289bc588SBarry Smith for ( i=0; i<m; i++ ) { 312d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v += n; continue;} 313289bc588SBarry Smith for ( j=0; j<n; j++ ) { 314d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v++; continue;} 315289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 316289bc588SBarry Smith } 317289bc588SBarry Smith } 318289bc588SBarry Smith } 319e8d4e0b9SBarry Smith } 320289bc588SBarry Smith return 0; 321289bc588SBarry Smith } 322e8d4e0b9SBarry Smith 323289bc588SBarry Smith /* -----------------------------------------------------------------*/ 32407a0e7adSLois Curfman McInnes static int MatCopy_Dense_Private(Mat matin,Mat *newmat) 325289bc588SBarry Smith { 32644a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 327289bc588SBarry Smith int ierr; 328289bc588SBarry Smith Mat newi; 32944a69424SLois Curfman McInnes Mat_Dense *l; 3306b5873e3SBarry Smith if ((ierr = MatCreateSequentialDense(matin->comm,mat->m,mat->n,&newi))) 3316b5873e3SBarry Smith SETERR(ierr,0); 33244a69424SLois Curfman McInnes l = (Mat_Dense *) newi->data; 333289bc588SBarry Smith MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 334289bc588SBarry Smith *newmat = newi; 335289bc588SBarry Smith return 0; 336289bc588SBarry Smith } 337da3a660dSBarry Smith #include "viewer.h" 338289bc588SBarry Smith 33944a69424SLois Curfman McInnes int MatView_Dense(PetscObject obj,Viewer ptr) 340289bc588SBarry Smith { 341289bc588SBarry Smith Mat matin = (Mat) obj; 34244a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 343289bc588SBarry Smith Scalar *v; 344289bc588SBarry Smith int i,j; 34523242f5aSBarry Smith PetscObject vobj = (PetscObject) ptr; 346da3a660dSBarry Smith 34723242f5aSBarry Smith if (ptr == 0) { 34823242f5aSBarry Smith ptr = STDOUT_VIEWER; vobj = (PetscObject) ptr; 34923242f5aSBarry Smith } 35023242f5aSBarry Smith if (vobj->cookie == VIEWER_COOKIE && vobj->type == MATLAB_VIEWER) { 3516f75cfa4SBarry Smith return ViewerMatlabPutArray_Private(ptr,mat->m,mat->n,mat->v); 352da3a660dSBarry Smith } 353da3a660dSBarry Smith else { 3544a0ce102SLois Curfman McInnes FILE *fd = ViewerFileGetPointer_Private(ptr); 355289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 356289bc588SBarry Smith v = mat->v + i; 357289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 358289bc588SBarry Smith #if defined(PETSC_COMPLEX) 3598e37dc09SBarry Smith fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 360289bc588SBarry Smith #else 3618e37dc09SBarry Smith fprintf(fd,"%6.4e ",*v); v += mat->m; 362289bc588SBarry Smith #endif 363289bc588SBarry Smith } 3648e37dc09SBarry Smith fprintf(fd,"\n"); 365289bc588SBarry Smith } 366da3a660dSBarry Smith } 367289bc588SBarry Smith return 0; 368289bc588SBarry Smith } 369289bc588SBarry Smith 370289bc588SBarry Smith 37144a69424SLois Curfman McInnes static int MatDestroy_Dense(PetscObject obj) 372289bc588SBarry Smith { 373289bc588SBarry Smith Mat mat = (Mat) obj; 37444a69424SLois Curfman McInnes Mat_Dense *l = (Mat_Dense *) mat->data; 375a5a9c739SBarry Smith #if defined(PETSC_LOG) 376a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 377a5a9c739SBarry Smith #endif 378289bc588SBarry Smith if (l->pivots) FREE(l->pivots); 379289bc588SBarry Smith FREE(l); 380a5a9c739SBarry Smith PLogObjectDestroy(mat); 3819e25ed09SBarry Smith PETSCHEADERDESTROY(mat); 382289bc588SBarry Smith return 0; 383289bc588SBarry Smith } 384289bc588SBarry Smith 38544a69424SLois Curfman McInnes static int MatTrans_Dense(Mat matin) 386289bc588SBarry Smith { 38744a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 388289bc588SBarry Smith int k,j; 389289bc588SBarry Smith Scalar *v = mat->v, tmp; 390289bc588SBarry Smith if (mat->m != mat->n) { 391289bc588SBarry Smith SETERR(1,"Cannot transpose rectangular dense matrix"); 392289bc588SBarry Smith } 393289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 394289bc588SBarry Smith for ( k=0; k<j; k++ ) { 395289bc588SBarry Smith tmp = v[j + k*mat->n]; 396289bc588SBarry Smith v[j + k*mat->n] = v[k + j*mat->n]; 397289bc588SBarry Smith v[k + j*mat->n] = tmp; 398289bc588SBarry Smith } 399289bc588SBarry Smith } 400289bc588SBarry Smith return 0; 401289bc588SBarry Smith } 402289bc588SBarry Smith 40344a69424SLois Curfman McInnes static int MatEqual_Dense(Mat matin1,Mat matin2) 404289bc588SBarry Smith { 40544a69424SLois Curfman McInnes Mat_Dense *mat1 = (Mat_Dense *) matin1->data; 40644a69424SLois Curfman McInnes Mat_Dense *mat2 = (Mat_Dense *) matin2->data; 407289bc588SBarry Smith int i; 408289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 409289bc588SBarry Smith if (mat1->m != mat2->m) return 0; 410289bc588SBarry Smith if (mat1->n != mat2->n) return 0; 411289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 412289bc588SBarry Smith if (*v1 != *v2) return 0; 413289bc588SBarry Smith v1++; v2++; 414289bc588SBarry Smith } 415289bc588SBarry Smith return 1; 416289bc588SBarry Smith } 417289bc588SBarry Smith 41846fc4a8fSLois Curfman McInnes static int MatGetDiagonal_Dense(Mat matin,Vec v) 419289bc588SBarry Smith { 42044a69424SLois Curfman McInnes Mat_Dense *mat = (Mat_Dense *) matin->data; 4216abc6512SBarry Smith int i, n; 4226abc6512SBarry Smith Scalar *x; 423289bc588SBarry Smith CHKTYPE(v,SEQVECTOR); 424289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 425289bc588SBarry Smith if (n != mat->m) SETERR(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); 439da3a660dSBarry Smith if (m != mat->m) SETERR(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); 448da3a660dSBarry Smith if (n != mat->n) SETERR(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 { 505289bc588SBarry Smith SETERR(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; 5226f0a148fSBarry Smith MEMSET(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; 531260925b8SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERR(ierr); 532260925b8SBarry Smith ierr = ISGetIndices(is,&rows); CHKERR(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 { 5640754003eSLois Curfman McInnes SETERR(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; 5710754003eSLois Curfman McInnes int *irow, *icol, nrows, ncols, *cwork, jstart; 5720754003eSLois Curfman McInnes Scalar *vwork, *val; 5730754003eSLois Curfman McInnes Mat newmat; 5740754003eSLois Curfman McInnes 5750754003eSLois Curfman McInnes ierr = ISGetIndices(isrow,&irow); CHKERR(ierr); 5760754003eSLois Curfman McInnes ierr = ISGetIndices(iscol,&icol); CHKERR(ierr); 5770754003eSLois Curfman McInnes ierr = ISGetSize(isrow,&nrows); CHKERR(ierr); 5780754003eSLois Curfman McInnes ierr = ISGetSize(iscol,&ncols); CHKERR(ierr); 5790754003eSLois Curfman McInnes 5800754003eSLois Curfman McInnes smap = (int *) MALLOC(oldcols*sizeof(int)); CHKPTR(smap); 5810754003eSLois Curfman McInnes cwork = (int *) MALLOC(ncols*sizeof(int)); CHKPTR(cwork); 5820754003eSLois Curfman McInnes vwork = (Scalar *) MALLOC(ncols*sizeof(Scalar)); CHKPTR(vwork); 5830754003eSLois Curfman McInnes memset((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); 5880754003eSLois Curfman McInnes CHKERR(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); 5990754003eSLois Curfman McInnes CHKERR(ierr); 6000754003eSLois Curfman McInnes } 601ee50ffe9SBarry Smith ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERR(ierr); 602ee50ffe9SBarry Smith ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERR(ierr); 6030754003eSLois Curfman McInnes 6040754003eSLois Curfman McInnes /* Free work space */ 6050754003eSLois Curfman McInnes FREE(smap); FREE(cwork); FREE(vwork); 6060754003eSLois Curfman McInnes ierr = ISRestoreIndices(isrow,&irow); CHKERR(ierr); 6070754003eSLois Curfman McInnes ierr = ISRestoreIndices(iscol,&icol); CHKERR(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, 62144a69424SLois Curfman McInnes MatTrans_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, 63107a0e7adSLois Curfman McInnes MatCopy_Dense_Private}; 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 646*dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 647d65003e9SLois Curfman McInnes 648*dbd7a890SLois 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); 65844a69424SLois Curfman McInnes l = (Mat_Dense *) MALLOC(size); CHKPTR(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 671289bc588SBarry Smith MEMSET(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