1 #ifndef lint 2 static char vcid[] = "$Id: dense.c,v 1.118 1996/12/18 17:33:36 balay Exp balay $"; 3 #endif 4 /* 5 Defines the basic matrix operations for sequential dense. 6 */ 7 8 #include "src/mat/impls/dense/seq/dense.h" 9 #include "pinclude/plapack.h" 10 #include "pinclude/pviewer.h" 11 12 #undef __FUNCTION__ 13 #define __FUNCTION__ "MatAXPY_SeqDense" 14 int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y) 15 { 16 Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data; 17 int N = x->m*x->n, one = 1; 18 BLaxpy_( &N, alpha, x->v, &one, y->v, &one ); 19 PLogFlops(2*N-1); 20 return 0; 21 } 22 23 #undef __FUNCTION__ 24 #define __FUNCTION__ "MatGetInfo_SeqDense" 25 static int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 26 { 27 Mat_SeqDense *a = (Mat_SeqDense *) A->data; 28 int i,N = a->m*a->n,count = 0; 29 Scalar *v = a->v; 30 for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 31 32 info->rows_global = (double)a->m; 33 info->columns_global = (double)a->n; 34 info->rows_local = (double)a->m; 35 info->columns_local = (double)a->n; 36 info->block_size = 1.0; 37 info->nz_allocated = (double)N; 38 info->nz_used = (double)count; 39 info->nz_unneeded = (double)(N-count); 40 info->assemblies = (double)A->num_ass; 41 info->mallocs = 0; 42 info->memory = A->mem; 43 info->fill_ratio_given = 0; 44 info->fill_ratio_needed = 0; 45 info->factor_mallocs = 0; 46 47 return 0; 48 } 49 50 #undef __FUNCTION__ 51 #define __FUNCTION__ "MatScale_SeqDense" 52 static int MatScale_SeqDense(Scalar *alpha,Mat inA) 53 { 54 Mat_SeqDense *a = (Mat_SeqDense *) inA->data; 55 int one = 1, nz; 56 57 nz = a->m*a->n; 58 BLscal_( &nz, alpha, a->v, &one ); 59 PLogFlops(nz); 60 return 0; 61 } 62 63 /* ---------------------------------------------------------------*/ 64 /* COMMENT: I have chosen to hide column permutation in the pivots, 65 rather than put it in the Mat->col slot.*/ 66 #undef __FUNCTION__ 67 #define __FUNCTION__ "MatLUFactor_SeqDense" 68 static int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f) 69 { 70 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 71 int info; 72 73 if (!mat->pivots) { 74 mat->pivots = (int *) PetscMalloc((mat->m+1)*sizeof(int));CHKPTRQ(mat->pivots); 75 PLogObjectMemory(A,mat->m*sizeof(int)); 76 } 77 LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 78 if (info<0) SETERRQ(1,"Bad argument to LU factorization"); 79 if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 80 A->factor = FACTOR_LU; 81 PLogFlops((2*mat->n*mat->n*mat->n)/3); 82 return 0; 83 } 84 85 #undef __FUNCTION__ 86 #define __FUNCTION__ "MatConvertSameType_SeqDense" 87 static int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues) 88 { 89 Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l; 90 int ierr; 91 Mat newi; 92 93 ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi); CHKERRQ(ierr); 94 l = (Mat_SeqDense *) newi->data; 95 if (cpvalues == COPY_VALUES) { 96 PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 97 } 98 newi->assembled = PETSC_TRUE; 99 *newmat = newi; 100 return 0; 101 } 102 103 #undef __FUNCTION__ 104 #define __FUNCTION__ "MatLUFactorSymbolic_SeqDense" 105 static int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact) 106 { 107 return MatConvertSameType_SeqDense(A,fact,PETSC_FALSE); 108 } 109 110 #undef __FUNCTION__ 111 #define __FUNCTION__ "MatLUFactorNumeric_SeqDense" 112 static int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 113 { 114 Mat_SeqDense *mat = (Mat_SeqDense*) A->data, *l = (Mat_SeqDense*) (*fact)->data; 115 /* copy the numerical values */ 116 PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 117 (*fact)->factor = 0; 118 return MatLUFactor(*fact,0,0,1.0); 119 } 120 121 #undef __FUNCTION__ 122 #define __FUNCTION__ "MatCholeskyFactorSymbolic_SeqDense" 123 static int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact) 124 { 125 return MatConvert(A,MATSAME,fact); 126 } 127 128 #undef __FUNCTION__ 129 #define __FUNCTION__ "MatCholeskyFactorNumeric_SeqDense" 130 static int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 131 { 132 return MatCholeskyFactor(*fact,0,1.0); 133 } 134 135 #undef __FUNCTION__ 136 #define __FUNCTION__ "MatCholeskyFactor_SeqDense" 137 static int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f) 138 { 139 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 140 int info; 141 142 if (mat->pivots) { 143 PetscFree(mat->pivots); 144 PLogObjectMemory(A,-mat->m*sizeof(int)); 145 mat->pivots = 0; 146 } 147 LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 148 if (info) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization"); 149 A->factor = FACTOR_CHOLESKY; 150 PLogFlops((mat->n*mat->n*mat->n)/3); 151 return 0; 152 } 153 154 #undef __FUNCTION__ 155 #define __FUNCTION__ "MatSolve_SeqDense" 156 static int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 157 { 158 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 159 int one = 1, info, ierr; 160 Scalar *x, *y; 161 162 ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 163 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 164 if (A->factor == FACTOR_LU) { 165 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 166 } 167 else if (A->factor == FACTOR_CHOLESKY){ 168 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 169 } 170 else SETERRQ(1,"Matrix must be factored to solve"); 171 if (info) SETERRQ(1,"MBad solve"); 172 PLogFlops(mat->n*mat->n - mat->n); 173 return 0; 174 } 175 176 #undef __FUNCTION__ 177 #define __FUNCTION__ "MatSolveTrans_SeqDense" 178 static int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy) 179 { 180 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 181 int one = 1, info; 182 Scalar *x, *y; 183 184 VecGetArray(xx,&x); VecGetArray(yy,&y); 185 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 186 /* assume if pivots exist then use LU; else Cholesky */ 187 if (mat->pivots) { 188 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 189 } 190 else { 191 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 192 } 193 if (info) SETERRQ(1,"Bad solve"); 194 PLogFlops(mat->n*mat->n - mat->n); 195 return 0; 196 } 197 198 #undef __FUNCTION__ 199 #define __FUNCTION__ "MatSolveAdd_SeqDense" 200 static int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 201 { 202 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 203 int one = 1, info,ierr; 204 Scalar *x, *y, sone = 1.0; 205 Vec tmp = 0; 206 207 VecGetArray(xx,&x); VecGetArray(yy,&y); 208 if (yy == zz) { 209 ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 210 PLogObjectParent(A,tmp); 211 ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 212 } 213 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 214 /* assume if pivots exist then use LU; else Cholesky */ 215 if (mat->pivots) { 216 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 217 } 218 else { 219 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 220 } 221 if (info) SETERRQ(1,"Bad solve"); 222 if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 223 else VecAXPY(&sone,zz,yy); 224 PLogFlops(mat->n*mat->n - mat->n); 225 return 0; 226 } 227 228 #undef __FUNCTION__ 229 #define __FUNCTION__ "MatSolveTransAdd_SeqDense" 230 static int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 231 { 232 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 233 int one = 1, info,ierr; 234 Scalar *x, *y, sone = 1.0; 235 Vec tmp; 236 237 VecGetArray(xx,&x); VecGetArray(yy,&y); 238 if (yy == zz) { 239 ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 240 PLogObjectParent(A,tmp); 241 ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 242 } 243 PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 244 /* assume if pivots exist then use LU; else Cholesky */ 245 if (mat->pivots) { 246 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 247 } 248 else { 249 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 250 } 251 if (info) SETERRQ(1,"Bad solve"); 252 if (tmp) { 253 ierr = VecAXPY(&sone,tmp,yy); CHKERRQ(ierr); 254 ierr = VecDestroy(tmp); CHKERRQ(ierr); 255 } 256 else { 257 ierr = VecAXPY(&sone,zz,yy); CHKERRQ(ierr); 258 } 259 PLogFlops(mat->n*mat->n - mat->n); 260 return 0; 261 } 262 /* ------------------------------------------------------------------*/ 263 #undef __FUNCTION__ 264 #define __FUNCTION__ "MatRelax_SeqDense" 265 static int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 266 double shift,int its,Vec xx) 267 { 268 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 269 Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 270 int o = 1,ierr, m = mat->m, i; 271 272 if (flag & SOR_ZERO_INITIAL_GUESS) { 273 /* this is a hack fix, should have another version without 274 the second BLdot */ 275 ierr = VecSet(&zero,xx); CHKERRQ(ierr); 276 } 277 VecGetArray(xx,&x); VecGetArray(bb,&b); 278 while (its--) { 279 if (flag & SOR_FORWARD_SWEEP){ 280 for ( i=0; i<m; i++ ) { 281 #if defined(PETSC_COMPLEX) 282 /* cannot use BLAS dot for complex because compiler/linker is 283 not happy about returning a double complex */ 284 int _i; 285 Scalar sum = b[i]; 286 for ( _i=0; _i<m; _i++ ) { 287 sum -= conj(v[i+_i*m])*x[_i]; 288 } 289 xt = sum; 290 #else 291 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 292 #endif 293 x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 294 } 295 } 296 if (flag & SOR_BACKWARD_SWEEP) { 297 for ( i=m-1; i>=0; i-- ) { 298 #if defined(PETSC_COMPLEX) 299 /* cannot use BLAS dot for complex because compiler/linker is 300 not happy about returning a double complex */ 301 int _i; 302 Scalar sum = b[i]; 303 for ( _i=0; _i<m; _i++ ) { 304 sum -= conj(v[i+_i*m])*x[_i]; 305 } 306 xt = sum; 307 #else 308 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 309 #endif 310 x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 311 } 312 } 313 } 314 return 0; 315 } 316 317 /* -----------------------------------------------------------------*/ 318 #undef __FUNCTION__ 319 #define __FUNCTION__ "MatMultTrans_SeqDense" 320 int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy) 321 { 322 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 323 Scalar *v = mat->v, *x, *y; 324 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 325 VecGetArray(xx,&x), VecGetArray(yy,&y); 326 LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 327 PLogFlops(2*mat->m*mat->n - mat->n); 328 return 0; 329 } 330 331 #undef __FUNCTION__ 332 #define __FUNCTION__ "MatMult_SeqDense" 333 int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 334 { 335 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 336 Scalar *v = mat->v, *x, *y; 337 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 338 VecGetArray(xx,&x); VecGetArray(yy,&y); 339 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 340 PLogFlops(2*mat->m*mat->n - mat->m); 341 return 0; 342 } 343 344 #undef __FUNCTION__ 345 #define __FUNCTION__ "MatMultAdd_SeqDense" 346 int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 347 { 348 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 349 Scalar *v = mat->v, *x, *y, *z; 350 int _One=1; Scalar _DOne=1.0; 351 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 352 if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 353 LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 354 PLogFlops(2*mat->m*mat->n); 355 return 0; 356 } 357 358 #undef __FUNCTION__ 359 #define __FUNCTION__ "MatMultTransAdd_SeqDense" 360 int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 361 { 362 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 363 Scalar *v = mat->v, *x, *y, *z; 364 int _One=1; 365 Scalar _DOne=1.0; 366 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 367 if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar)); 368 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 369 PLogFlops(2*mat->m*mat->n); 370 return 0; 371 } 372 373 /* -----------------------------------------------------------------*/ 374 #undef __FUNCTION__ 375 #define __FUNCTION__ "MatGetRow_SeqDense" 376 static int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 377 { 378 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 379 Scalar *v; 380 int i; 381 382 *ncols = mat->n; 383 if (cols) { 384 *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int)); CHKPTRQ(*cols); 385 for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 386 } 387 if (vals) { 388 *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar)); CHKPTRQ(*vals); 389 v = mat->v + row; 390 for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 391 } 392 return 0; 393 } 394 395 #undef __FUNCTION__ 396 #define __FUNCTION__ "MatRestoreRow_SeqDense" 397 static int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 398 { 399 if (cols) { PetscFree(*cols); } 400 if (vals) { PetscFree(*vals); } 401 return 0; 402 } 403 /* ----------------------------------------------------------------*/ 404 #undef __FUNCTION__ 405 #define __FUNCTION__ "MatSetValues_SeqDense" 406 static int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 407 int *indexn,Scalar *v,InsertMode addv) 408 { 409 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 410 int i,j; 411 412 if (!mat->roworiented) { 413 if (addv == INSERT_VALUES) { 414 for ( j=0; j<n; j++ ) { 415 for ( i=0; i<m; i++ ) { 416 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 417 } 418 } 419 } 420 else { 421 for ( j=0; j<n; j++ ) { 422 for ( i=0; i<m; i++ ) { 423 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 424 } 425 } 426 } 427 } 428 else { 429 if (addv == INSERT_VALUES) { 430 for ( i=0; i<m; i++ ) { 431 for ( j=0; j<n; j++ ) { 432 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 433 } 434 } 435 } 436 else { 437 for ( i=0; i<m; i++ ) { 438 for ( j=0; j<n; j++ ) { 439 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 440 } 441 } 442 } 443 } 444 return 0; 445 } 446 447 #undef __FUNCTION__ 448 #define __FUNCTION__ "MatGetValues_SeqDense" 449 static int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v) 450 { 451 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 452 int i, j; 453 Scalar *vpt = v; 454 455 /* row-oriented output */ 456 for ( i=0; i<m; i++ ) { 457 for ( j=0; j<n; j++ ) { 458 *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]]; 459 } 460 } 461 return 0; 462 } 463 464 /* -----------------------------------------------------------------*/ 465 466 #include "sys.h" 467 468 #undef __FUNCTION__ 469 #define __FUNCTION__ "MatLoad_SeqDense" 470 int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A) 471 { 472 Mat_SeqDense *a; 473 Mat B; 474 int *scols, i, j, nz, ierr, fd, header[4], size; 475 int *rowlengths = 0, M, N, *cols; 476 Scalar *vals, *svals, *v; 477 MPI_Comm comm = ((PetscObject)viewer)->comm; 478 479 MPI_Comm_size(comm,&size); 480 if (size > 1) SETERRQ(1,"view must have one processor"); 481 ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 482 ierr = PetscBinaryRead(fd,header,4,BINARY_INT); CHKERRQ(ierr); 483 if (header[0] != MAT_COOKIE) SETERRQ(1,"Not matrix object"); 484 M = header[1]; N = header[2]; nz = header[3]; 485 486 if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 487 ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 488 B = *A; 489 a = (Mat_SeqDense *) B->data; 490 491 /* read in nonzero values */ 492 ierr = PetscBinaryRead(fd,a->v,M*N,BINARY_SCALAR); CHKERRQ(ierr); 493 494 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 495 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 496 ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr); 497 } else { 498 /* read row lengths */ 499 rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) ); CHKPTRQ(rowlengths); 500 ierr = PetscBinaryRead(fd,rowlengths,M,BINARY_INT); CHKERRQ(ierr); 501 502 /* create our matrix */ 503 ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 504 B = *A; 505 a = (Mat_SeqDense *) B->data; 506 v = a->v; 507 508 /* read column indices and nonzeros */ 509 cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cols); 510 ierr = PetscBinaryRead(fd,cols,nz,BINARY_INT); CHKERRQ(ierr); 511 vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) ); CHKPTRQ(vals); 512 ierr = PetscBinaryRead(fd,vals,nz,BINARY_SCALAR); CHKERRQ(ierr); 513 514 /* insert into matrix */ 515 for ( i=0; i<M; i++ ) { 516 for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j]; 517 svals += rowlengths[i]; scols += rowlengths[i]; 518 } 519 PetscFree(vals); PetscFree(cols); PetscFree(rowlengths); 520 521 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 522 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 523 } 524 return 0; 525 } 526 527 #include "pinclude/pviewer.h" 528 #include "sys.h" 529 530 #undef __FUNCTION__ 531 #define __FUNCTION__ "MatView_SeqDense_ASCII" 532 static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 533 { 534 Mat_SeqDense *a = (Mat_SeqDense *) A->data; 535 int ierr, i, j, format; 536 FILE *fd; 537 char *outputname; 538 Scalar *v; 539 540 ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 541 ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 542 ierr = ViewerGetFormat(viewer,&format); 543 if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 544 return 0; /* do nothing for now */ 545 } 546 else if (format == VIEWER_FORMAT_ASCII_COMMON) { 547 for ( i=0; i<a->m; i++ ) { 548 v = a->v + i; 549 fprintf(fd,"row %d:",i); 550 for ( j=0; j<a->n; j++ ) { 551 #if defined(PETSC_COMPLEX) 552 if (real(*v) != 0.0 && imag(*v) != 0.0) 553 fprintf(fd," %d %g + %g i",j,real(*v),imag(*v)); 554 else if (real(*v)) fprintf(fd," %d %g ",j,real(*v)); 555 v += a->m; 556 #else 557 if (*v) fprintf(fd," %d %g ",j,*v); 558 v += a->m; 559 #endif 560 } 561 fprintf(fd,"\n"); 562 } 563 } 564 else { 565 #if defined(PETSC_COMPLEX) 566 int allreal = 1; 567 /* determine if matrix has all real values */ 568 v = a->v; 569 for ( i=0; i<a->m*a->n; i++ ) { 570 if (imag(v[i])) { allreal = 0; break ;} 571 } 572 #endif 573 for ( i=0; i<a->m; i++ ) { 574 v = a->v + i; 575 for ( j=0; j<a->n; j++ ) { 576 #if defined(PETSC_COMPLEX) 577 if (allreal) { 578 fprintf(fd,"%6.4e ",real(*v)); v += a->m; 579 } else { 580 fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m; 581 } 582 #else 583 fprintf(fd,"%6.4e ",*v); v += a->m; 584 #endif 585 } 586 fprintf(fd,"\n"); 587 } 588 } 589 fflush(fd); 590 return 0; 591 } 592 593 #undef __FUNCTION__ 594 #define __FUNCTION__ "MatView_SeqDense_Binary" 595 static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 596 { 597 Mat_SeqDense *a = (Mat_SeqDense *) A->data; 598 int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 599 int format; 600 Scalar *v, *anonz,*vals; 601 602 ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 603 604 ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 605 if (format == VIEWER_FORMAT_BINARY_NATIVE) { 606 /* store the matrix as a dense matrix */ 607 col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens); 608 col_lens[0] = MAT_COOKIE; 609 col_lens[1] = m; 610 col_lens[2] = n; 611 col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 612 ierr = PetscBinaryWrite(fd,col_lens,4,BINARY_INT,1); CHKERRQ(ierr); 613 PetscFree(col_lens); 614 615 /* write out matrix, by rows */ 616 vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals); 617 v = a->v; 618 for ( i=0; i<m; i++ ) { 619 for ( j=0; j<n; j++ ) { 620 vals[i + j*m] = *v++; 621 } 622 } 623 ierr = PetscBinaryWrite(fd,vals,n*m,BINARY_SCALAR,0); CHKERRQ(ierr); 624 PetscFree(vals); 625 } else { 626 col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens); 627 col_lens[0] = MAT_COOKIE; 628 col_lens[1] = m; 629 col_lens[2] = n; 630 col_lens[3] = nz; 631 632 /* store lengths of each row and write (including header) to file */ 633 for ( i=0; i<m; i++ ) col_lens[4+i] = n; 634 ierr = PetscBinaryWrite(fd,col_lens,4+m,BINARY_INT,1); CHKERRQ(ierr); 635 636 /* Possibly should write in smaller increments, not whole matrix at once? */ 637 /* store column indices (zero start index) */ 638 ict = 0; 639 for ( i=0; i<m; i++ ) { 640 for ( j=0; j<n; j++ ) col_lens[ict++] = j; 641 } 642 ierr = PetscBinaryWrite(fd,col_lens,nz,BINARY_INT,0); CHKERRQ(ierr); 643 PetscFree(col_lens); 644 645 /* store nonzero values */ 646 anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz); 647 ict = 0; 648 for ( i=0; i<m; i++ ) { 649 v = a->v + i; 650 for ( j=0; j<n; j++ ) { 651 anonz[ict++] = *v; v += a->m; 652 } 653 } 654 ierr = PetscBinaryWrite(fd,anonz,nz,BINARY_SCALAR,0); CHKERRQ(ierr); 655 PetscFree(anonz); 656 } 657 return 0; 658 } 659 660 #undef __FUNCTION__ 661 #define __FUNCTION__ "MatView_SeqDense" 662 static int MatView_SeqDense(PetscObject obj,Viewer viewer) 663 { 664 Mat A = (Mat) obj; 665 Mat_SeqDense *a = (Mat_SeqDense*) A->data; 666 ViewerType vtype; 667 int ierr; 668 669 ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 670 671 if (vtype == MATLAB_VIEWER) { 672 return ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); 673 } 674 else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 675 return MatView_SeqDense_ASCII(A,viewer); 676 } 677 else if (vtype == BINARY_FILE_VIEWER) { 678 return MatView_SeqDense_Binary(A,viewer); 679 } 680 return 0; 681 } 682 683 #undef __FUNCTION__ 684 #define __FUNCTION__ "MatDestroy_SeqDense" 685 static int MatDestroy_SeqDense(PetscObject obj) 686 { 687 Mat mat = (Mat) obj; 688 Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 689 int ierr; 690 691 #if defined(PETSC_LOG) 692 PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 693 #endif 694 if (l->pivots) PetscFree(l->pivots); 695 if (!l->user_alloc) PetscFree(l->v); 696 PetscFree(l); 697 if (mat->mapping) { 698 ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr); 699 } 700 PLogObjectDestroy(mat); 701 PetscHeaderDestroy(mat); 702 return 0; 703 } 704 705 #undef __FUNCTION__ 706 #define __FUNCTION__ "MatTranspose_SeqDense" 707 static int MatTranspose_SeqDense(Mat A,Mat *matout) 708 { 709 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 710 int k, j, m, n; 711 Scalar *v, tmp; 712 713 v = mat->v; m = mat->m; n = mat->n; 714 if (matout == PETSC_NULL) { /* in place transpose */ 715 if (m != n) SETERRQ(PETSC_ERR_SUP,"Supports square matrix only in-place"); 716 for ( j=0; j<m; j++ ) { 717 for ( k=0; k<j; k++ ) { 718 tmp = v[j + k*n]; 719 v[j + k*n] = v[k + j*n]; 720 v[k + j*n] = tmp; 721 } 722 } 723 } 724 else { /* out-of-place transpose */ 725 int ierr; 726 Mat tmat; 727 Mat_SeqDense *tmatd; 728 Scalar *v2; 729 ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr); 730 tmatd = (Mat_SeqDense *) tmat->data; 731 v = mat->v; v2 = tmatd->v; 732 for ( j=0; j<n; j++ ) { 733 for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 734 } 735 ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 736 ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 737 *matout = tmat; 738 } 739 return 0; 740 } 741 742 #undef __FUNCTION__ 743 #define __FUNCTION__ "MatEqual_SeqDense" 744 static int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg) 745 { 746 Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 747 Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 748 int i; 749 Scalar *v1 = mat1->v, *v2 = mat2->v; 750 751 if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,"Both matrices should be of type MATSEQDENSE"); 752 if (mat1->m != mat2->m) {*flg = PETSC_FALSE; return 0;} 753 if (mat1->n != mat2->n) {*flg = PETSC_FALSE; return 0;} 754 for ( i=0; i<mat1->m*mat1->n; i++ ) { 755 if (*v1 != *v2) {*flg = PETSC_FALSE; return 0;} 756 v1++; v2++; 757 } 758 *flg = PETSC_TRUE; 759 return 0; 760 } 761 762 #undef __FUNCTION__ 763 #define __FUNCTION__ "MatGetDiagonal_SeqDense" 764 static int MatGetDiagonal_SeqDense(Mat A,Vec v) 765 { 766 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 767 int i, n, len; 768 Scalar *x, zero = 0.0; 769 770 VecSet(&zero,v); 771 VecGetArray(v,&x); VecGetSize(v,&n); 772 len = PetscMin(mat->m,mat->n); 773 if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 774 for ( i=0; i<len; i++ ) { 775 x[i] = mat->v[i*mat->m + i]; 776 } 777 return 0; 778 } 779 780 #undef __FUNCTION__ 781 #define __FUNCTION__ "MatDiagonalScale_SeqDense" 782 static int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 783 { 784 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 785 Scalar *l,*r,x,*v; 786 int i,j,m = mat->m, n = mat->n; 787 788 if (ll) { 789 VecGetArray(ll,&l); VecGetSize(ll,&m); 790 if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 791 for ( i=0; i<m; i++ ) { 792 x = l[i]; 793 v = mat->v + i; 794 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 795 } 796 PLogFlops(n*m); 797 } 798 if (rr) { 799 VecGetArray(rr,&r); VecGetSize(rr,&n); 800 if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 801 for ( i=0; i<n; i++ ) { 802 x = r[i]; 803 v = mat->v + i*m; 804 for ( j=0; j<m; j++ ) { (*v++) *= x;} 805 } 806 PLogFlops(n*m); 807 } 808 return 0; 809 } 810 811 #undef __FUNCTION__ 812 #define __FUNCTION__ "MatNorm_SeqDense" 813 static int MatNorm_SeqDense(Mat A,NormType type,double *norm) 814 { 815 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 816 Scalar *v = mat->v; 817 double sum = 0.0; 818 int i, j; 819 820 if (type == NORM_FROBENIUS) { 821 for (i=0; i<mat->n*mat->m; i++ ) { 822 #if defined(PETSC_COMPLEX) 823 sum += real(conj(*v)*(*v)); v++; 824 #else 825 sum += (*v)*(*v); v++; 826 #endif 827 } 828 *norm = sqrt(sum); 829 PLogFlops(2*mat->n*mat->m); 830 } 831 else if (type == NORM_1) { 832 *norm = 0.0; 833 for ( j=0; j<mat->n; j++ ) { 834 sum = 0.0; 835 for ( i=0; i<mat->m; i++ ) { 836 sum += PetscAbsScalar(*v); v++; 837 } 838 if (sum > *norm) *norm = sum; 839 } 840 PLogFlops(mat->n*mat->m); 841 } 842 else if (type == NORM_INFINITY) { 843 *norm = 0.0; 844 for ( j=0; j<mat->m; j++ ) { 845 v = mat->v + j; 846 sum = 0.0; 847 for ( i=0; i<mat->n; i++ ) { 848 sum += PetscAbsScalar(*v); v += mat->m; 849 } 850 if (sum > *norm) *norm = sum; 851 } 852 PLogFlops(mat->n*mat->m); 853 } 854 else { 855 SETERRQ(PETSC_ERR_SUP,"No two norm"); 856 } 857 return 0; 858 } 859 860 #undef __FUNCTION__ 861 #define __FUNCTION__ "MatSetOption_SeqDense" 862 static int MatSetOption_SeqDense(Mat A,MatOption op) 863 { 864 Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 865 866 if (op == MAT_ROW_ORIENTED) aij->roworiented = 1; 867 else if (op == MAT_COLUMN_ORIENTED) aij->roworiented = 0; 868 else if (op == MAT_ROWS_SORTED || 869 op == MAT_ROWS_UNSORTED || 870 op == MAT_COLUMNS_SORTED || 871 op == MAT_COLUMNS_UNSORTED || 872 op == MAT_SYMMETRIC || 873 op == MAT_STRUCTURALLY_SYMMETRIC || 874 op == MAT_NO_NEW_NONZERO_LOCATIONS || 875 op == MAT_YES_NEW_NONZERO_LOCATIONS || 876 op == MAT_NO_NEW_DIAGONALS || 877 op == MAT_YES_NEW_DIAGONALS || 878 op == MAT_IGNORE_OFF_PROCESSOR_ENTRIES) 879 PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n"); 880 else 881 {SETERRQ(PETSC_ERR_SUP,"unknown option");} 882 return 0; 883 } 884 885 #undef __FUNCTION__ 886 #define __FUNCTION__ "MatZeroEntries_SeqDense" 887 static int MatZeroEntries_SeqDense(Mat A) 888 { 889 Mat_SeqDense *l = (Mat_SeqDense *) A->data; 890 PetscMemzero(l->v,l->m*l->n*sizeof(Scalar)); 891 return 0; 892 } 893 894 #undef __FUNCTION__ 895 #define __FUNCTION__ "MatGetBlockSize_SeqDense" 896 static int MatGetBlockSize_SeqDense(Mat A,int *bs) 897 { 898 *bs = 1; 899 return 0; 900 } 901 902 #undef __FUNCTION__ 903 #define __FUNCTION__ "MatZeroRows_SeqDense" 904 static int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 905 { 906 Mat_SeqDense *l = (Mat_SeqDense *) A->data; 907 int n = l->n, i, j,ierr,N, *rows; 908 Scalar *slot; 909 910 ierr = ISGetSize(is,&N); CHKERRQ(ierr); 911 ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 912 for ( i=0; i<N; i++ ) { 913 slot = l->v + rows[i]; 914 for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 915 } 916 if (diag) { 917 for ( i=0; i<N; i++ ) { 918 slot = l->v + (n+1)*rows[i]; 919 *slot = *diag; 920 } 921 } 922 ISRestoreIndices(is,&rows); 923 return 0; 924 } 925 926 #undef __FUNCTION__ 927 #define __FUNCTION__ "MatGetSize_SeqDense" 928 static int MatGetSize_SeqDense(Mat A,int *m,int *n) 929 { 930 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 931 *m = mat->m; *n = mat->n; 932 return 0; 933 } 934 935 #undef __FUNCTION__ 936 #define __FUNCTION__ "MatGetOwnershipRange_SeqDense" 937 static int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 938 { 939 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 940 *m = 0; *n = mat->m; 941 return 0; 942 } 943 944 #undef __FUNCTION__ 945 #define __FUNCTION__ "MatGetArray_SeqDense" 946 static int MatGetArray_SeqDense(Mat A,Scalar **array) 947 { 948 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 949 *array = mat->v; 950 return 0; 951 } 952 953 #undef __FUNCTION__ 954 #define __FUNCTION__ "MatRestoreArray_SeqDense" 955 static int MatRestoreArray_SeqDense(Mat A,Scalar **array) 956 { 957 return 0; 958 } 959 960 #undef __FUNCTION__ 961 #define __FUNCTION__ "MatGetSubMatrix_SeqDense" 962 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall, 963 Mat *submat) 964 { 965 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 966 int nznew, *smap, i, j, ierr, oldcols = mat->n; 967 int *irow, *icol, nrows, ncols, *cwork; 968 Scalar *vwork, *val; 969 Mat newmat; 970 971 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 972 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 973 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 974 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 975 976 smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap); 977 cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork); 978 vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 979 PetscMemzero((char*)smap,oldcols*sizeof(int)); 980 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 981 982 /* Create and fill new matrix */ 983 ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr); 984 for (i=0; i<nrows; i++) { 985 nznew = 0; 986 val = mat->v + irow[i]; 987 for (j=0; j<oldcols; j++) { 988 if (smap[j]) { 989 cwork[nznew] = smap[j] - 1; 990 vwork[nznew++] = val[j * mat->m]; 991 } 992 } 993 ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES); 994 CHKERRQ(ierr); 995 } 996 ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 997 ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 998 999 /* Free work space */ 1000 PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 1001 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1002 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 1003 *submat = newmat; 1004 return 0; 1005 } 1006 1007 #undef __FUNCTION__ 1008 #define __FUNCTION__ "MatGetSubMatrices_SeqDense" 1009 static int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1010 Mat **B) 1011 { 1012 int ierr,i; 1013 1014 if (scall == MAT_INITIAL_MATRIX) { 1015 *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1016 } 1017 1018 for ( i=0; i<n; i++ ) { 1019 ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 1020 } 1021 return 0; 1022 } 1023 1024 #undef __FUNCTION__ 1025 #define __FUNCTION__ "MatCopy_SeqDense" 1026 static int MatCopy_SeqDense(Mat A, Mat B) 1027 { 1028 Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data; 1029 if (B->type != MATSEQDENSE) return MatCopy_Basic(A,B); 1030 if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1031 PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar)); 1032 return 0; 1033 } 1034 1035 /* -------------------------------------------------------------------*/ 1036 static struct _MatOps MatOps = {MatSetValues_SeqDense, 1037 MatGetRow_SeqDense, 1038 MatRestoreRow_SeqDense, 1039 MatMult_SeqDense, 1040 MatMultAdd_SeqDense, 1041 MatMultTrans_SeqDense, 1042 MatMultTransAdd_SeqDense, 1043 MatSolve_SeqDense, 1044 MatSolveAdd_SeqDense, 1045 MatSolveTrans_SeqDense, 1046 MatSolveTransAdd_SeqDense, 1047 MatLUFactor_SeqDense, 1048 MatCholeskyFactor_SeqDense, 1049 MatRelax_SeqDense, 1050 MatTranspose_SeqDense, 1051 MatGetInfo_SeqDense, 1052 MatEqual_SeqDense, 1053 MatGetDiagonal_SeqDense, 1054 MatDiagonalScale_SeqDense, 1055 MatNorm_SeqDense, 1056 0, 1057 0, 1058 0, 1059 MatSetOption_SeqDense, 1060 MatZeroEntries_SeqDense, 1061 MatZeroRows_SeqDense, 1062 MatLUFactorSymbolic_SeqDense, 1063 MatLUFactorNumeric_SeqDense, 1064 MatCholeskyFactorSymbolic_SeqDense, 1065 MatCholeskyFactorNumeric_SeqDense, 1066 MatGetSize_SeqDense, 1067 MatGetSize_SeqDense, 1068 MatGetOwnershipRange_SeqDense, 1069 0, 1070 0, 1071 MatGetArray_SeqDense, 1072 MatRestoreArray_SeqDense, 1073 0, 1074 MatConvertSameType_SeqDense,0,0,0,0, 1075 MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0, 1076 MatGetValues_SeqDense, 1077 MatCopy_SeqDense,0,MatScale_SeqDense, 1078 0,0,0,MatGetBlockSize_SeqDense}; 1079 1080 #undef __FUNCTION__ 1081 #define __FUNCTION__ "MatCreateSeqDense" 1082 /*@C 1083 MatCreateSeqDense - Creates a sequential dense matrix that 1084 is stored in column major order (the usual Fortran 77 manner). Many 1085 of the matrix operations use the BLAS and LAPACK routines. 1086 1087 Input Parameters: 1088 . comm - MPI communicator, set to MPI_COMM_SELF 1089 . m - number of rows 1090 . n - number of columns 1091 . data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1092 to control all matrix memory allocation. 1093 1094 Output Parameter: 1095 . A - the matrix 1096 1097 Notes: 1098 The data input variable is intended primarily for Fortran programmers 1099 who wish to allocate their own matrix memory space. Most users should 1100 set data=PETSC_NULL. 1101 1102 .keywords: dense, matrix, LAPACK, BLAS 1103 1104 .seealso: MatCreate(), MatSetValues() 1105 @*/ 1106 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A) 1107 { 1108 Mat B; 1109 Mat_SeqDense *b; 1110 int ierr,flg,size; 1111 1112 MPI_Comm_size(comm,&size); 1113 if (size > 1) SETERRQ(1,"Comm must be of size 1"); 1114 1115 *A = 0; 1116 PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQDENSE,comm); 1117 PLogObjectCreate(B); 1118 b = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b); 1119 PetscMemzero(b,sizeof(Mat_SeqDense)); 1120 PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1121 B->destroy = MatDestroy_SeqDense; 1122 B->view = MatView_SeqDense; 1123 B->factor = 0; 1124 B->mapping = 0; 1125 PLogObjectMemory(B,sizeof(struct _Mat)); 1126 B->data = (void *) b; 1127 1128 b->m = m; B->m = m; B->M = m; 1129 b->n = n; B->n = n; B->N = n; 1130 b->pivots = 0; 1131 b->roworiented = 1; 1132 if (data == PETSC_NULL) { 1133 b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v); 1134 PetscMemzero(b->v,m*n*sizeof(Scalar)); 1135 b->user_alloc = 0; 1136 PLogObjectMemory(B,n*m*sizeof(Scalar)); 1137 } 1138 else { /* user-allocated storage */ 1139 b->v = data; 1140 b->user_alloc = 1; 1141 } 1142 ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 1143 if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);} 1144 1145 *A = B; 1146 return 0; 1147 } 1148 1149 #undef __FUNCTION__ 1150 #define __FUNCTION__ "MatCreate_SeqDense" 1151 int MatCreate_SeqDense(Mat A,Mat *newmat) 1152 { 1153 Mat_SeqDense *m = (Mat_SeqDense *) A->data; 1154 return MatCreateSeqDense(A->comm,m->m,m->n,PETSC_NULL,newmat); 1155 } 1156 1157 1158 1159 1160 1161