1cb512458SBarry Smith #ifndef lint 2*9e25ed09SBarry Smith static char vcid[] = "$Id: dense.c,v 1.11 1995/03/10 04:44:45 bsmith 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 */ 19289bc588SBarry Smith } MatiSD; 20289bc588SBarry Smith 21289bc588SBarry Smith 22289bc588SBarry Smith static int MatiSDnz(Mat matin,int *nz) 23289bc588SBarry Smith { 24289bc588SBarry Smith MatiSD *mat = (MatiSD *) 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++;} 28289bc588SBarry Smith *nz = count; return 0; 29289bc588SBarry Smith } 30289bc588SBarry Smith static int MatiSDmemory(Mat matin,int *mem) 31289bc588SBarry Smith { 32289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 33289bc588SBarry Smith *mem = mat->m*mat->n*sizeof(Scalar); return 0; 34289bc588SBarry Smith } 35289bc588SBarry Smith 36289bc588SBarry Smith /* ---------------------------------------------------------------*/ 37289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots, 38289bc588SBarry Smith rather than put it in the Mat->col slot.*/ 39289bc588SBarry Smith static int MatiSDlufactor(Mat matin,IS row,IS col) 40289bc588SBarry Smith { 41289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 426abc6512SBarry Smith int info; 43289bc588SBarry Smith if (!mat->pivots) { 44289bc588SBarry Smith mat->pivots = (int *) MALLOC( mat->m*sizeof(int) ); 45289bc588SBarry Smith CHKPTR(mat->pivots); 46289bc588SBarry Smith } 47289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 48289bc588SBarry Smith if (info) SETERR(1,"Bad LU factorization"); 49289bc588SBarry Smith matin->factor = FACTOR_LU; 50289bc588SBarry Smith return 0; 51289bc588SBarry Smith } 52289bc588SBarry Smith static int MatiSDlufactorsymbolic(Mat matin,IS row,IS col,Mat *fact) 53289bc588SBarry Smith { 54289bc588SBarry Smith int ierr; 556abc6512SBarry Smith if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0); 56289bc588SBarry Smith return 0; 57289bc588SBarry Smith } 5820563c6bSBarry Smith static int MatiSDlufactornumeric(Mat matin,Mat *fact) 59289bc588SBarry Smith { 6020563c6bSBarry Smith return MatLUFactor(*fact,0,0); 61289bc588SBarry Smith } 62289bc588SBarry Smith static int MatiSDchfactorsymbolic(Mat matin,IS row,Mat *fact) 63289bc588SBarry Smith { 64289bc588SBarry Smith int ierr; 656abc6512SBarry Smith if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0); 66289bc588SBarry Smith return 0; 67289bc588SBarry Smith } 6820563c6bSBarry Smith static int MatiSDchfactornumeric(Mat matin,Mat *fact) 69289bc588SBarry Smith { 7020563c6bSBarry Smith return MatCholeskyFactor(*fact,0); 71289bc588SBarry Smith } 72289bc588SBarry Smith static int MatiSDchfactor(Mat matin,IS perm) 73289bc588SBarry Smith { 74289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 756abc6512SBarry Smith int info; 76289bc588SBarry Smith if (mat->pivots) {FREE(mat->pivots); mat->pivots = 0;} 77289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 78289bc588SBarry Smith if (info) SETERR(1,"Bad Cholesky factorization"); 79289bc588SBarry Smith matin->factor = FACTOR_CHOLESKY; 80289bc588SBarry Smith return 0; 81289bc588SBarry Smith } 82289bc588SBarry Smith 83289bc588SBarry Smith static int MatiSDsolve(Mat matin,Vec xx,Vec yy) 84289bc588SBarry Smith { 85289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 866abc6512SBarry Smith int one = 1, info; 876abc6512SBarry Smith Scalar *x, *y; 88289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 89289bc588SBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 90*9e25ed09SBarry 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 } 94*9e25ed09SBarry 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 } 98*9e25ed09SBarry Smith else SETERR(1,"Matrix must be factored to solve"); 99289bc588SBarry Smith if (info) SETERR(1,"Bad solve"); 100289bc588SBarry Smith return 0; 101289bc588SBarry Smith } 102289bc588SBarry Smith static int MatiSDsolvetrans(Mat matin,Vec xx,Vec yy) 103da3a660dSBarry Smith { 104da3a660dSBarry Smith MatiSD *mat = (MatiSD *) matin->data; 1056abc6512SBarry Smith int one = 1, info; 1066abc6512SBarry Smith Scalar *x, *y; 107da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 108da3a660dSBarry Smith MEMCPY(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 } 118da3a660dSBarry Smith if (info) SETERR(1,"Bad solve"); 119da3a660dSBarry Smith return 0; 120da3a660dSBarry Smith } 121da3a660dSBarry Smith static int MatiSDsolveadd(Mat matin,Vec xx,Vec zz,Vec yy) 122da3a660dSBarry Smith { 123da3a660dSBarry Smith MatiSD *mat = (MatiSD *) 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) { 129da3a660dSBarry Smith ierr = VecCreate(yy,&tmp); CHKERR(ierr); 130da3a660dSBarry Smith ierr = VecCopy(yy,tmp); CHKERR(ierr); 131da3a660dSBarry Smith } 132da3a660dSBarry Smith MEMCPY(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 } 142da3a660dSBarry Smith if (info) SETERR(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 } 147da3a660dSBarry Smith static int MatiSDsolvetransadd(Mat matin,Vec xx,Vec zz, Vec yy) 148da3a660dSBarry Smith { 149da3a660dSBarry Smith MatiSD *mat = (MatiSD *) 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) { 155da3a660dSBarry Smith ierr = VecCreate(yy,&tmp); CHKERR(ierr); 156da3a660dSBarry Smith ierr = VecCopy(yy,tmp); CHKERR(ierr); 157da3a660dSBarry Smith } 158da3a660dSBarry Smith MEMCPY(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 } 168da3a660dSBarry Smith if (info) SETERR(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 /* ------------------------------------------------------------------*/ 174d6dfbf8fSBarry Smith static int MatiSDrelax(Mat matin,Vec bb,double omega,int flag,double shift, 175289bc588SBarry Smith int its,Vec xx) 176289bc588SBarry Smith { 177289bc588SBarry Smith MatiSD *mat = (MatiSD *) 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 */ 1846abc6512SBarry Smith if ((ierr = VecSet(&zero,xx))) SETERR(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 /* -----------------------------------------------------------------*/ 205289bc588SBarry Smith static int MatiSDmulttrans(Mat matin,Vec xx,Vec yy) 206289bc588SBarry Smith { 207289bc588SBarry Smith MatiSD *mat = (MatiSD *) 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 } 215289bc588SBarry Smith static int MatiSDmult(Mat matin,Vec xx,Vec yy) 216289bc588SBarry Smith { 217289bc588SBarry Smith MatiSD *mat = (MatiSD *) 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 } 225289bc588SBarry Smith static int MatiSDmultadd(Mat matin,Vec xx,Vec zz,Vec yy) 226289bc588SBarry Smith { 227289bc588SBarry Smith MatiSD *mat = (MatiSD *) 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); 231289bc588SBarry Smith if (zz != yy) MEMCPY(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 } 236289bc588SBarry Smith static int MatiSDmulttransadd(Mat matin,Vec xx,Vec zz,Vec yy) 237289bc588SBarry Smith { 238289bc588SBarry Smith MatiSD *mat = (MatiSD *) 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); 244289bc588SBarry Smith if (zz != yy) MEMCPY(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 /* -----------------------------------------------------------------*/ 251289bc588SBarry Smith static int MatiSDgetrow(Mat matin,int row,int *ncols,int **cols, 252289bc588SBarry Smith Scalar **vals) 253289bc588SBarry Smith { 254289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 255289bc588SBarry Smith Scalar *v; 256289bc588SBarry Smith int i; 257289bc588SBarry Smith *ncols = mat->n; 258289bc588SBarry Smith if (cols) { 259289bc588SBarry Smith *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols); 260289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) *cols[i] = i; 261289bc588SBarry Smith } 262289bc588SBarry Smith if (vals) { 263289bc588SBarry Smith *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*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 } 269289bc588SBarry Smith static int MatiSDrestorerow(Mat matin,int row,int *ncols,int **cols, 270289bc588SBarry Smith Scalar **vals) 271289bc588SBarry Smith { 272289bc588SBarry Smith if (cols) { FREE(*cols); } 273289bc588SBarry Smith if (vals) { FREE(*vals); } 274289bc588SBarry Smith return 0; 275289bc588SBarry Smith } 276289bc588SBarry Smith /* ----------------------------------------------------------------*/ 277e8d4e0b9SBarry Smith static int MatiSDinsert(Mat matin,int m,int *indexm,int n, 278e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 279289bc588SBarry Smith { 280289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 281289bc588SBarry Smith int i,j; 282d6dfbf8fSBarry Smith 283289bc588SBarry Smith if (!mat->roworiented) { 284e8d4e0b9SBarry Smith if (addv == InsertValues) { 285289bc588SBarry Smith for ( j=0; j<n; j++ ) { 286d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v += m; continue;} 287289bc588SBarry Smith for ( i=0; i<m; i++ ) { 288d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v++; continue;} 289289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 290289bc588SBarry Smith } 291289bc588SBarry Smith } 292289bc588SBarry Smith } 293289bc588SBarry Smith else { 294289bc588SBarry Smith for ( j=0; j<n; j++ ) { 295d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v += m; continue;} 296289bc588SBarry Smith for ( i=0; i<m; i++ ) { 297d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v++; continue;} 298289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 299289bc588SBarry Smith } 300289bc588SBarry Smith } 301289bc588SBarry Smith } 302e8d4e0b9SBarry Smith } 303e8d4e0b9SBarry Smith else { 304e8d4e0b9SBarry Smith if (addv == InsertValues) { 305e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 306d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v += n; continue;} 307e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 308d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v++; continue;} 309e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 310e8d4e0b9SBarry Smith } 311e8d4e0b9SBarry Smith } 312e8d4e0b9SBarry Smith } 313289bc588SBarry Smith else { 314289bc588SBarry Smith for ( i=0; i<m; i++ ) { 315d6dfbf8fSBarry Smith if (indexm[i] < 0) {*v += n; continue;} 316289bc588SBarry Smith for ( j=0; j<n; j++ ) { 317d6dfbf8fSBarry Smith if (indexn[j] < 0) {*v++; continue;} 318289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 319289bc588SBarry Smith } 320289bc588SBarry Smith } 321289bc588SBarry Smith } 322e8d4e0b9SBarry Smith } 323289bc588SBarry Smith return 0; 324289bc588SBarry Smith } 325e8d4e0b9SBarry Smith 326289bc588SBarry Smith /* -----------------------------------------------------------------*/ 327289bc588SBarry Smith static int MatiSDcopy(Mat matin,Mat *newmat) 328289bc588SBarry Smith { 329289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 330289bc588SBarry Smith int ierr; 331289bc588SBarry Smith Mat newi; 332289bc588SBarry Smith MatiSD *l; 3336abc6512SBarry Smith if ((ierr = MatCreateSequentialDense(mat->m,mat->n,&newi))) SETERR(ierr,0); 334289bc588SBarry Smith l = (MatiSD *) newi->data; 335289bc588SBarry Smith MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 336289bc588SBarry Smith *newmat = newi; 337289bc588SBarry Smith return 0; 338289bc588SBarry Smith } 339da3a660dSBarry Smith #include "viewer.h" 340289bc588SBarry Smith 341289bc588SBarry Smith int MatiSDview(PetscObject obj,Viewer ptr) 342289bc588SBarry Smith { 343289bc588SBarry Smith Mat matin = (Mat) obj; 344289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 345289bc588SBarry Smith Scalar *v; 346289bc588SBarry Smith int i,j; 347da3a660dSBarry Smith PetscObject ojb = (PetscObject) ptr; 348da3a660dSBarry Smith 3496abc6512SBarry Smith if (ojb && ojb->cookie == VIEWER_COOKIE && ojb->type == MATLAB_VIEWER) { 350da3a660dSBarry Smith return ViewerMatlabPutArray(ptr,mat->m,mat->n,mat->v); 351da3a660dSBarry Smith } 352da3a660dSBarry Smith else { 353289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 354289bc588SBarry Smith v = mat->v + i; 355289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 356289bc588SBarry Smith #if defined(PETSC_COMPLEX) 357289bc588SBarry Smith printf("%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 358289bc588SBarry Smith #else 359289bc588SBarry Smith printf("%6.4e ",*v); v += mat->m; 360289bc588SBarry Smith #endif 361289bc588SBarry Smith } 362289bc588SBarry Smith printf("\n"); 363289bc588SBarry Smith } 364da3a660dSBarry Smith } 365289bc588SBarry Smith return 0; 366289bc588SBarry Smith } 367289bc588SBarry Smith 368289bc588SBarry Smith 369289bc588SBarry Smith static int MatiSDdestroy(PetscObject obj) 370289bc588SBarry Smith { 371289bc588SBarry Smith Mat mat = (Mat) obj; 372289bc588SBarry Smith MatiSD *l = (MatiSD *) mat->data; 373289bc588SBarry Smith if (l->pivots) FREE(l->pivots); 374289bc588SBarry Smith FREE(l); 375*9e25ed09SBarry Smith PETSCHEADERDESTROY(mat); 376289bc588SBarry Smith return 0; 377289bc588SBarry Smith } 378289bc588SBarry Smith 379289bc588SBarry Smith static int MatiSDtrans(Mat matin) 380289bc588SBarry Smith { 381289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 382289bc588SBarry Smith int k,j; 383289bc588SBarry Smith Scalar *v = mat->v, tmp; 384289bc588SBarry Smith if (mat->m != mat->n) { 385289bc588SBarry Smith SETERR(1,"Cannot transpose rectangular dense matrix"); 386289bc588SBarry Smith } 387289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 388289bc588SBarry Smith for ( k=0; k<j; k++ ) { 389289bc588SBarry Smith tmp = v[j + k*mat->n]; 390289bc588SBarry Smith v[j + k*mat->n] = v[k + j*mat->n]; 391289bc588SBarry Smith v[k + j*mat->n] = tmp; 392289bc588SBarry Smith } 393289bc588SBarry Smith } 394289bc588SBarry Smith return 0; 395289bc588SBarry Smith } 396289bc588SBarry Smith 397289bc588SBarry Smith static int MatiSDequal(Mat matin1,Mat matin2) 398289bc588SBarry Smith { 399289bc588SBarry Smith MatiSD *mat1 = (MatiSD *) matin1->data; 400289bc588SBarry Smith MatiSD *mat2 = (MatiSD *) matin2->data; 401289bc588SBarry Smith int i; 402289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 403289bc588SBarry Smith if (mat1->m != mat2->m) return 0; 404289bc588SBarry Smith if (mat1->n != mat2->n) return 0; 405289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 406289bc588SBarry Smith if (*v1 != *v2) return 0; 407289bc588SBarry Smith v1++; v2++; 408289bc588SBarry Smith } 409289bc588SBarry Smith return 1; 410289bc588SBarry Smith } 411289bc588SBarry Smith 412289bc588SBarry Smith static int MatiSDgetdiag(Mat matin,Vec v) 413289bc588SBarry Smith { 414289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 4156abc6512SBarry Smith int i, n; 4166abc6512SBarry Smith Scalar *x; 417289bc588SBarry Smith CHKTYPE(v,SEQVECTOR); 418289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 419289bc588SBarry Smith if (n != mat->m) SETERR(1,"Nonconforming matrix and vector"); 420289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 421289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 422289bc588SBarry Smith } 423289bc588SBarry Smith return 0; 424289bc588SBarry Smith } 425289bc588SBarry Smith 426da3a660dSBarry Smith static int MatiSDscale(Mat matin,Vec ll,Vec rr) 427289bc588SBarry Smith { 428da3a660dSBarry Smith MatiSD *mat = (MatiSD *) matin->data; 429da3a660dSBarry Smith Scalar *l,*r,x,*v; 430da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 43128988994SBarry Smith if (ll) { 432da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 433da3a660dSBarry Smith if (m != mat->m) SETERR(1,"Left scaling vector wrong length"); 434da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 435da3a660dSBarry Smith x = l[i]; 436da3a660dSBarry Smith v = mat->v + i; 437da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 438da3a660dSBarry Smith } 439da3a660dSBarry Smith } 44028988994SBarry Smith if (rr) { 441da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 442da3a660dSBarry Smith if (n != mat->n) SETERR(1,"Right scaling vector wrong length"); 443da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 444da3a660dSBarry Smith x = r[i]; 445da3a660dSBarry Smith v = mat->v + i*m; 446da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 447da3a660dSBarry Smith } 448da3a660dSBarry Smith } 449289bc588SBarry Smith return 0; 450289bc588SBarry Smith } 451289bc588SBarry Smith 452da3a660dSBarry Smith 453289bc588SBarry Smith static int MatiSDnorm(Mat matin,int type,double *norm) 454289bc588SBarry Smith { 455289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 456289bc588SBarry Smith Scalar *v = mat->v; 457289bc588SBarry Smith double sum = 0.0; 458289bc588SBarry Smith int i, j; 459289bc588SBarry Smith if (type == NORM_FROBENIUS) { 460289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 461289bc588SBarry Smith #if defined(PETSC_COMPLEX) 462289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 463289bc588SBarry Smith #else 464289bc588SBarry Smith sum += (*v)*(*v); v++; 465289bc588SBarry Smith #endif 466289bc588SBarry Smith } 467289bc588SBarry Smith *norm = sqrt(sum); 468289bc588SBarry Smith } 469289bc588SBarry Smith else if (type == NORM_1) { 470289bc588SBarry Smith *norm = 0.0; 471289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 472289bc588SBarry Smith sum = 0.0; 473289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 474289bc588SBarry Smith #if defined(PETSC_COMPLEX) 475289bc588SBarry Smith sum += abs(*v++); 476289bc588SBarry Smith #else 477289bc588SBarry Smith sum += fabs(*v++); 478289bc588SBarry Smith #endif 479289bc588SBarry Smith } 480289bc588SBarry Smith if (sum > *norm) *norm = sum; 481289bc588SBarry Smith } 482289bc588SBarry Smith } 483289bc588SBarry Smith else if (type == NORM_INFINITY) { 484289bc588SBarry Smith *norm = 0.0; 485289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 486289bc588SBarry Smith v = mat->v + j; 487289bc588SBarry Smith sum = 0.0; 488289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 489289bc588SBarry Smith #if defined(PETSC_COMPLEX) 490289bc588SBarry Smith sum += abs(*v); v += mat->m; 491289bc588SBarry Smith #else 492289bc588SBarry Smith sum += fabs(*v); v += mat->m; 493289bc588SBarry Smith #endif 494289bc588SBarry Smith } 495289bc588SBarry Smith if (sum > *norm) *norm = sum; 496289bc588SBarry Smith } 497289bc588SBarry Smith } 498289bc588SBarry Smith else { 499289bc588SBarry Smith SETERR(1,"No support for the two norm yet"); 500289bc588SBarry Smith } 501289bc588SBarry Smith return 0; 502289bc588SBarry Smith } 503289bc588SBarry Smith 504289bc588SBarry Smith static int MatiDenseinsopt(Mat aijin,int op) 505289bc588SBarry Smith { 506289bc588SBarry Smith MatiSD *aij = (MatiSD *) aijin->data; 507289bc588SBarry Smith if (op == ROW_ORIENTED) aij->roworiented = 1; 508289bc588SBarry Smith else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 509289bc588SBarry Smith /* doesn't care about sorted rows or columns */ 510289bc588SBarry Smith return 0; 511289bc588SBarry Smith } 512289bc588SBarry Smith 5136f0a148fSBarry Smith static int MatiZero(Mat A) 5146f0a148fSBarry Smith { 5156f0a148fSBarry Smith MatiSD *l = (MatiSD *) A->data; 5166f0a148fSBarry Smith MEMSET(l->v,0,l->m*l->n*sizeof(Scalar)); 5176f0a148fSBarry Smith return 0; 5186f0a148fSBarry Smith } 5196f0a148fSBarry Smith 520260925b8SBarry Smith static int MatiZerorows(Mat A,IS is,Scalar *diag) 5216f0a148fSBarry Smith { 5226f0a148fSBarry Smith MatiSD *l = (MatiSD *) A->data; 5236abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 5246f0a148fSBarry Smith Scalar *slot; 525260925b8SBarry Smith ierr = ISGetLocalSize(is,&N); CHKERR(ierr); 526260925b8SBarry Smith ierr = ISGetIndices(is,&rows); CHKERR(ierr); 5276f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5286f0a148fSBarry Smith slot = l->v + rows[i]; 5296f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 5306f0a148fSBarry Smith } 5316f0a148fSBarry Smith if (diag) { 5326f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 5336f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 534260925b8SBarry Smith *slot = *diag; 5356f0a148fSBarry Smith } 5366f0a148fSBarry Smith } 537260925b8SBarry Smith ISRestoreIndices(is,&rows); 5386f0a148fSBarry Smith return 0; 5396f0a148fSBarry Smith } 540289bc588SBarry Smith /* -------------------------------------------------------------------*/ 541e8d4e0b9SBarry Smith static struct _MatOps MatOps = {MatiSDinsert, 542289bc588SBarry Smith MatiSDgetrow, MatiSDrestorerow, 543289bc588SBarry Smith MatiSDmult, MatiSDmultadd, MatiSDmulttrans, MatiSDmulttransadd, 544da3a660dSBarry Smith MatiSDsolve,MatiSDsolveadd,MatiSDsolvetrans,MatiSDsolvetransadd, 545289bc588SBarry Smith MatiSDlufactor,MatiSDchfactor, 546289bc588SBarry Smith MatiSDrelax, 547289bc588SBarry Smith MatiSDtrans, 548289bc588SBarry Smith MatiSDnz,MatiSDmemory,MatiSDequal, 549289bc588SBarry Smith MatiSDcopy, 550289bc588SBarry Smith MatiSDgetdiag,MatiSDscale,MatiSDnorm, 551289bc588SBarry Smith 0,0, 5526f0a148fSBarry Smith 0, MatiDenseinsopt,MatiZero,MatiZerorows,0, 553289bc588SBarry Smith MatiSDlufactorsymbolic,MatiSDlufactornumeric, 554289bc588SBarry Smith MatiSDchfactorsymbolic,MatiSDchfactornumeric 555289bc588SBarry Smith }; 55620563c6bSBarry Smith /*@ 55720563c6bSBarry Smith MatCreateSequentialDense - Creates a sequential dense matrix that 55820563c6bSBarry Smith is stored in the usual Fortran 77 manner. Many of the matrix 55920563c6bSBarry Smith operations use the BLAS and LAPACK routines. 560289bc588SBarry Smith 56120563c6bSBarry Smith Input Parameters: 56220563c6bSBarry Smith . m, n - the number of rows and columns in the matrix. 56320563c6bSBarry Smith 56420563c6bSBarry Smith Output Parameter: 56520563c6bSBarry Smith . newmat - the matrix created. 56620563c6bSBarry Smith 56720563c6bSBarry Smith Keywords: dense matrix, lapack, blas 56820563c6bSBarry Smith @*/ 569289bc588SBarry Smith int MatCreateSequentialDense(int m,int n,Mat *newmat) 570289bc588SBarry Smith { 571289bc588SBarry Smith int size = sizeof(MatiSD) + m*n*sizeof(Scalar); 572289bc588SBarry Smith Mat mat; 573289bc588SBarry Smith MatiSD *l; 574289bc588SBarry Smith *newmat = 0; 575*9e25ed09SBarry Smith PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSESEQ,MPI_COMM_SELF); 576289bc588SBarry Smith l = (MatiSD *) MALLOC(size); CHKPTR(l); 577289bc588SBarry Smith mat->ops = &MatOps; 578289bc588SBarry Smith mat->destroy = MatiSDdestroy; 579289bc588SBarry Smith mat->view = MatiSDview; 580289bc588SBarry Smith mat->data = (void *) l; 581289bc588SBarry Smith mat->factor = 0; 582289bc588SBarry Smith mat->col = 0; 583289bc588SBarry Smith mat->row = 0; 584d6dfbf8fSBarry Smith 585289bc588SBarry Smith l->m = m; 586289bc588SBarry Smith l->n = n; 587289bc588SBarry Smith l->v = (Scalar *) (l + 1); 588289bc588SBarry Smith l->pivots = 0; 589289bc588SBarry Smith l->roworiented = 1; 590d6dfbf8fSBarry Smith 591289bc588SBarry Smith MEMSET(l->v,0,m*n*sizeof(Scalar)); 592289bc588SBarry Smith *newmat = mat; 593289bc588SBarry Smith return 0; 594289bc588SBarry Smith } 595289bc588SBarry Smith 596289bc588SBarry Smith int MatiSDCreate(Mat matin,Mat *newmat) 597289bc588SBarry Smith { 598289bc588SBarry Smith MatiSD *m = (MatiSD *) matin->data; 599289bc588SBarry Smith return MatCreateSequentialDense(m->m,m->n,newmat); 600289bc588SBarry Smith } 601