1289bc588SBarry Smith 2289bc588SBarry Smith /* 3289bc588SBarry Smith Standard Fortran style matrices 4289bc588SBarry Smith */ 5289bc588SBarry Smith #include "ptscimpl.h" 6289bc588SBarry Smith #include "plapack.h" 7289bc588SBarry Smith #include "matimpl.h" 8289bc588SBarry Smith #include "math.h" 9289bc588SBarry Smith #include "vec/vecimpl.h" 10289bc588SBarry Smith 11289bc588SBarry Smith typedef struct { 12289bc588SBarry Smith Scalar *v; 13289bc588SBarry Smith int roworiented; 14289bc588SBarry Smith int m,n,pad; 15289bc588SBarry Smith int *pivots; /* pivots in LU factorization */ 16289bc588SBarry Smith } MatiSD; 17289bc588SBarry Smith 18289bc588SBarry Smith 19289bc588SBarry Smith static int MatiSDnz(Mat matin,int *nz) 20289bc588SBarry Smith { 21289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 22289bc588SBarry Smith int i,N = mat->m*mat->n,count = 0; 23289bc588SBarry Smith Scalar *v = mat->v; 24289bc588SBarry Smith for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 25289bc588SBarry Smith *nz = count; return 0; 26289bc588SBarry Smith } 27289bc588SBarry Smith static int MatiSDmemory(Mat matin,int *mem) 28289bc588SBarry Smith { 29289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 30289bc588SBarry Smith *mem = mat->m*mat->n*sizeof(Scalar); 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.*/ 36289bc588SBarry Smith static int MatiSDlufactor(Mat matin,IS row,IS col) 37289bc588SBarry Smith { 38289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 39289bc588SBarry Smith int ierr, one = 1, 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 } 49289bc588SBarry Smith static int MatiSDlufactorsymbolic(Mat matin,IS row,IS col,Mat *fact) 50289bc588SBarry Smith { 51289bc588SBarry Smith int ierr; 52289bc588SBarry Smith if (ierr = MatCopy(matin,fact)) SETERR(ierr,0); 53289bc588SBarry Smith return 0; 54289bc588SBarry Smith } 55289bc588SBarry Smith static int MatiSDlufactornumeric(Mat matin,Mat fact) 56289bc588SBarry Smith { 57289bc588SBarry Smith return MatLUFactor(fact,0,0); 58289bc588SBarry Smith } 59289bc588SBarry Smith static int MatiSDchfactorsymbolic(Mat matin,IS row,Mat *fact) 60289bc588SBarry Smith { 61289bc588SBarry Smith int ierr; 62289bc588SBarry Smith if (ierr = MatCopy(matin,fact)) SETERR(ierr,0); 63289bc588SBarry Smith return 0; 64289bc588SBarry Smith } 65289bc588SBarry Smith static int MatiSDchfactornumeric(Mat matin,Mat fact) 66289bc588SBarry Smith { 67289bc588SBarry Smith return MatCholeskyFactor(fact,0); 68289bc588SBarry Smith } 69289bc588SBarry Smith static int MatiSDchfactor(Mat matin,IS perm) 70289bc588SBarry Smith { 71289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 72289bc588SBarry Smith int ierr, one = 1, 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 80289bc588SBarry Smith static int MatiSDsolve(Mat matin,Vec xx,Vec yy) 81289bc588SBarry Smith { 82289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 83289bc588SBarry Smith int i,j, one = 1, info; 84289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 85289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 86289bc588SBarry Smith MEMCPY(y,x,mat->m*sizeof(Scalar)); 87289bc588SBarry Smith /* assume if pivots exist then LU else Cholesky */ 88289bc588SBarry Smith if (mat->pivots) { 89289bc588SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 90289bc588SBarry Smith y, &mat->m, &info ); 91289bc588SBarry Smith } 92289bc588SBarry Smith else { 93289bc588SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 94289bc588SBarry Smith y, &mat->m, &info ); 95289bc588SBarry Smith } 96289bc588SBarry Smith if (info) SETERR(1,"Bad solve"); 97289bc588SBarry Smith return 0; 98289bc588SBarry Smith } 99289bc588SBarry Smith static int MatiSDsolvetrans(Mat matin,Vec xx,Vec yy) 100289bc588SBarry Smith {return 0;} 101289bc588SBarry Smith 102289bc588SBarry Smith /* ------------------------------------------------------------------*/ 103289bc588SBarry Smith static int MatiSDrelax(Mat matin,Vec bb,double omega,int flag,IS is, 104289bc588SBarry Smith int its,Vec xx) 105289bc588SBarry Smith { 106289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 107289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 108289bc588SBarry Smith int o = 1, tmp,n = mat->n,ierr, m = mat->m, i, j; 109289bc588SBarry Smith 110289bc588SBarry Smith if (is) SETERR(1,"No support for ordering in relaxation"); 111289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 112289bc588SBarry Smith /* this is a hack fix, should have another version without 113289bc588SBarry Smith the second BLdot */ 114289bc588SBarry Smith if (ierr = VecSet(&zero,xx)) SETERR(ierr,0); 115289bc588SBarry Smith } 116289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 117289bc588SBarry Smith while (its--) { 118289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 119289bc588SBarry Smith for ( i=0; i<m; i++ ) { 120289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 121289bc588SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/v[i + i*m] + x[i]); 122289bc588SBarry Smith } 123289bc588SBarry Smith } 124289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 125289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 126289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 127289bc588SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt/v[i + i*m] + x[i]); 128289bc588SBarry Smith } 129289bc588SBarry Smith } 130289bc588SBarry Smith } 131289bc588SBarry Smith return 0; 132289bc588SBarry Smith } 133289bc588SBarry Smith 134289bc588SBarry Smith /* -----------------------------------------------------------------*/ 135289bc588SBarry Smith static int MatiSDmulttrans(Mat matin,Vec xx,Vec yy) 136289bc588SBarry Smith { 137289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 138289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 139289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 140289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 141289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 142289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 143289bc588SBarry Smith return 0; 144289bc588SBarry Smith } 145289bc588SBarry Smith static int MatiSDmult(Mat matin,Vec xx,Vec yy) 146289bc588SBarry Smith { 147289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 148289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 149289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 150289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 151289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 152289bc588SBarry Smith x, &_One, &_DZero, y, &_One ); 153289bc588SBarry Smith return 0; 154289bc588SBarry Smith } 155289bc588SBarry Smith static int MatiSDmultadd(Mat matin,Vec xx,Vec zz,Vec yy) 156289bc588SBarry Smith { 157289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 158289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 159289bc588SBarry Smith int _One=1; Scalar _DOne=1.0, _DZero=0.0; 160289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 161289bc588SBarry Smith if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 162289bc588SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 163289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 164289bc588SBarry Smith return 0; 165289bc588SBarry Smith } 166289bc588SBarry Smith static int MatiSDmulttransadd(Mat matin,Vec xx,Vec zz,Vec yy) 167289bc588SBarry Smith { 168289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 169289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 170289bc588SBarry Smith int _One=1; 171289bc588SBarry Smith Scalar _DOne=1.0, _DZero=0.0; 172289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 173289bc588SBarry Smith VecGetArray(zz,&z); 174289bc588SBarry Smith if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 175289bc588SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 176289bc588SBarry Smith x, &_One, &_DOne, y, &_One ); 177289bc588SBarry Smith return 0; 178289bc588SBarry Smith } 179289bc588SBarry Smith 180289bc588SBarry Smith /* -----------------------------------------------------------------*/ 181289bc588SBarry Smith static int MatiSDgetrow(Mat matin,int row,int *ncols,int **cols, 182289bc588SBarry Smith Scalar **vals) 183289bc588SBarry Smith { 184289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 185289bc588SBarry Smith Scalar *v; 186289bc588SBarry Smith int i; 187289bc588SBarry Smith *ncols = mat->n; 188289bc588SBarry Smith if (cols) { 189289bc588SBarry Smith *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols); 190289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) *cols[i] = i; 191289bc588SBarry Smith } 192289bc588SBarry Smith if (vals) { 193289bc588SBarry Smith *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals); 194289bc588SBarry Smith v = mat->v + row; 195289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 196289bc588SBarry Smith } 197289bc588SBarry Smith return 0; 198289bc588SBarry Smith } 199289bc588SBarry Smith static int MatiSDrestorerow(Mat matin,int row,int *ncols,int **cols, 200289bc588SBarry Smith Scalar **vals) 201289bc588SBarry Smith { 202289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 203289bc588SBarry Smith if (cols) { FREE(*cols); } 204289bc588SBarry Smith if (vals) { FREE(*vals); } 205289bc588SBarry Smith return 0; 206289bc588SBarry Smith } 207289bc588SBarry Smith /* ----------------------------------------------------------------*/ 208289bc588SBarry Smith static int MatiSDinsert(Mat matin,Scalar *v,int m,int *indexm,int n, 209289bc588SBarry Smith int *indexn) 210289bc588SBarry Smith { 211289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 212289bc588SBarry Smith int i,j; 213289bc588SBarry Smith if (!mat->roworiented) { 214289bc588SBarry Smith for ( j=0; j<n; j++ ) { 215289bc588SBarry Smith for ( i=0; i<m; i++ ) { 216289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 217289bc588SBarry Smith } 218289bc588SBarry Smith } 219289bc588SBarry Smith } 220289bc588SBarry Smith else { 221289bc588SBarry Smith for ( i=0; i<m; i++ ) { 222289bc588SBarry Smith for ( j=0; j<n; j++ ) { 223289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 224289bc588SBarry Smith } 225289bc588SBarry Smith } 226289bc588SBarry Smith } 227289bc588SBarry Smith return 0; 228289bc588SBarry Smith } 229289bc588SBarry Smith static int MatiSDinsertadd(Mat matin,Scalar *v,int m,int *indexm, 230289bc588SBarry Smith int n,int *indexn) 231289bc588SBarry Smith { 232289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 233289bc588SBarry Smith int i,j; 234289bc588SBarry Smith if (!mat->roworiented) { 235289bc588SBarry Smith for ( j=0; j<n; j++ ) { 236289bc588SBarry Smith for ( i=0; i<m; i++ ) { 237289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 238289bc588SBarry Smith } 239289bc588SBarry Smith } 240289bc588SBarry Smith } 241289bc588SBarry Smith else { 242289bc588SBarry Smith for ( i=0; i<m; i++ ) { 243289bc588SBarry Smith for ( j=0; j<n; j++ ) { 244289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 245289bc588SBarry Smith } 246289bc588SBarry Smith } 247289bc588SBarry Smith } 248289bc588SBarry Smith return 0; 249289bc588SBarry Smith } 250289bc588SBarry Smith /* -----------------------------------------------------------------*/ 251289bc588SBarry Smith static int MatiSDcopy(Mat matin,Mat *newmat) 252289bc588SBarry Smith { 253289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 254289bc588SBarry Smith int ierr; 255289bc588SBarry Smith Mat newi; 256289bc588SBarry Smith MatiSD *l; 257289bc588SBarry Smith if (ierr = MatCreateSequentialDense(mat->m,mat->n,&newi)) SETERR(ierr,0); 258289bc588SBarry Smith l = (MatiSD *) newi->data; 259289bc588SBarry Smith MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 260289bc588SBarry Smith *newmat = newi; 261289bc588SBarry Smith return 0; 262289bc588SBarry Smith } 263289bc588SBarry Smith 264289bc588SBarry Smith int MatiSDview(PetscObject obj,Viewer ptr) 265289bc588SBarry Smith { 266289bc588SBarry Smith Mat matin = (Mat) obj; 267289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 268289bc588SBarry Smith Scalar *v; 269289bc588SBarry Smith int i,j; 270289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 271289bc588SBarry Smith v = mat->v + i; 272289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 273289bc588SBarry Smith #if defined(PETSC_COMPLEX) 274289bc588SBarry Smith printf("%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 275289bc588SBarry Smith #else 276289bc588SBarry Smith printf("%6.4e ",*v); v += mat->m; 277289bc588SBarry Smith #endif 278289bc588SBarry Smith } 279289bc588SBarry Smith printf("\n"); 280289bc588SBarry Smith } 281289bc588SBarry Smith return 0; 282289bc588SBarry Smith } 283289bc588SBarry Smith 284289bc588SBarry Smith 285289bc588SBarry Smith static int MatiSDdestroy(PetscObject obj) 286289bc588SBarry Smith { 287289bc588SBarry Smith Mat mat = (Mat) obj; 288289bc588SBarry Smith MatiSD *l = (MatiSD *) mat->data; 289289bc588SBarry Smith if (l->pivots) FREE(l->pivots); 290289bc588SBarry Smith FREE(l); 291289bc588SBarry Smith FREE(mat); 292289bc588SBarry Smith return 0; 293289bc588SBarry Smith } 294289bc588SBarry Smith 295289bc588SBarry Smith static int MatiSDtrans(Mat matin) 296289bc588SBarry Smith { 297289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 298289bc588SBarry Smith int k,j; 299289bc588SBarry Smith Scalar *v = mat->v, tmp; 300289bc588SBarry Smith if (mat->m != mat->n) { 301289bc588SBarry Smith SETERR(1,"Cannot transpose rectangular dense matrix"); 302289bc588SBarry Smith } 303289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 304289bc588SBarry Smith for ( k=0; k<j; k++ ) { 305289bc588SBarry Smith tmp = v[j + k*mat->n]; 306289bc588SBarry Smith v[j + k*mat->n] = v[k + j*mat->n]; 307289bc588SBarry Smith v[k + j*mat->n] = tmp; 308289bc588SBarry Smith } 309289bc588SBarry Smith } 310289bc588SBarry Smith return 0; 311289bc588SBarry Smith } 312289bc588SBarry Smith 313289bc588SBarry Smith static int MatiSDequal(Mat matin1,Mat matin2) 314289bc588SBarry Smith { 315289bc588SBarry Smith MatiSD *mat1 = (MatiSD *) matin1->data; 316289bc588SBarry Smith MatiSD *mat2 = (MatiSD *) matin2->data; 317289bc588SBarry Smith int i; 318289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 319289bc588SBarry Smith if (mat1->m != mat2->m) return 0; 320289bc588SBarry Smith if (mat1->n != mat2->n) return 0; 321289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 322289bc588SBarry Smith if (*v1 != *v2) return 0; 323289bc588SBarry Smith v1++; v2++; 324289bc588SBarry Smith } 325289bc588SBarry Smith return 1; 326289bc588SBarry Smith } 327289bc588SBarry Smith 328289bc588SBarry Smith static int MatiSDgetdiag(Mat matin,Vec v) 329289bc588SBarry Smith { 330289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 331289bc588SBarry Smith int i,j, n; 332289bc588SBarry Smith Scalar *x, zero = 0.0; 333289bc588SBarry Smith CHKTYPE(v,SEQVECTOR); 334289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 335289bc588SBarry Smith if (n != mat->m) SETERR(1,"Nonconforming matrix and vector"); 336289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 337289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 338289bc588SBarry Smith } 339289bc588SBarry Smith return 0; 340289bc588SBarry Smith } 341289bc588SBarry Smith 342289bc588SBarry Smith static int MatiSDscale(Mat matin,Vec l,Vec r) 343289bc588SBarry Smith { 344289bc588SBarry Smith return 0; 345289bc588SBarry Smith } 346289bc588SBarry Smith 347289bc588SBarry Smith static int MatiSDnorm(Mat matin,int type,double *norm) 348289bc588SBarry Smith { 349289bc588SBarry Smith MatiSD *mat = (MatiSD *) matin->data; 350289bc588SBarry Smith Scalar *v = mat->v; 351289bc588SBarry Smith double sum = 0.0; 352289bc588SBarry Smith int i, j; 353289bc588SBarry Smith if (type == NORM_FROBENIUS) { 354289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 355289bc588SBarry Smith #if defined(PETSC_COMPLEX) 356289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 357289bc588SBarry Smith #else 358289bc588SBarry Smith sum += (*v)*(*v); v++; 359289bc588SBarry Smith #endif 360289bc588SBarry Smith } 361289bc588SBarry Smith *norm = sqrt(sum); 362289bc588SBarry Smith } 363289bc588SBarry Smith else if (type == NORM_1) { 364289bc588SBarry Smith *norm = 0.0; 365289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 366289bc588SBarry Smith sum = 0.0; 367289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 368289bc588SBarry Smith #if defined(PETSC_COMPLEX) 369289bc588SBarry Smith sum += abs(*v++); 370289bc588SBarry Smith #else 371289bc588SBarry Smith sum += fabs(*v++); 372289bc588SBarry Smith #endif 373289bc588SBarry Smith } 374289bc588SBarry Smith if (sum > *norm) *norm = sum; 375289bc588SBarry Smith } 376289bc588SBarry Smith } 377289bc588SBarry Smith else if (type == NORM_INFINITY) { 378289bc588SBarry Smith *norm = 0.0; 379289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 380289bc588SBarry Smith v = mat->v + j; 381289bc588SBarry Smith sum = 0.0; 382289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 383289bc588SBarry Smith #if defined(PETSC_COMPLEX) 384289bc588SBarry Smith sum += abs(*v); v += mat->m; 385289bc588SBarry Smith #else 386289bc588SBarry Smith sum += fabs(*v); v += mat->m; 387289bc588SBarry Smith #endif 388289bc588SBarry Smith } 389289bc588SBarry Smith if (sum > *norm) *norm = sum; 390289bc588SBarry Smith } 391289bc588SBarry Smith } 392289bc588SBarry Smith else { 393289bc588SBarry Smith SETERR(1,"No support for the two norm yet"); 394289bc588SBarry Smith } 395289bc588SBarry Smith return 0; 396289bc588SBarry Smith } 397289bc588SBarry Smith 398289bc588SBarry Smith static int MatiDenseinsopt(Mat aijin,int op) 399289bc588SBarry Smith { 400289bc588SBarry Smith MatiSD *aij = (MatiSD *) aijin->data; 401289bc588SBarry Smith if (op == ROW_ORIENTED) aij->roworiented = 1; 402289bc588SBarry Smith else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 403289bc588SBarry Smith /* doesn't care about sorted rows or columns */ 404289bc588SBarry Smith return 0; 405289bc588SBarry Smith } 406289bc588SBarry Smith 407289bc588SBarry Smith /* -------------------------------------------------------------------*/ 408289bc588SBarry Smith static struct _MatOps MatOps = {MatiSDinsert,MatiSDinsert, 409289bc588SBarry Smith MatiSDgetrow, MatiSDrestorerow, 410289bc588SBarry Smith MatiSDmult, MatiSDmultadd, MatiSDmulttrans, MatiSDmulttransadd, 411*8c37ef55SBarry Smith MatiSDsolve,0,MatiSDsolvetrans,0, 412289bc588SBarry Smith MatiSDlufactor,MatiSDchfactor, 413289bc588SBarry Smith MatiSDrelax, 414289bc588SBarry Smith MatiSDtrans, 415289bc588SBarry Smith MatiSDnz,MatiSDmemory,MatiSDequal, 416289bc588SBarry Smith MatiSDcopy, 417289bc588SBarry Smith MatiSDgetdiag,MatiSDscale,MatiSDnorm, 418289bc588SBarry Smith 0,0, 419289bc588SBarry Smith 0, MatiDenseinsopt,0,0,0, 420289bc588SBarry Smith MatiSDlufactorsymbolic,MatiSDlufactornumeric, 421289bc588SBarry Smith MatiSDchfactorsymbolic,MatiSDchfactornumeric 422289bc588SBarry Smith }; 423289bc588SBarry Smith 424289bc588SBarry Smith int MatCreateSequentialDense(int m,int n,Mat *newmat) 425289bc588SBarry Smith { 426289bc588SBarry Smith int size = sizeof(MatiSD) + m*n*sizeof(Scalar); 427289bc588SBarry Smith Mat mat; 428289bc588SBarry Smith MatiSD *l; 429289bc588SBarry Smith *newmat = 0; 430289bc588SBarry Smith CREATEHEADER(mat,_Mat); 431289bc588SBarry Smith l = (MatiSD *) MALLOC(size); CHKPTR(l); 432289bc588SBarry Smith mat->cookie = MAT_COOKIE; 433289bc588SBarry Smith mat->ops = &MatOps; 434289bc588SBarry Smith mat->destroy = MatiSDdestroy; 435289bc588SBarry Smith mat->view = MatiSDview; 436289bc588SBarry Smith mat->data = (void *) l; 437289bc588SBarry Smith mat->type = MATDENSESEQ; 438289bc588SBarry Smith mat->factor = 0; 439289bc588SBarry Smith mat->col = 0; 440289bc588SBarry Smith mat->row = 0; 441289bc588SBarry Smith l->m = m; 442289bc588SBarry Smith l->n = n; 443289bc588SBarry Smith l->v = (Scalar *) (l + 1); 444289bc588SBarry Smith l->pivots = 0; 445289bc588SBarry Smith l->roworiented = 1; 446289bc588SBarry Smith MEMSET(l->v,0,m*n*sizeof(Scalar)); 447289bc588SBarry Smith *newmat = mat; 448289bc588SBarry Smith return 0; 449289bc588SBarry Smith } 450289bc588SBarry Smith 451289bc588SBarry Smith int MatiSDCreate(Mat matin,Mat *newmat) 452289bc588SBarry Smith { 453289bc588SBarry Smith MatiSD *m = (MatiSD *) matin->data; 454289bc588SBarry Smith return MatCreateSequentialDense(m->m,m->n,newmat); 455289bc588SBarry Smith } 456