1 #ifndef lint 2 static char vcid[] = "$Id: dense.c,v 1.66 1995/10/18 03:18:20 curfman Exp curfman $"; 3 #endif 4 5 #include "dense.h" 6 #include "pinclude/plapack.h" 7 #include "pinclude/pviewer.h" 8 9 int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y) 10 { 11 Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data; 12 int N = x->m*x->n, one = 1; 13 BLaxpy_( &N, alpha, x->v, &one, y->v, &one ); 14 PLogFlops(2*N-1); 15 return 0; 16 } 17 18 static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 19 { 20 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 21 int i,N = mat->m*mat->n,count = 0; 22 Scalar *v = mat->v; 23 for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 24 *nz = count; *nzalloc = N; *mem = (int)A->mem; 25 return 0; 26 } 27 28 /* ---------------------------------------------------------------*/ 29 /* COMMENT: I have chosen to hide column permutation in the pivots, 30 rather than put it in the Mat->col slot.*/ 31 static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f) 32 { 33 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 34 int info; 35 if (!mat->pivots) { 36 mat->pivots = (int *) PETSCMALLOC(mat->m*sizeof(int));CHKPTRQ(mat->pivots); 37 PLogObjectMemory(A,mat->m*sizeof(int)); 38 } 39 LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 40 if (info) SETERRQ(1,"MatLUFactor_SeqDense:Bad LU factorization"); 41 A->factor = FACTOR_LU; 42 return 0; 43 } 44 static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact) 45 { 46 int ierr; 47 ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr); 48 return 0; 49 } 50 static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 51 { 52 return MatLUFactor(*fact,0,0,1.0); 53 } 54 static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact) 55 { 56 int ierr; 57 ierr = MatConvert(A,MATSAME,fact); CHKERRQ(ierr); 58 return 0; 59 } 60 static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 61 { 62 return MatCholeskyFactor(*fact,0,1.0); 63 } 64 static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f) 65 { 66 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 67 int info; 68 if (mat->pivots) { 69 PETSCFREE(mat->pivots); 70 PLogObjectMemory(A,-mat->m*sizeof(int)); 71 mat->pivots = 0; 72 } 73 LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 74 if (info) SETERRQ(1,"MatCholeskyFactor_SeqDense:Bad factorization"); 75 A->factor = FACTOR_CHOLESKY; 76 return 0; 77 } 78 79 static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 80 { 81 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 82 int one = 1, info; 83 Scalar *x, *y; 84 VecGetArray(xx,&x); VecGetArray(yy,&y); 85 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 86 if (A->factor == FACTOR_LU) { 87 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 88 } 89 else if (A->factor == FACTOR_CHOLESKY){ 90 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 91 } 92 else SETERRQ(1,"MatSolve_SeqDense:Matrix must be factored to solve"); 93 if (info) SETERRQ(1,"MatSolve_SeqDense:Bad solve"); 94 return 0; 95 } 96 static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy) 97 { 98 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 99 int one = 1, info; 100 Scalar *x, *y; 101 VecGetArray(xx,&x); VecGetArray(yy,&y); 102 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 103 /* assume if pivots exist then use LU; else Cholesky */ 104 if (mat->pivots) { 105 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 106 } 107 else { 108 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 109 } 110 if (info) SETERRQ(1,"MatSolveTrans_SeqDense:Bad solve"); 111 return 0; 112 } 113 static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 114 { 115 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 116 int one = 1, info,ierr; 117 Scalar *x, *y, sone = 1.0; 118 Vec tmp = 0; 119 VecGetArray(xx,&x); VecGetArray(yy,&y); 120 if (yy == zz) { 121 ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 122 PLogObjectParent(A,tmp); 123 ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 124 } 125 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 126 /* assume if pivots exist then use LU; else Cholesky */ 127 if (mat->pivots) { 128 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 129 } 130 else { 131 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 132 } 133 if (info) SETERRQ(1,"MatSolveAdd_SeqDense:Bad solve"); 134 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 135 else VecAXPY(&sone,zz,yy); 136 return 0; 137 } 138 static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 139 { 140 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 141 int one = 1, info,ierr; 142 Scalar *x, *y, sone = 1.0; 143 Vec tmp; 144 VecGetArray(xx,&x); VecGetArray(yy,&y); 145 if (yy == zz) { 146 ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 147 PLogObjectParent(A,tmp); 148 ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 149 } 150 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 151 /* assume if pivots exist then use LU; else Cholesky */ 152 if (mat->pivots) { 153 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 154 } 155 else { 156 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 157 } 158 if (info) SETERRQ(1,"MatSolveTransAdd_SeqDense:Bad solve"); 159 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 160 else VecAXPY(&sone,zz,yy); 161 return 0; 162 } 163 /* ------------------------------------------------------------------*/ 164 static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 165 double shift,int its,Vec xx) 166 { 167 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 168 Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 169 int o = 1,ierr, m = mat->m, i; 170 171 if (flag & SOR_ZERO_INITIAL_GUESS) { 172 /* this is a hack fix, should have another version without 173 the second BLdot */ 174 ierr = VecSet(&zero,xx); CHKERRQ(ierr); 175 } 176 VecGetArray(xx,&x); VecGetArray(bb,&b); 177 while (its--) { 178 if (flag & SOR_FORWARD_SWEEP){ 179 for ( i=0; i<m; i++ ) { 180 #if defined(PETSC_COMPLEX) 181 /* cannot use BLAS dot for complex because compiler/linker is 182 not happy about returning a double complex */ 183 int _i; 184 Scalar sum = b[i]; 185 for ( _i=0; _i<m; _i++ ) { 186 sum -= conj(v[i+_i*m])*x[_i]; 187 } 188 xt = sum; 189 #else 190 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 191 #endif 192 x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 193 } 194 } 195 if (flag & SOR_BACKWARD_SWEEP) { 196 for ( i=m-1; i>=0; i-- ) { 197 #if defined(PETSC_COMPLEX) 198 /* cannot use BLAS dot for complex because compiler/linker is 199 not happy about returning a double complex */ 200 int _i; 201 Scalar sum = b[i]; 202 for ( _i=0; _i<m; _i++ ) { 203 sum -= conj(v[i+_i*m])*x[_i]; 204 } 205 xt = sum; 206 #else 207 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 208 #endif 209 x[i] = (1. - omega)*x[i] + omega*(xt/(v[i + i*m]+shift) + x[i]); 210 } 211 } 212 } 213 return 0; 214 } 215 216 /* -----------------------------------------------------------------*/ 217 static int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy) 218 { 219 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 220 Scalar *v = mat->v, *x, *y; 221 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 222 VecGetArray(xx,&x), VecGetArray(yy,&y); 223 LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 224 return 0; 225 } 226 static int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 227 { 228 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 229 Scalar *v = mat->v, *x, *y; 230 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 231 VecGetArray(xx,&x); VecGetArray(yy,&y); 232 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 233 return 0; 234 } 235 static int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 236 { 237 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 238 Scalar *v = mat->v, *x, *y, *z; 239 int _One=1; Scalar _DOne=1.0; 240 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 241 if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 242 LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 243 return 0; 244 } 245 static int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 246 { 247 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 248 Scalar *v = mat->v, *x, *y, *z; 249 int _One=1; 250 Scalar _DOne=1.0; 251 VecGetArray(xx,&x); VecGetArray(yy,&y); 252 VecGetArray(zz,&z); 253 if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 254 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 255 return 0; 256 } 257 258 /* -----------------------------------------------------------------*/ 259 static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 260 { 261 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 262 Scalar *v; 263 int i; 264 *ncols = mat->n; 265 if (cols) { 266 *cols = (int *) PETSCMALLOC(mat->n*sizeof(int)); CHKPTRQ(*cols); 267 for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 268 } 269 if (vals) { 270 *vals = (Scalar *) PETSCMALLOC(mat->n*sizeof(Scalar)); CHKPTRQ(*vals); 271 v = mat->v + row; 272 for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 273 } 274 return 0; 275 } 276 static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 277 { 278 if (cols) { PETSCFREE(*cols); } 279 if (vals) { PETSCFREE(*vals); } 280 return 0; 281 } 282 /* ----------------------------------------------------------------*/ 283 static int MatInsert_SeqDense(Mat A,int m,int *indexm,int n, 284 int *indexn,Scalar *v,InsertMode addv) 285 { 286 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 287 int i,j; 288 289 if (!mat->roworiented) { 290 if (addv == INSERT_VALUES) { 291 for ( j=0; j<n; j++ ) { 292 for ( i=0; i<m; i++ ) { 293 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 294 } 295 } 296 } 297 else { 298 for ( j=0; j<n; j++ ) { 299 for ( i=0; i<m; i++ ) { 300 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 301 } 302 } 303 } 304 } 305 else { 306 if (addv == INSERT_VALUES) { 307 for ( i=0; i<m; i++ ) { 308 for ( j=0; j<n; j++ ) { 309 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 310 } 311 } 312 } 313 else { 314 for ( i=0; i<m; i++ ) { 315 for ( j=0; j<n; j++ ) { 316 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 317 } 318 } 319 } 320 } 321 return 0; 322 } 323 324 /* -----------------------------------------------------------------*/ 325 static int MatCopyPrivate_SeqDense(Mat A,Mat *newmat) 326 { 327 Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l; 328 int ierr; 329 Mat newi; 330 331 ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,&newi); CHKERRQ(ierr); 332 l = (Mat_SeqDense *) newi->data; 333 PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 334 *newmat = newi; 335 return 0; 336 } 337 338 #include "pinclude/pviewer.h" 339 #include "sysio.h" 340 341 static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 342 { 343 Mat_SeqDense *a = (Mat_SeqDense *) A->data; 344 int ierr, i, j, format; 345 FILE *fd; 346 char *outputname; 347 Scalar *v; 348 349 ierr = ViewerFileGetPointer_Private(viewer,&fd); CHKERRQ(ierr); 350 ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 351 ierr = ViewerFileGetFormat_Private(viewer,&format); 352 if (format == FILE_FORMAT_INFO) { 353 ; /* do nothing for now */ 354 } 355 else { 356 for ( i=0; i<a->m; i++ ) { 357 v = a->v + i; 358 for ( j=0; j<a->n; j++ ) { 359 #if defined(PETSC_COMPLEX) 360 fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m; 361 #else 362 fprintf(fd,"%6.4e ",*v); v += a->m; 363 #endif 364 } 365 fprintf(fd,"\n"); 366 } 367 } 368 fflush(fd); 369 return 0; 370 } 371 372 static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 373 { 374 Mat_SeqDense *a = (Mat_SeqDense *) A->data; 375 int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 376 Scalar *v, *anonz; 377 378 ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 379 col_lens = (int *) PETSCMALLOC( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens); 380 col_lens[0] = MAT_COOKIE; 381 col_lens[1] = m; 382 col_lens[2] = n; 383 col_lens[3] = nz; 384 385 /* store lengths of each row and write (including header) to file */ 386 for ( i=0; i<m; i++ ) col_lens[4+i] = n; 387 ierr = SYWrite(fd,col_lens,4+m,SYINT,1); CHKERRQ(ierr); 388 389 /* Possibly should write in smaller increments, not whole matrix at once? */ 390 /* store column indices (zero start index) */ 391 ict = 0; 392 for ( i=0; i<m; i++ ) { 393 for ( j=0; j<n; j++ ) col_lens[ict++] = j; 394 } 395 ierr = SYWrite(fd,col_lens,nz,SYINT,0); CHKERRQ(ierr); 396 PETSCFREE(col_lens); 397 398 /* store nonzero values */ 399 anonz = (Scalar *) PETSCMALLOC((nz)*sizeof(Scalar)); CHKPTRQ(anonz); 400 ict = 0; 401 for ( i=0; i<m; i++ ) { 402 v = a->v + i; 403 for ( j=0; j<n; j++ ) { 404 anonz[ict++] = *v; v += a->m; 405 } 406 } 407 ierr = SYWrite(fd,anonz,nz,SYSCALAR,0); CHKERRQ(ierr); 408 PETSCFREE(anonz); 409 return 0; 410 } 411 412 static int MatView_SeqDense(PetscObject obj,Viewer viewer) 413 { 414 Mat A = (Mat) obj; 415 Mat_SeqDense *a = (Mat_SeqDense*) A->data; 416 PetscObject vobj = (PetscObject) viewer; 417 418 if (!viewer) { 419 viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 420 } 421 if (vobj->cookie == VIEWER_COOKIE) { 422 if (vobj->type == MATLAB_VIEWER) { 423 return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); 424 } 425 else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER) { 426 return MatView_SeqDense_ASCII(A,viewer); 427 } 428 else if (vobj->type == BINARY_FILE_VIEWER) { 429 return MatView_SeqDense_Binary(A,viewer); 430 } 431 } 432 return 0; 433 } 434 435 static int MatDestroy_SeqDense(PetscObject obj) 436 { 437 Mat mat = (Mat) obj; 438 Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 439 #if defined(PETSC_LOG) 440 PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 441 #endif 442 if (l->pivots) PETSCFREE(l->pivots); 443 PETSCFREE(l); 444 PLogObjectDestroy(mat); 445 PETSCHEADERDESTROY(mat); 446 return 0; 447 } 448 449 static int MatTranspose_SeqDense(Mat A,Mat *matout) 450 { 451 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 452 int k, j, m, n; 453 Scalar *v, tmp; 454 455 v = mat->v; m = mat->m; n = mat->n; 456 if (!matout) { /* in place transpose */ 457 if (m != n) SETERRQ(1,"MatTranspose_SeqDense:Not for rectangular matrix in place"); 458 for ( j=0; j<m; j++ ) { 459 for ( k=0; k<j; k++ ) { 460 tmp = v[j + k*n]; 461 v[j + k*n] = v[k + j*n]; 462 v[k + j*n] = tmp; 463 } 464 } 465 } 466 else { /* out-of-place transpose */ 467 int ierr; 468 Mat tmat; 469 Mat_SeqDense *tmatd; 470 Scalar *v2; 471 ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,&tmat); CHKERRQ(ierr); 472 tmatd = (Mat_SeqDense *) tmat->data; 473 v = mat->v; v2 = tmatd->v; 474 for ( j=0; j<n; j++ ) { 475 for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 476 } 477 ierr = MatAssemblyBegin(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 478 ierr = MatAssemblyEnd(tmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 479 *matout = tmat; 480 } 481 return 0; 482 } 483 484 static int MatEqual_SeqDense(Mat A1,Mat A2) 485 { 486 Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 487 Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 488 int i; 489 Scalar *v1 = mat1->v, *v2 = mat2->v; 490 if (mat1->m != mat2->m) return 0; 491 if (mat1->n != mat2->n) return 0; 492 for ( i=0; i<mat1->m*mat1->n; i++ ) { 493 if (*v1 != *v2) return 0; 494 v1++; v2++; 495 } 496 return 1; 497 } 498 499 static int MatGetDiagonal_SeqDense(Mat A,Vec v) 500 { 501 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 502 int i, n; 503 Scalar *x; 504 VecGetArray(v,&x); VecGetSize(v,&n); 505 if (n != mat->m) SETERRQ(1,"MatGetDiagonal_SeqDense:Nonconforming mat and vec"); 506 for ( i=0; i<mat->m; i++ ) { 507 x[i] = mat->v[i*mat->m + i]; 508 } 509 return 0; 510 } 511 512 static int MatScale_SeqDense(Mat A,Vec ll,Vec rr) 513 { 514 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 515 Scalar *l,*r,x,*v; 516 int i,j,m = mat->m, n = mat->n; 517 if (ll) { 518 VecGetArray(ll,&l); VecGetSize(ll,&m); 519 if (m != mat->m) SETERRQ(1,"MatScale_SeqDense:Left scaling vec wrong size"); 520 for ( i=0; i<m; i++ ) { 521 x = l[i]; 522 v = mat->v + i; 523 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 524 } 525 } 526 if (rr) { 527 VecGetArray(rr,&r); VecGetSize(rr,&n); 528 if (n != mat->n) SETERRQ(1,"MatScale_SeqDense:Right scaling vec wrong size"); 529 for ( i=0; i<n; i++ ) { 530 x = r[i]; 531 v = mat->v + i*m; 532 for ( j=0; j<m; j++ ) { (*v++) *= x;} 533 } 534 } 535 return 0; 536 } 537 538 static int MatNorm_SeqDense(Mat A,MatNormType type,double *norm) 539 { 540 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 541 Scalar *v = mat->v; 542 double sum = 0.0; 543 int i, j; 544 if (type == NORM_FROBENIUS) { 545 for (i=0; i<mat->n*mat->m; i++ ) { 546 #if defined(PETSC_COMPLEX) 547 sum += real(conj(*v)*(*v)); v++; 548 #else 549 sum += (*v)*(*v); v++; 550 #endif 551 } 552 *norm = sqrt(sum); 553 } 554 else if (type == NORM_1) { 555 *norm = 0.0; 556 for ( j=0; j<mat->n; j++ ) { 557 sum = 0.0; 558 for ( i=0; i<mat->m; i++ ) { 559 #if defined(PETSC_COMPLEX) 560 sum += abs(*v++); 561 #else 562 sum += fabs(*v++); 563 #endif 564 } 565 if (sum > *norm) *norm = sum; 566 } 567 } 568 else if (type == NORM_INFINITY) { 569 *norm = 0.0; 570 for ( j=0; j<mat->m; j++ ) { 571 v = mat->v + j; 572 sum = 0.0; 573 for ( i=0; i<mat->n; i++ ) { 574 #if defined(PETSC_COMPLEX) 575 sum += abs(*v); v += mat->m; 576 #else 577 sum += fabs(*v); v += mat->m; 578 #endif 579 } 580 if (sum > *norm) *norm = sum; 581 } 582 } 583 else { 584 SETERRQ(1,"MatNorm_SeqDense:No two norm"); 585 } 586 return 0; 587 } 588 589 static int MatSetOption_SeqDense(Mat A,MatOption op) 590 { 591 Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 592 if (op == ROW_ORIENTED) aij->roworiented = 1; 593 else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 594 else if (op == ROWS_SORTED || 595 op == SYMMETRIC_MATRIX || 596 op == STRUCTURALLY_SYMMETRIC_MATRIX || 597 op == NO_NEW_NONZERO_LOCATIONS || 598 op == YES_NEW_NONZERO_LOCATIONS || 599 op == NO_NEW_DIAGONALS || 600 op == YES_NEW_DIAGONALS) 601 PLogInfo((PetscObject)A,"Info:MatSetOption_SeqDense:Option ignored\n"); 602 else 603 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqDense:unknown option");} 604 return 0; 605 } 606 607 static int MatZeroEntries_SeqDense(Mat A) 608 { 609 Mat_SeqDense *l = (Mat_SeqDense *) A->data; 610 PetscZero(l->v,l->m*l->n*sizeof(Scalar)); 611 return 0; 612 } 613 614 static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 615 { 616 Mat_SeqDense *l = (Mat_SeqDense *) A->data; 617 int n = l->n, i, j,ierr,N, *rows; 618 Scalar *slot; 619 ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 620 ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 621 for ( i=0; i<N; i++ ) { 622 slot = l->v + rows[i]; 623 for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 624 } 625 if (diag) { 626 for ( i=0; i<N; i++ ) { 627 slot = l->v + (n+1)*rows[i]; 628 *slot = *diag; 629 } 630 } 631 ISRestoreIndices(is,&rows); 632 return 0; 633 } 634 635 static int MatGetSize_SeqDense(Mat A,int *m,int *n) 636 { 637 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 638 *m = mat->m; *n = mat->n; 639 return 0; 640 } 641 642 static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 643 { 644 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 645 *m = 0; *n = mat->m; 646 return 0; 647 } 648 649 static int MatGetArray_SeqDense(Mat A,Scalar **array) 650 { 651 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 652 *array = mat->v; 653 return 0; 654 } 655 656 657 static int MatGetSubMatrixInPlace_SeqDense(Mat A,IS isrow,IS iscol) 658 { 659 SETERRQ(1,"MatGetSubMatrixInPlace_SeqDense:not done"); 660 } 661 662 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,Mat *submat) 663 { 664 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 665 int nznew, *smap, i, j, ierr, oldcols = mat->n; 666 int *irow, *icol, nrows, ncols, *cwork; 667 Scalar *vwork, *val; 668 Mat newmat; 669 670 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 671 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 672 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 673 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 674 675 smap = (int *) PETSCMALLOC(oldcols*sizeof(int)); CHKPTRQ(smap); 676 cwork = (int *) PETSCMALLOC(ncols*sizeof(int)); CHKPTRQ(cwork); 677 vwork = (Scalar *) PETSCMALLOC(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 678 PetscZero((char*)smap,oldcols*sizeof(int)); 679 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 680 681 /* Create and fill new matrix */ 682 ierr = MatCreateSeqDense(A->comm,nrows,ncols,&newmat);CHKERRQ(ierr); 683 for (i=0; i<nrows; i++) { 684 nznew = 0; 685 val = mat->v + irow[i]; 686 for (j=0; j<oldcols; j++) { 687 if (smap[j]) { 688 cwork[nznew] = smap[j] - 1; 689 vwork[nznew++] = val[j * mat->m]; 690 } 691 } 692 ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES); 693 CHKERRQ(ierr); 694 } 695 ierr = MatAssemblyBegin(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 696 ierr = MatAssemblyEnd(newmat,FINAL_ASSEMBLY); CHKERRQ(ierr); 697 698 /* Free work space */ 699 PETSCFREE(smap); PETSCFREE(cwork); PETSCFREE(vwork); 700 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 701 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 702 *submat = newmat; 703 return 0; 704 } 705 706 /* -------------------------------------------------------------------*/ 707 static struct _MatOps MatOps = {MatInsert_SeqDense, 708 MatGetRow_SeqDense, MatRestoreRow_SeqDense, 709 MatMult_SeqDense, MatMultAdd_SeqDense, 710 MatMultTrans_SeqDense, MatMultTransAdd_SeqDense, 711 MatSolve_SeqDense,MatSolveAdd_SeqDense, 712 MatSolveTrans_SeqDense,MatSolveTransAdd_SeqDense, 713 MatLUFactor_SeqDense,MatCholeskyFactor_SeqDense, 714 MatRelax_SeqDense, 715 MatTranspose_SeqDense, 716 MatGetInfo_SeqDense,MatEqual_SeqDense, 717 MatGetDiagonal_SeqDense,MatScale_SeqDense,MatNorm_SeqDense, 718 0,0, 719 0, MatSetOption_SeqDense,MatZeroEntries_SeqDense,MatZeroRows_SeqDense,0, 720 MatLUFactorSymbolic_SeqDense,MatLUFactorNumeric_SeqDense, 721 MatCholeskyFactorSymbolic_SeqDense,MatCholeskyFactorNumeric_SeqDense, 722 MatGetSize_SeqDense,MatGetSize_SeqDense,MatGetOwnershipRange_SeqDense, 723 0,0,MatGetArray_SeqDense,0,0, 724 MatGetSubMatrix_SeqDense,MatGetSubMatrixInPlace_SeqDense, 725 MatCopyPrivate_SeqDense,0,0,0,0, 726 MatAXPY_SeqDense}; 727 728 /*@C 729 MatCreateSeqDense - Creates a sequential dense matrix that 730 is stored in column major order (the usual Fortran 77 manner). Many 731 of the matrix operations use the BLAS and LAPACK routines. 732 733 Input Parameters: 734 . comm - MPI communicator, set to MPI_COMM_SELF 735 . m - number of rows 736 . n - number of column 737 738 Output Parameter: 739 . newmat - the matrix 740 741 .keywords: dense, matrix, LAPACK, BLAS 742 743 .seealso: MatCreate(), MatSetValues() 744 @*/ 745 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Mat *newmat) 746 { 747 int size = sizeof(Mat_SeqDense) + m*n*sizeof(Scalar); 748 Mat mat; 749 Mat_SeqDense *l; 750 *newmat = 0; 751 PETSCHEADERCREATE(mat,_Mat,MAT_COOKIE,MATSEQDENSE,comm); 752 PLogObjectCreate(mat); 753 l = (Mat_SeqDense *) PETSCMALLOC(size); CHKPTRQ(l); 754 PetscMemcpy(&mat->ops,&MatOps,sizeof(struct _MatOps)); 755 mat->destroy = MatDestroy_SeqDense; 756 mat->view = MatView_SeqDense; 757 mat->data = (void *) l; 758 mat->factor = 0; 759 PLogObjectMemory(mat,sizeof(struct _Mat) + size); 760 761 l->m = m; 762 l->n = n; 763 l->v = (Scalar *) (l + 1); 764 l->pivots = 0; 765 l->roworiented = 1; 766 767 PetscZero(l->v,m*n*sizeof(Scalar)); 768 *newmat = mat; 769 return 0; 770 } 771 772 int MatCreate_SeqDense(Mat A,Mat *newmat) 773 { 774 Mat_SeqDense *m = (Mat_SeqDense *) A->data; 775 return MatCreateSeqDense(A->comm,m->m,m->n,newmat); 776 } 777