1 #ifndef lint 2 static char vcid[] = "$Id: dense.c,v 1.18 1995/03/23 22:57:37 curfman Exp bsmith $"; 3 #endif 4 5 /* 6 Standard Fortran style matrices 7 */ 8 #include "ptscimpl.h" 9 #include "plapack.h" 10 #include "matimpl.h" 11 #include "math.h" 12 #include "vec/vecimpl.h" 13 14 typedef struct { 15 Scalar *v; 16 int roworiented; 17 int m,n,pad; 18 int *pivots; /* pivots in LU factorization */ 19 } Mat_Dense; 20 21 22 static int MatNz_Dense(Mat matin,int *nz) 23 { 24 Mat_Dense *mat = (Mat_Dense *) matin->data; 25 int i,N = mat->m*mat->n,count = 0; 26 Scalar *v = mat->v; 27 for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 28 *nz = count; return 0; 29 } 30 static int MatMemory_Dense(Mat matin,int *mem) 31 { 32 Mat_Dense *mat = (Mat_Dense *) matin->data; 33 *mem = mat->m*mat->n*sizeof(Scalar); return 0; 34 } 35 36 /* ---------------------------------------------------------------*/ 37 /* COMMENT: I have chosen to hide column permutation in the pivots, 38 rather than put it in the Mat->col slot.*/ 39 static int MatLUFactor_Dense(Mat matin,IS row,IS col) 40 { 41 Mat_Dense *mat = (Mat_Dense *) matin->data; 42 int info; 43 if (!mat->pivots) { 44 mat->pivots = (int *) MALLOC( mat->m*sizeof(int) ); 45 CHKPTR(mat->pivots); 46 } 47 LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 48 if (info) SETERR(1,"Bad LU factorization"); 49 matin->factor = FACTOR_LU; 50 return 0; 51 } 52 static int MatLUFactorSymbolic_Dense(Mat matin,IS row,IS col,Mat *fact) 53 { 54 int ierr; 55 if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0); 56 return 0; 57 } 58 static int MatLUFactorNumeric_Dense(Mat matin,Mat *fact) 59 { 60 return MatLUFactor(*fact,0,0); 61 } 62 static int MatChFactorSymbolic_Dense(Mat matin,IS row,Mat *fact) 63 { 64 int ierr; 65 if ((ierr = MatCopy(matin,fact))) SETERR(ierr,0); 66 return 0; 67 } 68 static int MatChFactorNumeric_Dense(Mat matin,Mat *fact) 69 { 70 return MatCholeskyFactor(*fact,0); 71 } 72 static int MatChFactor_Dense(Mat matin,IS perm) 73 { 74 Mat_Dense *mat = (Mat_Dense *) matin->data; 75 int info; 76 if (mat->pivots) {FREE(mat->pivots); mat->pivots = 0;} 77 LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 78 if (info) SETERR(1,"Bad Cholesky factorization"); 79 matin->factor = FACTOR_CHOLESKY; 80 return 0; 81 } 82 83 static int MatSolve_Dense(Mat matin,Vec xx,Vec yy) 84 { 85 Mat_Dense *mat = (Mat_Dense *) matin->data; 86 int one = 1, info; 87 Scalar *x, *y; 88 VecGetArray(xx,&x); VecGetArray(yy,&y); 89 MEMCPY(y,x,mat->m*sizeof(Scalar)); 90 if (matin->factor == FACTOR_LU) { 91 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 92 y, &mat->m, &info ); 93 } 94 else if (matin->factor == FACTOR_CHOLESKY){ 95 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 96 y, &mat->m, &info ); 97 } 98 else SETERR(1,"Matrix must be factored to solve"); 99 if (info) SETERR(1,"Bad solve"); 100 return 0; 101 } 102 static int MatSolveTrans_Dense(Mat matin,Vec xx,Vec yy) 103 { 104 Mat_Dense *mat = (Mat_Dense *) matin->data; 105 int one = 1, info; 106 Scalar *x, *y; 107 VecGetArray(xx,&x); VecGetArray(yy,&y); 108 MEMCPY(y,x,mat->m*sizeof(Scalar)); 109 /* assume if pivots exist then LU else Cholesky */ 110 if (mat->pivots) { 111 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 112 y, &mat->m, &info ); 113 } 114 else { 115 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 116 y, &mat->m, &info ); 117 } 118 if (info) SETERR(1,"Bad solve"); 119 return 0; 120 } 121 static int MatSolveAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 122 { 123 Mat_Dense *mat = (Mat_Dense *) matin->data; 124 int one = 1, info,ierr; 125 Scalar *x, *y, sone = 1.0; 126 Vec tmp = 0; 127 VecGetArray(xx,&x); VecGetArray(yy,&y); 128 if (yy == zz) { 129 ierr = VecCreate(yy,&tmp); CHKERR(ierr); 130 ierr = VecCopy(yy,tmp); CHKERR(ierr); 131 } 132 MEMCPY(y,x,mat->m*sizeof(Scalar)); 133 /* assume if pivots exist then LU else Cholesky */ 134 if (mat->pivots) { 135 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 136 y, &mat->m, &info ); 137 } 138 else { 139 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 140 y, &mat->m, &info ); 141 } 142 if (info) SETERR(1,"Bad solve"); 143 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 144 else VecAXPY(&sone,zz,yy); 145 return 0; 146 } 147 static int MatSolveTransAdd_Dense(Mat matin,Vec xx,Vec zz, Vec yy) 148 { 149 Mat_Dense *mat = (Mat_Dense *) matin->data; 150 int one = 1, info,ierr; 151 Scalar *x, *y, sone = 1.0; 152 Vec tmp; 153 VecGetArray(xx,&x); VecGetArray(yy,&y); 154 if (yy == zz) { 155 ierr = VecCreate(yy,&tmp); CHKERR(ierr); 156 ierr = VecCopy(yy,tmp); CHKERR(ierr); 157 } 158 MEMCPY(y,x,mat->m*sizeof(Scalar)); 159 /* assume if pivots exist then LU else Cholesky */ 160 if (mat->pivots) { 161 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots, 162 y, &mat->m, &info ); 163 } 164 else { 165 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 166 y, &mat->m, &info ); 167 } 168 if (info) SETERR(1,"Bad solve"); 169 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 170 else VecAXPY(&sone,zz,yy); 171 return 0; 172 } 173 /* ------------------------------------------------------------------*/ 174 static int MatRelax_Dense(Mat matin,Vec bb,double omega,int flag,double shift, 175 int its,Vec xx) 176 { 177 Mat_Dense *mat = (Mat_Dense *) matin->data; 178 Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 179 int o = 1,ierr, m = mat->m, i; 180 181 if (flag & SOR_ZERO_INITIAL_GUESS) { 182 /* this is a hack fix, should have another version without 183 the second BLdot */ 184 if ((ierr = VecSet(&zero,xx))) SETERR(ierr,0); 185 } 186 VecGetArray(xx,&x); VecGetArray(bb,&b); 187 while (its--) { 188 if (flag & SOR_FORWARD_SWEEP){ 189 for ( i=0; i<m; i++ ) { 190 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 191 x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 192 } 193 } 194 if (flag & SOR_BACKWARD_SWEEP) { 195 for ( i=m-1; i>=0; i-- ) { 196 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 197 x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 198 } 199 } 200 } 201 return 0; 202 } 203 204 /* -----------------------------------------------------------------*/ 205 static int MatMultTrans_Dense(Mat matin,Vec xx,Vec yy) 206 { 207 Mat_Dense *mat = (Mat_Dense *) matin->data; 208 Scalar *v = mat->v, *x, *y; 209 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 210 VecGetArray(xx,&x), VecGetArray(yy,&y); 211 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 212 x, &_One, &_DZero, y, &_One ); 213 return 0; 214 } 215 static int MatMult_Dense(Mat matin,Vec xx,Vec yy) 216 { 217 Mat_Dense *mat = (Mat_Dense *) matin->data; 218 Scalar *v = mat->v, *x, *y; 219 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 220 VecGetArray(xx,&x); VecGetArray(yy,&y); 221 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 222 x, &_One, &_DZero, y, &_One ); 223 return 0; 224 } 225 static int MatMultAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 226 { 227 Mat_Dense *mat = (Mat_Dense *) matin->data; 228 Scalar *v = mat->v, *x, *y, *z; 229 int _One=1; Scalar _DOne=1.0; 230 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 231 if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 232 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 233 x, &_One, &_DOne, y, &_One ); 234 return 0; 235 } 236 static int MatMultTransAdd_Dense(Mat matin,Vec xx,Vec zz,Vec yy) 237 { 238 Mat_Dense *mat = (Mat_Dense *) matin->data; 239 Scalar *v = mat->v, *x, *y, *z; 240 int _One=1; 241 Scalar _DOne=1.0; 242 VecGetArray(xx,&x); VecGetArray(yy,&y); 243 VecGetArray(zz,&z); 244 if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 245 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 246 x, &_One, &_DOne, y, &_One ); 247 return 0; 248 } 249 250 /* -----------------------------------------------------------------*/ 251 static int MatGetRow_Dense(Mat matin,int row,int *ncols,int **cols, 252 Scalar **vals) 253 { 254 Mat_Dense *mat = (Mat_Dense *) matin->data; 255 Scalar *v; 256 int i; 257 *ncols = mat->n; 258 if (cols) { 259 *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols); 260 for ( i=0; i<mat->n; i++ ) *cols[i] = i; 261 } 262 if (vals) { 263 *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals); 264 v = mat->v + row; 265 for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 266 } 267 return 0; 268 } 269 static int MatRestoreRow_Dense(Mat matin,int row,int *ncols,int **cols, 270 Scalar **vals) 271 { 272 if (cols) { FREE(*cols); } 273 if (vals) { FREE(*vals); } 274 return 0; 275 } 276 /* ----------------------------------------------------------------*/ 277 static int MatInsert_Dense(Mat matin,int m,int *indexm,int n, 278 int *indexn,Scalar *v,InsertMode addv) 279 { 280 Mat_Dense *mat = (Mat_Dense *) matin->data; 281 int i,j; 282 283 if (!mat->roworiented) { 284 if (addv == InsertValues) { 285 for ( j=0; j<n; j++ ) { 286 if (indexn[j] < 0) {*v += m; continue;} 287 for ( i=0; i<m; i++ ) { 288 if (indexm[i] < 0) {*v++; continue;} 289 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 290 } 291 } 292 } 293 else { 294 for ( j=0; j<n; j++ ) { 295 if (indexn[j] < 0) {*v += m; continue;} 296 for ( i=0; i<m; i++ ) { 297 if (indexm[i] < 0) {*v++; continue;} 298 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 299 } 300 } 301 } 302 } 303 else { 304 if (addv == InsertValues) { 305 for ( i=0; i<m; i++ ) { 306 if (indexm[i] < 0) {*v += n; continue;} 307 for ( j=0; j<n; j++ ) { 308 if (indexn[j] < 0) {*v++; continue;} 309 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 310 } 311 } 312 } 313 else { 314 for ( i=0; i<m; i++ ) { 315 if (indexm[i] < 0) {*v += n; continue;} 316 for ( j=0; j<n; j++ ) { 317 if (indexn[j] < 0) {*v++; continue;} 318 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 319 } 320 } 321 } 322 } 323 return 0; 324 } 325 326 /* -----------------------------------------------------------------*/ 327 static int MatCopy_Dense(Mat matin,Mat *newmat) 328 { 329 Mat_Dense *mat = (Mat_Dense *) matin->data; 330 int ierr; 331 Mat newi; 332 Mat_Dense *l; 333 if ((ierr = MatCreateSequentialDense(mat->m,mat->n,&newi))) SETERR(ierr,0); 334 l = (Mat_Dense *) newi->data; 335 MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 336 *newmat = newi; 337 return 0; 338 } 339 #include "viewer.h" 340 341 int MatView_Dense(PetscObject obj,Viewer ptr) 342 { 343 Mat matin = (Mat) obj; 344 Mat_Dense *mat = (Mat_Dense *) matin->data; 345 Scalar *v; 346 int i,j; 347 PetscObject ojb = (PetscObject) ptr; 348 349 if (ojb && ojb->cookie == VIEWER_COOKIE && ojb->type == MATLAB_VIEWER) { 350 return ViewerMatlabPutArray(ptr,mat->m,mat->n,mat->v); 351 } 352 else { 353 FILE *fd = ViewerFileGetPointer(ptr); 354 for ( i=0; i<mat->m; i++ ) { 355 v = mat->v + i; 356 for ( j=0; j<mat->n; j++ ) { 357 #if defined(PETSC_COMPLEX) 358 fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 359 #else 360 fprintf(fd,"%6.4e ",*v); v += mat->m; 361 #endif 362 } 363 fprintf(fd,"\n"); 364 } 365 } 366 return 0; 367 } 368 369 370 static int MatDestroy_Dense(PetscObject obj) 371 { 372 Mat mat = (Mat) obj; 373 Mat_Dense *l = (Mat_Dense *) mat->data; 374 #if defined(PETSC_LOG) 375 PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 376 #endif 377 if (l->pivots) FREE(l->pivots); 378 FREE(l); 379 PLogObjectDestroy(mat); 380 PETSCHEADERDESTROY(mat); 381 return 0; 382 } 383 384 static int MatTrans_Dense(Mat matin) 385 { 386 Mat_Dense *mat = (Mat_Dense *) matin->data; 387 int k,j; 388 Scalar *v = mat->v, tmp; 389 if (mat->m != mat->n) { 390 SETERR(1,"Cannot transpose rectangular dense matrix"); 391 } 392 for ( j=0; j<mat->m; j++ ) { 393 for ( k=0; k<j; k++ ) { 394 tmp = v[j + k*mat->n]; 395 v[j + k*mat->n] = v[k + j*mat->n]; 396 v[k + j*mat->n] = tmp; 397 } 398 } 399 return 0; 400 } 401 402 static int MatEqual_Dense(Mat matin1,Mat matin2) 403 { 404 Mat_Dense *mat1 = (Mat_Dense *) matin1->data; 405 Mat_Dense *mat2 = (Mat_Dense *) matin2->data; 406 int i; 407 Scalar *v1 = mat1->v, *v2 = mat2->v; 408 if (mat1->m != mat2->m) return 0; 409 if (mat1->n != mat2->n) return 0; 410 for ( i=0; i<mat1->m*mat1->n; i++ ) { 411 if (*v1 != *v2) return 0; 412 v1++; v2++; 413 } 414 return 1; 415 } 416 417 static int MatGetDiag_Dense(Mat matin,Vec v) 418 { 419 Mat_Dense *mat = (Mat_Dense *) matin->data; 420 int i, n; 421 Scalar *x; 422 CHKTYPE(v,SEQVECTOR); 423 VecGetArray(v,&x); VecGetSize(v,&n); 424 if (n != mat->m) SETERR(1,"Nonconforming matrix and vector"); 425 for ( i=0; i<mat->m; i++ ) { 426 x[i] = mat->v[i*mat->m + i]; 427 } 428 return 0; 429 } 430 431 static int MatScale_Dense(Mat matin,Vec ll,Vec rr) 432 { 433 Mat_Dense *mat = (Mat_Dense *) matin->data; 434 Scalar *l,*r,x,*v; 435 int i,j,m = mat->m, n = mat->n; 436 if (ll) { 437 VecGetArray(ll,&l); VecGetSize(ll,&m); 438 if (m != mat->m) SETERR(1,"Left scaling vector wrong length"); 439 for ( i=0; i<m; i++ ) { 440 x = l[i]; 441 v = mat->v + i; 442 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 443 } 444 } 445 if (rr) { 446 VecGetArray(rr,&r); VecGetSize(rr,&n); 447 if (n != mat->n) SETERR(1,"Right scaling vector wrong length"); 448 for ( i=0; i<n; i++ ) { 449 x = r[i]; 450 v = mat->v + i*m; 451 for ( j=0; j<m; j++ ) { (*v++) *= x;} 452 } 453 } 454 return 0; 455 } 456 457 458 static int MatNorm_Dense(Mat matin,int type,double *norm) 459 { 460 Mat_Dense *mat = (Mat_Dense *) matin->data; 461 Scalar *v = mat->v; 462 double sum = 0.0; 463 int i, j; 464 if (type == NORM_FROBENIUS) { 465 for (i=0; i<mat->n*mat->m; i++ ) { 466 #if defined(PETSC_COMPLEX) 467 sum += real(conj(*v)*(*v)); v++; 468 #else 469 sum += (*v)*(*v); v++; 470 #endif 471 } 472 *norm = sqrt(sum); 473 } 474 else if (type == NORM_1) { 475 *norm = 0.0; 476 for ( j=0; j<mat->n; j++ ) { 477 sum = 0.0; 478 for ( i=0; i<mat->m; i++ ) { 479 #if defined(PETSC_COMPLEX) 480 sum += abs(*v++); 481 #else 482 sum += fabs(*v++); 483 #endif 484 } 485 if (sum > *norm) *norm = sum; 486 } 487 } 488 else if (type == NORM_INFINITY) { 489 *norm = 0.0; 490 for ( j=0; j<mat->m; j++ ) { 491 v = mat->v + j; 492 sum = 0.0; 493 for ( i=0; i<mat->n; i++ ) { 494 #if defined(PETSC_COMPLEX) 495 sum += abs(*v); v += mat->m; 496 #else 497 sum += fabs(*v); v += mat->m; 498 #endif 499 } 500 if (sum > *norm) *norm = sum; 501 } 502 } 503 else { 504 SETERR(1,"No support for the two norm yet"); 505 } 506 return 0; 507 } 508 509 static int MatInsOpt_Dense(Mat aijin,int op) 510 { 511 Mat_Dense *aij = (Mat_Dense *) aijin->data; 512 if (op == ROW_ORIENTED) aij->roworiented = 1; 513 else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 514 /* doesn't care about sorted rows or columns */ 515 return 0; 516 } 517 518 static int MatZero_Dense(Mat A) 519 { 520 Mat_Dense *l = (Mat_Dense *) A->data; 521 MEMSET(l->v,0,l->m*l->n*sizeof(Scalar)); 522 return 0; 523 } 524 525 static int MatZeroRows_Dense(Mat A,IS is,Scalar *diag) 526 { 527 Mat_Dense *l = (Mat_Dense *) A->data; 528 int n = l->n, i, j,ierr,N, *rows; 529 Scalar *slot; 530 ierr = ISGetLocalSize(is,&N); CHKERR(ierr); 531 ierr = ISGetIndices(is,&rows); CHKERR(ierr); 532 for ( i=0; i<N; i++ ) { 533 slot = l->v + rows[i]; 534 for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 535 } 536 if (diag) { 537 for ( i=0; i<N; i++ ) { 538 slot = l->v + (n+1)*rows[i]; 539 *slot = *diag; 540 } 541 } 542 ISRestoreIndices(is,&rows); 543 return 0; 544 } 545 546 static int MatSize_Dense(Mat matin,int *m,int *n) 547 { 548 Mat_Dense *mat = (Mat_Dense *) matin->data; 549 *m = mat->m; *n = mat->n; 550 return 0; 551 } 552 553 static int MatGetArray_Dense(Mat matin,Scalar **array) 554 { 555 Mat_Dense *mat = (Mat_Dense *) matin->data; 556 *array = mat->v; 557 return 0; 558 } 559 /* -------------------------------------------------------------------*/ 560 static struct _MatOps MatOps = {MatInsert_Dense, 561 MatGetRow_Dense, MatRestoreRow_Dense, 562 MatMult_Dense, MatMultAdd_Dense, 563 MatMultTrans_Dense, MatMultTransAdd_Dense, 564 MatSolve_Dense,MatSolveAdd_Dense, 565 MatSolveTrans_Dense,MatSolveTransAdd_Dense, 566 MatLUFactor_Dense,MatChFactor_Dense, 567 MatRelax_Dense, 568 MatTrans_Dense, 569 MatNz_Dense,MatMemory_Dense,MatEqual_Dense, 570 MatCopy_Dense, 571 MatGetDiag_Dense,MatScale_Dense,MatNorm_Dense, 572 0,0, 573 0, MatInsOpt_Dense,MatZero_Dense,MatZeroRows_Dense,0, 574 MatLUFactorSymbolic_Dense,MatLUFactorNumeric_Dense, 575 MatChFactorSymbolic_Dense,MatChFactorNumeric_Dense, 576 MatSize_Dense,MatSize_Dense,0, 577 0,0,MatGetArray_Dense 578 }; 579 /*@ 580 MatCreateSequentialDense - Creates a sequential dense matrix that 581 is stored in the usual Fortran 77 manner. Many of the matrix 582 operations use the BLAS and LAPACK routines. 583 584 Input Parameters: 585 . m, n - the number of rows and columns in the matrix. 586 587 Output Parameter: 588 . newmat - the matrix created. 589 590 Keywords: dense matrix, lapack, blas 591 @*/ 592 int MatCreateSequentialDense(int m,int n,Mat *newmat) 593 { 594 int size = sizeof(Mat_Dense) + m*n*sizeof(Scalar); 595 Mat mat; 596 Mat_Dense *l; 597 *newmat = 0; 598 PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATDENSE,MPI_COMM_SELF); 599 PLogObjectCreate(mat); 600 l = (Mat_Dense *) MALLOC(size); CHKPTR(l); 601 mat->ops = &MatOps; 602 mat->destroy = MatDestroy_Dense; 603 mat->view = MatView_Dense; 604 mat->data = (void *) l; 605 mat->factor = 0; 606 mat->col = 0; 607 mat->row = 0; 608 609 l->m = m; 610 l->n = n; 611 l->v = (Scalar *) (l + 1); 612 l->pivots = 0; 613 l->roworiented = 1; 614 615 MEMSET(l->v,0,m*n*sizeof(Scalar)); 616 *newmat = mat; 617 return 0; 618 } 619 620 int MatCreate_Dense(Mat matin,Mat *newmat) 621 { 622 Mat_Dense *m = (Mat_Dense *) matin->data; 623 return MatCreateSequentialDense(m->m,m->n,newmat); 624 } 625