1 /*$Id: dense.c,v 1.208 2001/09/07 20:09:16 bsmith Exp $*/ 2 /* 3 Defines the basic matrix operations for sequential dense. 4 */ 5 6 #include "src/mat/impls/dense/seq/dense.h" 7 #include "petscblaslapack.h" 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "MatAXPY_SeqDense" 11 int MatAXPY_SeqDense(const PetscScalar *alpha,Mat X,Mat Y,MatStructure str) 12 { 13 Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data; 14 int N = X->m*X->n,m=X->m,ldax=x->lda,lday=y->lda, j,one = 1; 15 16 PetscFunctionBegin; 17 if (X->m != Y->m || X->n != Y->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 18 if (ldax>m || lday>m) { 19 for (j=0; j<X->n; j++) { 20 BLaxpy_(&m,(PetscScalar*)alpha,x->v+j*ldax,&one,y->v+j*lday,&one); 21 } 22 } else { 23 BLaxpy_(&N,(PetscScalar*)alpha,x->v,&one,y->v,&one); 24 } 25 PetscLogFlops(2*N-1); 26 PetscFunctionReturn(0); 27 } 28 29 #undef __FUNCT__ 30 #define __FUNCT__ "MatGetInfo_SeqDense" 31 int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 32 { 33 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 34 int i,N = A->m*A->n,count = 0; 35 PetscScalar *v = a->v; 36 37 PetscFunctionBegin; 38 for (i=0; i<N; i++) {if (*v != 0.0) count++; v++;} 39 40 info->rows_global = (double)A->m; 41 info->columns_global = (double)A->n; 42 info->rows_local = (double)A->m; 43 info->columns_local = (double)A->n; 44 info->block_size = 1.0; 45 info->nz_allocated = (double)N; 46 info->nz_used = (double)count; 47 info->nz_unneeded = (double)(N-count); 48 info->assemblies = (double)A->num_ass; 49 info->mallocs = 0; 50 info->memory = A->mem; 51 info->fill_ratio_given = 0; 52 info->fill_ratio_needed = 0; 53 info->factor_mallocs = 0; 54 55 PetscFunctionReturn(0); 56 } 57 58 #undef __FUNCT__ 59 #define __FUNCT__ "MatScale_SeqDense" 60 int MatScale_SeqDense(const PetscScalar *alpha,Mat A) 61 { 62 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 63 int one = 1,lda = a->lda,j,nz; 64 65 PetscFunctionBegin; 66 if (lda>A->m) { 67 nz = A->m; 68 for (j=0; j<A->n; j++) { 69 BLscal_(&nz,(PetscScalar*)alpha,a->v+j*lda,&one); 70 } 71 } else { 72 nz = A->m*A->n; 73 BLscal_(&nz,(PetscScalar*)alpha,a->v,&one); 74 } 75 PetscLogFlops(nz); 76 PetscFunctionReturn(0); 77 } 78 79 /* ---------------------------------------------------------------*/ 80 /* COMMENT: I have chosen to hide row permutation in the pivots, 81 rather than put it in the Mat->row slot.*/ 82 #undef __FUNCT__ 83 #define __FUNCT__ "MatLUFactor_SeqDense" 84 int MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo) 85 { 86 #if defined(PETSC_MISSING_LAPACK_GETRF) 87 PetscFunctionBegin; 88 SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 89 #else 90 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 91 int info,ierr; 92 93 PetscFunctionBegin; 94 if (!mat->pivots) { 95 ierr = PetscMalloc((A->m+1)*sizeof(int),&mat->pivots);CHKERRQ(ierr); 96 PetscLogObjectMemory(A,A->m*sizeof(int)); 97 } 98 A->factor = FACTOR_LU; 99 if (!A->m || !A->n) PetscFunctionReturn(0); 100 LAgetrf_(&A->m,&A->n,mat->v,&mat->lda,mat->pivots,&info); 101 if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization"); 102 if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 103 PetscLogFlops((2*A->n*A->n*A->n)/3); 104 #endif 105 PetscFunctionReturn(0); 106 } 107 108 #undef __FUNCT__ 109 #define __FUNCT__ "MatDuplicate_SeqDense" 110 int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 111 { 112 Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 113 int lda = mat->lda,j,m,ierr; 114 Mat newi; 115 116 PetscFunctionBegin; 117 ierr = MatCreateSeqDense(A->comm,A->m,A->n,PETSC_NULL,&newi);CHKERRQ(ierr); 118 if (cpvalues == MAT_COPY_VALUES) { 119 l = (Mat_SeqDense*)newi->data; 120 if (lda>A->m) { 121 m = A->m; 122 for (j=0; j<A->n; j++) { 123 ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 124 } 125 } else { 126 ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 127 } 128 } 129 newi->assembled = PETSC_TRUE; 130 *newmat = newi; 131 PetscFunctionReturn(0); 132 } 133 134 #undef __FUNCT__ 135 #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 136 int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact) 137 { 138 int ierr; 139 140 PetscFunctionBegin; 141 ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 142 PetscFunctionReturn(0); 143 } 144 145 #undef __FUNCT__ 146 #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 147 int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 148 { 149 Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data; 150 int lda1=mat->lda,lda2=l->lda, m=A->m,n=A->n, j,ierr; 151 MatFactorInfo info; 152 153 PetscFunctionBegin; 154 /* copy the numerical values */ 155 if (lda1>m || lda2>m ) { 156 for (j=0; j<n; j++) { 157 ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar)); CHKERRQ(ierr); 158 } 159 } else { 160 ierr = PetscMemcpy(l->v,mat->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 161 } 162 (*fact)->factor = 0; 163 ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr); 164 PetscFunctionReturn(0); 165 } 166 167 #undef __FUNCT__ 168 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 169 int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact) 170 { 171 int ierr; 172 173 PetscFunctionBegin; 174 ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr); 175 PetscFunctionReturn(0); 176 } 177 178 #undef __FUNCT__ 179 #define __FUNCT__ "MatCholeskyFactor_SeqDense" 180 int MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo) 181 { 182 #if defined(PETSC_MISSING_LAPACK_POTRF) 183 PetscFunctionBegin; 184 SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 185 #else 186 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 187 int info,ierr; 188 189 PetscFunctionBegin; 190 if (mat->pivots) { 191 ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 192 PetscLogObjectMemory(A,-A->m*sizeof(int)); 193 mat->pivots = 0; 194 } 195 if (!A->m || !A->n) PetscFunctionReturn(0); 196 LApotrf_("L",&A->n,mat->v,&mat->lda,&info); 197 if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %d",info-1); 198 A->factor = FACTOR_CHOLESKY; 199 PetscLogFlops((A->n*A->n*A->n)/3); 200 #endif 201 PetscFunctionReturn(0); 202 } 203 204 #undef __FUNCT__ 205 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 206 int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 207 { 208 int ierr; 209 MatFactorInfo info; 210 211 PetscFunctionBegin; 212 info.fill = 1.0; 213 ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr); 214 PetscFunctionReturn(0); 215 } 216 217 #undef __FUNCT__ 218 #define __FUNCT__ "MatSolve_SeqDense" 219 int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 220 { 221 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 222 int one = 1,info,ierr; 223 PetscScalar *x,*y; 224 225 PetscFunctionBegin; 226 if (!A->m || !A->n) PetscFunctionReturn(0); 227 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 228 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 229 ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 230 if (A->factor == FACTOR_LU) { 231 #if defined(PETSC_MISSING_LAPACK_GETRS) 232 SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 233 #else 234 LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 235 if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve"); 236 #endif 237 } else if (A->factor == FACTOR_CHOLESKY){ 238 #if defined(PETSC_MISSING_LAPACK_POTRS) 239 SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 240 #else 241 LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 242 if (info) SETERRQ(PETSC_ERR_LIB,"MBad solve"); 243 #endif 244 } 245 else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 246 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 247 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 248 PetscLogFlops(2*A->n*A->n - A->n); 249 PetscFunctionReturn(0); 250 } 251 252 #undef __FUNCT__ 253 #define __FUNCT__ "MatSolveTranspose_SeqDense" 254 int MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 255 { 256 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 257 int ierr,one = 1,info; 258 PetscScalar *x,*y; 259 260 PetscFunctionBegin; 261 if (!A->m || !A->n) PetscFunctionReturn(0); 262 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 263 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 264 ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 265 /* assume if pivots exist then use LU; else Cholesky */ 266 if (mat->pivots) { 267 #if defined(PETSC_MISSING_LAPACK_GETRS) 268 SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 269 #else 270 LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 271 if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 272 #endif 273 } else { 274 #if defined(PETSC_MISSING_LAPACK_POTRS) 275 SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 276 #else 277 LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 278 if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 279 #endif 280 } 281 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 282 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 283 PetscLogFlops(2*A->n*A->n - A->n); 284 PetscFunctionReturn(0); 285 } 286 287 #undef __FUNCT__ 288 #define __FUNCT__ "MatSolveAdd_SeqDense" 289 int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 290 { 291 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 292 int ierr,one = 1,info; 293 PetscScalar *x,*y,sone = 1.0; 294 Vec tmp = 0; 295 296 PetscFunctionBegin; 297 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 298 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 299 if (!A->m || !A->n) PetscFunctionReturn(0); 300 if (yy == zz) { 301 ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 302 PetscLogObjectParent(A,tmp); 303 ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 304 } 305 ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 306 /* assume if pivots exist then use LU; else Cholesky */ 307 if (mat->pivots) { 308 #if defined(PETSC_MISSING_LAPACK_GETRS) 309 SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 310 #else 311 LAgetrs_("N",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 312 if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 313 #endif 314 } else { 315 #if defined(PETSC_MISSING_LAPACK_POTRS) 316 SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 317 #else 318 LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 319 if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 320 #endif 321 } 322 if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 323 else {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);} 324 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 325 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 326 PetscLogFlops(2*A->n*A->n); 327 PetscFunctionReturn(0); 328 } 329 330 #undef __FUNCT__ 331 #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 332 int MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 333 { 334 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 335 int one = 1,info,ierr; 336 PetscScalar *x,*y,sone = 1.0; 337 Vec tmp; 338 339 PetscFunctionBegin; 340 if (!A->m || !A->n) PetscFunctionReturn(0); 341 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 342 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 343 if (yy == zz) { 344 ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 345 PetscLogObjectParent(A,tmp); 346 ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 347 } 348 ierr = PetscMemcpy(y,x,A->m*sizeof(PetscScalar));CHKERRQ(ierr); 349 /* assume if pivots exist then use LU; else Cholesky */ 350 if (mat->pivots) { 351 #if defined(PETSC_MISSING_LAPACK_GETRS) 352 SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 353 #else 354 LAgetrs_("T",&A->m,&one,mat->v,&mat->lda,mat->pivots,y,&A->m,&info); 355 if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 356 #endif 357 } else { 358 #if defined(PETSC_MISSING_LAPACK_POTRS) 359 SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 360 #else 361 LApotrs_("L",&A->m,&one,mat->v,&mat->lda,y,&A->m,&info); 362 if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 363 #endif 364 } 365 if (tmp) { 366 ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); 367 ierr = VecDestroy(tmp);CHKERRQ(ierr); 368 } else { 369 ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr); 370 } 371 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 372 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 373 PetscLogFlops(2*A->n*A->n); 374 PetscFunctionReturn(0); 375 } 376 /* ------------------------------------------------------------------*/ 377 #undef __FUNCT__ 378 #define __FUNCT__ "MatRelax_SeqDense" 379 int MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag, 380 PetscReal shift,int its,int lits,Vec xx) 381 { 382 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 383 PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 384 int ierr,m = A->m,i; 385 #if !defined(PETSC_USE_COMPLEX) 386 int o = 1; 387 #endif 388 389 PetscFunctionBegin; 390 if (flag & SOR_ZERO_INITIAL_GUESS) { 391 /* this is a hack fix, should have another version without the second BLdot */ 392 ierr = VecSet(&zero,xx);CHKERRQ(ierr); 393 } 394 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 395 ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 396 its = its*lits; 397 if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %d and local its %d both positive",its,lits); 398 while (its--) { 399 if (flag & SOR_FORWARD_SWEEP){ 400 for (i=0; i<m; i++) { 401 #if defined(PETSC_USE_COMPLEX) 402 /* cannot use BLAS dot for complex because compiler/linker is 403 not happy about returning a double complex */ 404 int _i; 405 PetscScalar sum = b[i]; 406 for (_i=0; _i<m; _i++) { 407 sum -= PetscConj(v[i+_i*m])*x[_i]; 408 } 409 xt = sum; 410 #else 411 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 412 #endif 413 x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 414 } 415 } 416 if (flag & SOR_BACKWARD_SWEEP) { 417 for (i=m-1; i>=0; i--) { 418 #if defined(PETSC_USE_COMPLEX) 419 /* cannot use BLAS dot for complex because compiler/linker is 420 not happy about returning a double complex */ 421 int _i; 422 PetscScalar sum = b[i]; 423 for (_i=0; _i<m; _i++) { 424 sum -= PetscConj(v[i+_i*m])*x[_i]; 425 } 426 xt = sum; 427 #else 428 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 429 #endif 430 x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 431 } 432 } 433 } 434 ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 435 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 436 PetscFunctionReturn(0); 437 } 438 439 /* -----------------------------------------------------------------*/ 440 #undef __FUNCT__ 441 #define __FUNCT__ "MatMultTranspose_SeqDense" 442 int MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 443 { 444 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 445 PetscScalar *v = mat->v,*x,*y; 446 int ierr,_One=1; 447 PetscScalar _DOne=1.0,_DZero=0.0; 448 449 PetscFunctionBegin; 450 if (!A->m || !A->n) PetscFunctionReturn(0); 451 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 452 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 453 LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 454 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 455 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 456 PetscLogFlops(2*A->m*A->n - A->n); 457 PetscFunctionReturn(0); 458 } 459 460 #undef __FUNCT__ 461 #define __FUNCT__ "MatMult_SeqDense" 462 int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 463 { 464 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 465 PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 466 int ierr,_One=1; 467 468 PetscFunctionBegin; 469 if (!A->m || !A->n) PetscFunctionReturn(0); 470 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 471 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 472 LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 473 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 474 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 475 PetscLogFlops(2*A->m*A->n - A->m); 476 PetscFunctionReturn(0); 477 } 478 479 #undef __FUNCT__ 480 #define __FUNCT__ "MatMultAdd_SeqDense" 481 int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 482 { 483 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 484 PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 485 int ierr,_One=1; 486 487 PetscFunctionBegin; 488 if (!A->m || !A->n) PetscFunctionReturn(0); 489 if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 490 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 491 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 492 LAgemv_("N",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 493 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 494 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 495 PetscLogFlops(2*A->m*A->n); 496 PetscFunctionReturn(0); 497 } 498 499 #undef __FUNCT__ 500 #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 501 int MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 502 { 503 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 504 PetscScalar *v = mat->v,*x,*y; 505 int ierr,_One=1; 506 PetscScalar _DOne=1.0; 507 508 PetscFunctionBegin; 509 if (!A->m || !A->n) PetscFunctionReturn(0); 510 if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 511 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 512 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 513 LAgemv_("T",&(A->m),&(A->n),&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 514 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 515 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 516 PetscLogFlops(2*A->m*A->n); 517 PetscFunctionReturn(0); 518 } 519 520 /* -----------------------------------------------------------------*/ 521 #undef __FUNCT__ 522 #define __FUNCT__ "MatGetRow_SeqDense" 523 int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals) 524 { 525 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 526 PetscScalar *v; 527 int i,ierr; 528 529 PetscFunctionBegin; 530 *ncols = A->n; 531 if (cols) { 532 ierr = PetscMalloc((A->n+1)*sizeof(int),cols);CHKERRQ(ierr); 533 for (i=0; i<A->n; i++) (*cols)[i] = i; 534 } 535 if (vals) { 536 ierr = PetscMalloc((A->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 537 v = mat->v + row; 538 for (i=0; i<A->n; i++) {(*vals)[i] = *v; v += mat->lda;} 539 } 540 PetscFunctionReturn(0); 541 } 542 543 #undef __FUNCT__ 544 #define __FUNCT__ "MatRestoreRow_SeqDense" 545 int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,PetscScalar **vals) 546 { 547 int ierr; 548 PetscFunctionBegin; 549 if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 550 if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 551 PetscFunctionReturn(0); 552 } 553 /* ----------------------------------------------------------------*/ 554 #undef __FUNCT__ 555 #define __FUNCT__ "MatSetValues_SeqDense" 556 int MatSetValues_SeqDense(Mat A,int m,const int indexm[],int n,const int indexn[],const PetscScalar v[],InsertMode addv) 557 { 558 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 559 int i,j,idx=0; 560 561 PetscFunctionBegin; 562 if (!mat->roworiented) { 563 if (addv == INSERT_VALUES) { 564 for (j=0; j<n; j++) { 565 if (indexn[j] < 0) {idx += m; continue;} 566 #if defined(PETSC_USE_BOPT_g) 567 if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[j],A->n-1); 568 #endif 569 for (i=0; i<m; i++) { 570 if (indexm[i] < 0) {idx++; continue;} 571 #if defined(PETSC_USE_BOPT_g) 572 if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[i],A->m-1); 573 #endif 574 mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 575 } 576 } 577 } else { 578 for (j=0; j<n; j++) { 579 if (indexn[j] < 0) {idx += m; continue;} 580 #if defined(PETSC_USE_BOPT_g) 581 if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[j],A->n-1); 582 #endif 583 for (i=0; i<m; i++) { 584 if (indexm[i] < 0) {idx++; continue;} 585 #if defined(PETSC_USE_BOPT_g) 586 if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[i],A->m-1); 587 #endif 588 mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 589 } 590 } 591 } 592 } else { 593 if (addv == INSERT_VALUES) { 594 for (i=0; i<m; i++) { 595 if (indexm[i] < 0) { idx += n; continue;} 596 #if defined(PETSC_USE_BOPT_g) 597 if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[j],A->m-1); 598 #endif 599 for (j=0; j<n; j++) { 600 if (indexn[j] < 0) { idx++; continue;} 601 #if defined(PETSC_USE_BOPT_g) 602 if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[i],A->n-1); 603 #endif 604 mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 605 } 606 } 607 } else { 608 for (i=0; i<m; i++) { 609 if (indexm[i] < 0) { idx += n; continue;} 610 #if defined(PETSC_USE_BOPT_g) 611 if (indexm[i] >= A->m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %d max %d",indexm[j],A->m-1); 612 #endif 613 for (j=0; j<n; j++) { 614 if (indexn[j] < 0) { idx++; continue;} 615 #if defined(PETSC_USE_BOPT_g) 616 if (indexn[j] >= A->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %d max %d",indexn[i],A->n-1); 617 #endif 618 mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 619 } 620 } 621 } 622 } 623 PetscFunctionReturn(0); 624 } 625 626 #undef __FUNCT__ 627 #define __FUNCT__ "MatGetValues_SeqDense" 628 int MatGetValues_SeqDense(Mat A,int m,const int indexm[],int n,const int indexn[],PetscScalar v[]) 629 { 630 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 631 int i,j; 632 PetscScalar *vpt = v; 633 634 PetscFunctionBegin; 635 /* row-oriented output */ 636 for (i=0; i<m; i++) { 637 for (j=0; j<n; j++) { 638 *vpt++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 639 } 640 } 641 PetscFunctionReturn(0); 642 } 643 644 /* -----------------------------------------------------------------*/ 645 646 #include "petscsys.h" 647 648 #undef __FUNCT__ 649 #define __FUNCT__ "MatLoad_SeqDense" 650 int MatLoad_SeqDense(PetscViewer viewer,const MatType type,Mat *A) 651 { 652 Mat_SeqDense *a; 653 Mat B; 654 int *scols,i,j,nz,ierr,fd,header[4],size; 655 int *rowlengths = 0,M,N,*cols; 656 PetscScalar *vals,*svals,*v,*w; 657 MPI_Comm comm = ((PetscObject)viewer)->comm; 658 659 PetscFunctionBegin; 660 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 661 if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 662 ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 663 ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 664 if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 665 M = header[1]; N = header[2]; nz = header[3]; 666 667 if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 668 ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 669 B = *A; 670 a = (Mat_SeqDense*)B->data; 671 v = a->v; 672 /* Allocate some temp space to read in the values and then flip them 673 from row major to column major */ 674 ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 675 /* read in nonzero values */ 676 ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 677 /* now flip the values and store them in the matrix*/ 678 for (j=0; j<N; j++) { 679 for (i=0; i<M; i++) { 680 *v++ =w[i*N+j]; 681 } 682 } 683 ierr = PetscFree(w);CHKERRQ(ierr); 684 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 685 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 686 } else { 687 /* read row lengths */ 688 ierr = PetscMalloc((M+1)*sizeof(int),&rowlengths);CHKERRQ(ierr); 689 ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 690 691 /* create our matrix */ 692 ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 693 B = *A; 694 a = (Mat_SeqDense*)B->data; 695 v = a->v; 696 697 /* read column indices and nonzeros */ 698 ierr = PetscMalloc((nz+1)*sizeof(int),&scols);CHKERRQ(ierr); 699 cols = scols; 700 ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 701 ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 702 vals = svals; 703 ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 704 705 /* insert into matrix */ 706 for (i=0; i<M; i++) { 707 for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 708 svals += rowlengths[i]; scols += rowlengths[i]; 709 } 710 ierr = PetscFree(vals);CHKERRQ(ierr); 711 ierr = PetscFree(cols);CHKERRQ(ierr); 712 ierr = PetscFree(rowlengths);CHKERRQ(ierr); 713 714 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 715 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 716 } 717 PetscFunctionReturn(0); 718 } 719 720 #include "petscsys.h" 721 722 #undef __FUNCT__ 723 #define __FUNCT__ "MatView_SeqDense_ASCII" 724 static int MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 725 { 726 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 727 int ierr,i,j; 728 char *name; 729 PetscScalar *v; 730 PetscViewerFormat format; 731 732 PetscFunctionBegin; 733 ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 734 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 735 if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 736 PetscFunctionReturn(0); /* do nothing for now */ 737 } else if (format == PETSC_VIEWER_ASCII_COMMON) { 738 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 739 for (i=0; i<A->m; i++) { 740 v = a->v + i; 741 ierr = PetscViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr); 742 for (j=0; j<A->n; j++) { 743 #if defined(PETSC_USE_COMPLEX) 744 if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 745 ierr = PetscViewerASCIIPrintf(viewer," (%d, %g + %g i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 746 } else if (PetscRealPart(*v)) { 747 ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,PetscRealPart(*v));CHKERRQ(ierr); 748 } 749 #else 750 if (*v) { 751 ierr = PetscViewerASCIIPrintf(viewer," (%d, %g) ",j,*v);CHKERRQ(ierr); 752 } 753 #endif 754 v += a->lda; 755 } 756 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 757 } 758 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 759 } else { 760 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 761 #if defined(PETSC_USE_COMPLEX) 762 PetscTruth allreal = PETSC_TRUE; 763 /* determine if matrix has all real values */ 764 v = a->v; 765 for (i=0; i<A->m*A->n; i++) { 766 if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 767 } 768 #endif 769 if (format == PETSC_VIEWER_ASCII_MATLAB) { 770 ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 771 ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %d %d \n",A->m,A->n);CHKERRQ(ierr); 772 ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%d,%d);\n",name,A->m,A->n);CHKERRQ(ierr); 773 ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 774 } 775 776 for (i=0; i<A->m; i++) { 777 v = a->v + i; 778 for (j=0; j<A->n; j++) { 779 #if defined(PETSC_USE_COMPLEX) 780 if (allreal) { 781 ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",PetscRealPart(*v));CHKERRQ(ierr); 782 } else { 783 ierr = PetscViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 784 } 785 #else 786 ierr = PetscViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); 787 #endif 788 v += a->lda; 789 } 790 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 791 } 792 if (format == PETSC_VIEWER_ASCII_MATLAB) { 793 ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 794 } 795 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 796 } 797 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 798 PetscFunctionReturn(0); 799 } 800 801 #undef __FUNCT__ 802 #define __FUNCT__ "MatView_SeqDense_Binary" 803 static int MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 804 { 805 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 806 int ict,j,n = A->n,m = A->m,i,fd,*col_lens,ierr,nz = m*n; 807 PetscScalar *v,*anonz,*vals; 808 PetscViewerFormat format; 809 810 PetscFunctionBegin; 811 ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 812 813 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 814 if (format == PETSC_VIEWER_BINARY_NATIVE) { 815 /* store the matrix as a dense matrix */ 816 ierr = PetscMalloc(4*sizeof(int),&col_lens);CHKERRQ(ierr); 817 col_lens[0] = MAT_FILE_COOKIE; 818 col_lens[1] = m; 819 col_lens[2] = n; 820 col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 821 ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr); 822 ierr = PetscFree(col_lens);CHKERRQ(ierr); 823 824 /* write out matrix, by rows */ 825 ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 826 v = a->v; 827 for (i=0; i<m; i++) { 828 for (j=0; j<n; j++) { 829 vals[i + j*m] = *v++; 830 } 831 } 832 ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr); 833 ierr = PetscFree(vals);CHKERRQ(ierr); 834 } else { 835 ierr = PetscMalloc((4+nz)*sizeof(int),&col_lens);CHKERRQ(ierr); 836 col_lens[0] = MAT_FILE_COOKIE; 837 col_lens[1] = m; 838 col_lens[2] = n; 839 col_lens[3] = nz; 840 841 /* store lengths of each row and write (including header) to file */ 842 for (i=0; i<m; i++) col_lens[4+i] = n; 843 ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr); 844 845 /* Possibly should write in smaller increments, not whole matrix at once? */ 846 /* store column indices (zero start index) */ 847 ict = 0; 848 for (i=0; i<m; i++) { 849 for (j=0; j<n; j++) col_lens[ict++] = j; 850 } 851 ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr); 852 ierr = PetscFree(col_lens);CHKERRQ(ierr); 853 854 /* store nonzero values */ 855 ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 856 ict = 0; 857 for (i=0; i<m; i++) { 858 v = a->v + i; 859 for (j=0; j<n; j++) { 860 anonz[ict++] = *v; v += a->lda; 861 } 862 } 863 ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr); 864 ierr = PetscFree(anonz);CHKERRQ(ierr); 865 } 866 PetscFunctionReturn(0); 867 } 868 869 #undef __FUNCT__ 870 #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 871 int MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 872 { 873 Mat A = (Mat) Aa; 874 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 875 int m = A->m,n = A->n,color,i,j,ierr; 876 PetscScalar *v = a->v; 877 PetscViewer viewer; 878 PetscDraw popup; 879 PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 880 PetscViewerFormat format; 881 882 PetscFunctionBegin; 883 884 ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 885 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 886 ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 887 888 /* Loop over matrix elements drawing boxes */ 889 if (format != PETSC_VIEWER_DRAW_CONTOUR) { 890 /* Blue for negative and Red for positive */ 891 color = PETSC_DRAW_BLUE; 892 for(j = 0; j < n; j++) { 893 x_l = j; 894 x_r = x_l + 1.0; 895 for(i = 0; i < m; i++) { 896 y_l = m - i - 1.0; 897 y_r = y_l + 1.0; 898 #if defined(PETSC_USE_COMPLEX) 899 if (PetscRealPart(v[j*m+i]) > 0.) { 900 color = PETSC_DRAW_RED; 901 } else if (PetscRealPart(v[j*m+i]) < 0.) { 902 color = PETSC_DRAW_BLUE; 903 } else { 904 continue; 905 } 906 #else 907 if (v[j*m+i] > 0.) { 908 color = PETSC_DRAW_RED; 909 } else if (v[j*m+i] < 0.) { 910 color = PETSC_DRAW_BLUE; 911 } else { 912 continue; 913 } 914 #endif 915 ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 916 } 917 } 918 } else { 919 /* use contour shading to indicate magnitude of values */ 920 /* first determine max of all nonzero values */ 921 for(i = 0; i < m*n; i++) { 922 if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 923 } 924 scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 925 ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 926 if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 927 for(j = 0; j < n; j++) { 928 x_l = j; 929 x_r = x_l + 1.0; 930 for(i = 0; i < m; i++) { 931 y_l = m - i - 1.0; 932 y_r = y_l + 1.0; 933 color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 934 ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 935 } 936 } 937 } 938 PetscFunctionReturn(0); 939 } 940 941 #undef __FUNCT__ 942 #define __FUNCT__ "MatView_SeqDense_Draw" 943 int MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 944 { 945 PetscDraw draw; 946 PetscTruth isnull; 947 PetscReal xr,yr,xl,yl,h,w; 948 int ierr; 949 950 PetscFunctionBegin; 951 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 952 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 953 if (isnull == PETSC_TRUE) PetscFunctionReturn(0); 954 955 ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 956 xr = A->n; yr = A->m; h = yr/10.0; w = xr/10.0; 957 xr += w; yr += h; xl = -w; yl = -h; 958 ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 959 ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 960 ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 961 PetscFunctionReturn(0); 962 } 963 964 #undef __FUNCT__ 965 #define __FUNCT__ "MatView_SeqDense" 966 int MatView_SeqDense(Mat A,PetscViewer viewer) 967 { 968 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 969 int ierr; 970 PetscTruth issocket,isascii,isbinary,isdraw; 971 972 PetscFunctionBegin; 973 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 974 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 975 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 976 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 977 978 if (issocket) { 979 if (a->lda>A->m) SETERRQ(1,"Case can not handle LDA"); 980 ierr = PetscViewerSocketPutScalar(viewer,A->m,A->n,a->v);CHKERRQ(ierr); 981 } else if (isascii) { 982 ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 983 } else if (isbinary) { 984 ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 985 } else if (isdraw) { 986 ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 987 } else { 988 SETERRQ1(1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 989 } 990 PetscFunctionReturn(0); 991 } 992 993 #undef __FUNCT__ 994 #define __FUNCT__ "MatDestroy_SeqDense" 995 int MatDestroy_SeqDense(Mat mat) 996 { 997 Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 998 int ierr; 999 1000 PetscFunctionBegin; 1001 #if defined(PETSC_USE_LOG) 1002 PetscLogObjectState((PetscObject)mat,"Rows %d Cols %d",mat->m,mat->n); 1003 #endif 1004 if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);} 1005 if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1006 ierr = PetscFree(l);CHKERRQ(ierr); 1007 PetscFunctionReturn(0); 1008 } 1009 1010 #undef __FUNCT__ 1011 #define __FUNCT__ "MatTranspose_SeqDense" 1012 int MatTranspose_SeqDense(Mat A,Mat *matout) 1013 { 1014 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1015 int k,j,m,n,M,ierr; 1016 PetscScalar *v,tmp; 1017 1018 PetscFunctionBegin; 1019 v = mat->v; m = A->m; M = mat->lda; n = A->n; 1020 if (!matout) { /* in place transpose */ 1021 if (m != n) { 1022 SETERRQ(1,"Can not transpose non-square matrix in place"); 1023 } else { 1024 for (j=0; j<m; j++) { 1025 for (k=0; k<j; k++) { 1026 tmp = v[j + k*M]; 1027 v[j + k*M] = v[k + j*M]; 1028 v[k + j*M] = tmp; 1029 } 1030 } 1031 } 1032 } else { /* out-of-place transpose */ 1033 Mat tmat; 1034 Mat_SeqDense *tmatd; 1035 PetscScalar *v2; 1036 1037 ierr = MatCreateSeqDense(A->comm,A->n,A->m,PETSC_NULL,&tmat);CHKERRQ(ierr); 1038 tmatd = (Mat_SeqDense*)tmat->data; 1039 v = mat->v; v2 = tmatd->v; 1040 for (j=0; j<n; j++) { 1041 for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1042 } 1043 ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1044 ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1045 *matout = tmat; 1046 } 1047 PetscFunctionReturn(0); 1048 } 1049 1050 #undef __FUNCT__ 1051 #define __FUNCT__ "MatEqual_SeqDense" 1052 int MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1053 { 1054 Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1055 Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 1056 int i,j; 1057 PetscScalar *v1 = mat1->v,*v2 = mat2->v; 1058 1059 PetscFunctionBegin; 1060 if (A1->m != A2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1061 if (A1->n != A2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1062 for (i=0; i<A1->m; i++) { 1063 v1 = mat1->v+i; v2 = mat2->v+i; 1064 for (j=0; j<A1->n; j++) { 1065 if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1066 v1 += mat1->lda; v2 += mat2->lda; 1067 } 1068 } 1069 *flg = PETSC_TRUE; 1070 PetscFunctionReturn(0); 1071 } 1072 1073 #undef __FUNCT__ 1074 #define __FUNCT__ "MatGetDiagonal_SeqDense" 1075 int MatGetDiagonal_SeqDense(Mat A,Vec v) 1076 { 1077 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1078 int ierr,i,n,len; 1079 PetscScalar *x,zero = 0.0; 1080 1081 PetscFunctionBegin; 1082 ierr = VecSet(&zero,v);CHKERRQ(ierr); 1083 ierr = VecGetSize(v,&n);CHKERRQ(ierr); 1084 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1085 len = PetscMin(A->m,A->n); 1086 if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 1087 for (i=0; i<len; i++) { 1088 x[i] = mat->v[i*mat->lda + i]; 1089 } 1090 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1091 PetscFunctionReturn(0); 1092 } 1093 1094 #undef __FUNCT__ 1095 #define __FUNCT__ "MatDiagonalScale_SeqDense" 1096 int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1097 { 1098 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1099 PetscScalar *l,*r,x,*v; 1100 int ierr,i,j,m = A->m,n = A->n; 1101 1102 PetscFunctionBegin; 1103 if (ll) { 1104 ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 1105 ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1106 if (m != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1107 for (i=0; i<m; i++) { 1108 x = l[i]; 1109 v = mat->v + i; 1110 for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1111 } 1112 ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1113 PetscLogFlops(n*m); 1114 } 1115 if (rr) { 1116 ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 1117 ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1118 if (n != A->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1119 for (i=0; i<n; i++) { 1120 x = r[i]; 1121 v = mat->v + i*m; 1122 for (j=0; j<m; j++) { (*v++) *= x;} 1123 } 1124 ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1125 PetscLogFlops(n*m); 1126 } 1127 PetscFunctionReturn(0); 1128 } 1129 1130 #undef __FUNCT__ 1131 #define __FUNCT__ "MatNorm_SeqDense" 1132 int MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1133 { 1134 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1135 PetscScalar *v = mat->v; 1136 PetscReal sum = 0.0; 1137 int lda=mat->lda,m=A->m,i,j; 1138 1139 PetscFunctionBegin; 1140 if (type == NORM_FROBENIUS) { 1141 if (lda>m) { 1142 for (j=0; j<A->n; j++) { 1143 v = mat->v+j*lda; 1144 for (i=0; i<m; i++) { 1145 #if defined(PETSC_USE_COMPLEX) 1146 sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1147 #else 1148 sum += (*v)*(*v); v++; 1149 #endif 1150 } 1151 } 1152 } else { 1153 for (i=0; i<A->n*A->m; i++) { 1154 #if defined(PETSC_USE_COMPLEX) 1155 sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1156 #else 1157 sum += (*v)*(*v); v++; 1158 #endif 1159 } 1160 } 1161 *nrm = sqrt(sum); 1162 PetscLogFlops(2*A->n*A->m); 1163 } else if (type == NORM_1) { 1164 *nrm = 0.0; 1165 for (j=0; j<A->n; j++) { 1166 v = mat->v + j*mat->lda; 1167 sum = 0.0; 1168 for (i=0; i<A->m; i++) { 1169 sum += PetscAbsScalar(*v); v++; 1170 } 1171 if (sum > *nrm) *nrm = sum; 1172 } 1173 PetscLogFlops(A->n*A->m); 1174 } else if (type == NORM_INFINITY) { 1175 *nrm = 0.0; 1176 for (j=0; j<A->m; j++) { 1177 v = mat->v + j; 1178 sum = 0.0; 1179 for (i=0; i<A->n; i++) { 1180 sum += PetscAbsScalar(*v); v += mat->lda; 1181 } 1182 if (sum > *nrm) *nrm = sum; 1183 } 1184 PetscLogFlops(A->n*A->m); 1185 } else { 1186 SETERRQ(PETSC_ERR_SUP,"No two norm"); 1187 } 1188 PetscFunctionReturn(0); 1189 } 1190 1191 #undef __FUNCT__ 1192 #define __FUNCT__ "MatSetOption_SeqDense" 1193 int MatSetOption_SeqDense(Mat A,MatOption op) 1194 { 1195 Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 1196 1197 PetscFunctionBegin; 1198 switch (op) { 1199 case MAT_ROW_ORIENTED: 1200 aij->roworiented = PETSC_TRUE; 1201 break; 1202 case MAT_COLUMN_ORIENTED: 1203 aij->roworiented = PETSC_FALSE; 1204 break; 1205 case MAT_ROWS_SORTED: 1206 case MAT_ROWS_UNSORTED: 1207 case MAT_COLUMNS_SORTED: 1208 case MAT_COLUMNS_UNSORTED: 1209 case MAT_NO_NEW_NONZERO_LOCATIONS: 1210 case MAT_YES_NEW_NONZERO_LOCATIONS: 1211 case MAT_NEW_NONZERO_LOCATION_ERR: 1212 case MAT_NO_NEW_DIAGONALS: 1213 case MAT_YES_NEW_DIAGONALS: 1214 case MAT_IGNORE_OFF_PROC_ENTRIES: 1215 case MAT_USE_HASH_TABLE: 1216 PetscLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 1217 break; 1218 case MAT_SYMMETRIC: 1219 case MAT_STRUCTURALLY_SYMMETRIC: 1220 case MAT_NOT_SYMMETRIC: 1221 case MAT_NOT_STRUCTURALLY_SYMMETRIC: 1222 case MAT_HERMITIAN: 1223 case MAT_NOT_HERMITIAN: 1224 case MAT_SYMMETRY_ETERNAL: 1225 case MAT_NOT_SYMMETRY_ETERNAL: 1226 break; 1227 default: 1228 SETERRQ(PETSC_ERR_SUP,"unknown option"); 1229 } 1230 PetscFunctionReturn(0); 1231 } 1232 1233 #undef __FUNCT__ 1234 #define __FUNCT__ "MatZeroEntries_SeqDense" 1235 int MatZeroEntries_SeqDense(Mat A) 1236 { 1237 Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1238 int lda=l->lda,m=A->m,j, ierr; 1239 1240 PetscFunctionBegin; 1241 if (lda>m) { 1242 for (j=0; j<A->n; j++) { 1243 ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar)); CHKERRQ(ierr); 1244 } 1245 } else { 1246 ierr = PetscMemzero(l->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1247 } 1248 PetscFunctionReturn(0); 1249 } 1250 1251 #undef __FUNCT__ 1252 #define __FUNCT__ "MatGetBlockSize_SeqDense" 1253 int MatGetBlockSize_SeqDense(Mat A,int *bs) 1254 { 1255 PetscFunctionBegin; 1256 *bs = 1; 1257 PetscFunctionReturn(0); 1258 } 1259 1260 #undef __FUNCT__ 1261 #define __FUNCT__ "MatZeroRows_SeqDense" 1262 int MatZeroRows_SeqDense(Mat A,IS is,const PetscScalar *diag) 1263 { 1264 Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1265 int n = A->n,i,j,ierr,N,*rows; 1266 PetscScalar *slot; 1267 1268 PetscFunctionBegin; 1269 ierr = ISGetLocalSize(is,&N);CHKERRQ(ierr); 1270 ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 1271 for (i=0; i<N; i++) { 1272 slot = l->v + rows[i]; 1273 for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 1274 } 1275 if (diag) { 1276 for (i=0; i<N; i++) { 1277 slot = l->v + (n+1)*rows[i]; 1278 *slot = *diag; 1279 } 1280 } 1281 ISRestoreIndices(is,&rows); 1282 PetscFunctionReturn(0); 1283 } 1284 1285 #undef __FUNCT__ 1286 #define __FUNCT__ "MatGetArray_SeqDense" 1287 int MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 1288 { 1289 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1290 1291 PetscFunctionBegin; 1292 *array = mat->v; 1293 PetscFunctionReturn(0); 1294 } 1295 1296 #undef __FUNCT__ 1297 #define __FUNCT__ "MatRestoreArray_SeqDense" 1298 int MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1299 { 1300 PetscFunctionBegin; 1301 *array = 0; /* user cannot accidently use the array later */ 1302 PetscFunctionReturn(0); 1303 } 1304 1305 #undef __FUNCT__ 1306 #define __FUNCT__ "MatGetSubMatrix_SeqDense" 1307 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 1308 { 1309 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1310 int i,j,ierr,m = A->m,*irow,*icol,nrows,ncols; 1311 PetscScalar *av,*bv,*v = mat->v; 1312 Mat newmat; 1313 1314 PetscFunctionBegin; 1315 ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 1316 ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1317 ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1318 ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 1319 1320 /* Check submatrixcall */ 1321 if (scall == MAT_REUSE_MATRIX) { 1322 int n_cols,n_rows; 1323 ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 1324 if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); 1325 newmat = *B; 1326 } else { 1327 /* Create and fill new matrix */ 1328 ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr); 1329 } 1330 1331 /* Now extract the data pointers and do the copy,column at a time */ 1332 bv = ((Mat_SeqDense*)newmat->data)->v; 1333 1334 for (i=0; i<ncols; i++) { 1335 av = v + m*icol[i]; 1336 for (j=0; j<nrows; j++) { 1337 *bv++ = av[irow[j]]; 1338 } 1339 } 1340 1341 /* Assemble the matrices so that the correct flags are set */ 1342 ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1343 ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1344 1345 /* Free work space */ 1346 ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 1347 ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1348 *B = newmat; 1349 PetscFunctionReturn(0); 1350 } 1351 1352 #undef __FUNCT__ 1353 #define __FUNCT__ "MatGetSubMatrices_SeqDense" 1354 int MatGetSubMatrices_SeqDense(Mat A,int n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1355 { 1356 int ierr,i; 1357 1358 PetscFunctionBegin; 1359 if (scall == MAT_INITIAL_MATRIX) { 1360 ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1361 } 1362 1363 for (i=0; i<n; i++) { 1364 ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1365 } 1366 PetscFunctionReturn(0); 1367 } 1368 1369 #undef __FUNCT__ 1370 #define __FUNCT__ "MatCopy_SeqDense" 1371 int MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 1372 { 1373 Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 1374 int lda1=a->lda,lda2=b->lda, m=A->m,n=A->n, j,ierr; 1375 1376 PetscFunctionBegin; 1377 /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 1378 if (A->ops->copy != B->ops->copy) { 1379 ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1380 PetscFunctionReturn(0); 1381 } 1382 if (m != B->m || n != B->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1383 if (lda1>m || lda2>m) { 1384 for (j=0; j<n; j++) { 1385 ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1386 } 1387 } else { 1388 ierr = PetscMemcpy(b->v,a->v,A->m*A->n*sizeof(PetscScalar));CHKERRQ(ierr); 1389 } 1390 PetscFunctionReturn(0); 1391 } 1392 1393 #undef __FUNCT__ 1394 #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1395 int MatSetUpPreallocation_SeqDense(Mat A) 1396 { 1397 int ierr; 1398 1399 PetscFunctionBegin; 1400 ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 1401 PetscFunctionReturn(0); 1402 } 1403 1404 /* -------------------------------------------------------------------*/ 1405 static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1406 MatGetRow_SeqDense, 1407 MatRestoreRow_SeqDense, 1408 MatMult_SeqDense, 1409 /* 4*/ MatMultAdd_SeqDense, 1410 MatMultTranspose_SeqDense, 1411 MatMultTransposeAdd_SeqDense, 1412 MatSolve_SeqDense, 1413 MatSolveAdd_SeqDense, 1414 MatSolveTranspose_SeqDense, 1415 /*10*/ MatSolveTransposeAdd_SeqDense, 1416 MatLUFactor_SeqDense, 1417 MatCholeskyFactor_SeqDense, 1418 MatRelax_SeqDense, 1419 MatTranspose_SeqDense, 1420 /*15*/ MatGetInfo_SeqDense, 1421 MatEqual_SeqDense, 1422 MatGetDiagonal_SeqDense, 1423 MatDiagonalScale_SeqDense, 1424 MatNorm_SeqDense, 1425 /*20*/ 0, 1426 0, 1427 0, 1428 MatSetOption_SeqDense, 1429 MatZeroEntries_SeqDense, 1430 /*25*/ MatZeroRows_SeqDense, 1431 MatLUFactorSymbolic_SeqDense, 1432 MatLUFactorNumeric_SeqDense, 1433 MatCholeskyFactorSymbolic_SeqDense, 1434 MatCholeskyFactorNumeric_SeqDense, 1435 /*30*/ MatSetUpPreallocation_SeqDense, 1436 0, 1437 0, 1438 MatGetArray_SeqDense, 1439 MatRestoreArray_SeqDense, 1440 /*35*/ MatDuplicate_SeqDense, 1441 0, 1442 0, 1443 0, 1444 0, 1445 /*40*/ MatAXPY_SeqDense, 1446 MatGetSubMatrices_SeqDense, 1447 0, 1448 MatGetValues_SeqDense, 1449 MatCopy_SeqDense, 1450 /*45*/ 0, 1451 MatScale_SeqDense, 1452 0, 1453 0, 1454 0, 1455 /*50*/ MatGetBlockSize_SeqDense, 1456 0, 1457 0, 1458 0, 1459 0, 1460 /*55*/ 0, 1461 0, 1462 0, 1463 0, 1464 0, 1465 /*60*/ 0, 1466 MatDestroy_SeqDense, 1467 MatView_SeqDense, 1468 MatGetPetscMaps_Petsc, 1469 0, 1470 /*65*/ 0, 1471 0, 1472 0, 1473 0, 1474 0, 1475 /*70*/ 0, 1476 0, 1477 0, 1478 0, 1479 0, 1480 /*75*/ 0, 1481 0, 1482 0, 1483 0, 1484 0, 1485 /*80*/ 0, 1486 0, 1487 0, 1488 0, 1489 /*85*/ MatLoad_SeqDense}; 1490 1491 #undef __FUNCT__ 1492 #define __FUNCT__ "MatCreateSeqDense" 1493 /*@C 1494 MatCreateSeqDense - Creates a sequential dense matrix that 1495 is stored in column major order (the usual Fortran 77 manner). Many 1496 of the matrix operations use the BLAS and LAPACK routines. 1497 1498 Collective on MPI_Comm 1499 1500 Input Parameters: 1501 + comm - MPI communicator, set to PETSC_COMM_SELF 1502 . m - number of rows 1503 . n - number of columns 1504 - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1505 to control all matrix memory allocation. 1506 1507 Output Parameter: 1508 . A - the matrix 1509 1510 Notes: 1511 The data input variable is intended primarily for Fortran programmers 1512 who wish to allocate their own matrix memory space. Most users should 1513 set data=PETSC_NULL. 1514 1515 Level: intermediate 1516 1517 .keywords: dense, matrix, LAPACK, BLAS 1518 1519 .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1520 @*/ 1521 int MatCreateSeqDense(MPI_Comm comm,int m,int n,PetscScalar *data,Mat *A) 1522 { 1523 int ierr; 1524 1525 PetscFunctionBegin; 1526 ierr = MatCreate(comm,m,n,m,n,A);CHKERRQ(ierr); 1527 ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1528 ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1529 PetscFunctionReturn(0); 1530 } 1531 1532 #undef __FUNCT__ 1533 #define __FUNCT__ "MatSeqDensePreallocation" 1534 /*@C 1535 MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1536 1537 Collective on MPI_Comm 1538 1539 Input Parameters: 1540 + A - the matrix 1541 - data - the array (or PETSC_NULL) 1542 1543 Notes: 1544 The data input variable is intended primarily for Fortran programmers 1545 who wish to allocate their own matrix memory space. Most users should 1546 set data=PETSC_NULL. 1547 1548 Level: intermediate 1549 1550 .keywords: dense, matrix, LAPACK, BLAS 1551 1552 .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1553 @*/ 1554 int MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 1555 { 1556 int ierr,(*f)(Mat,PetscScalar[]); 1557 1558 PetscFunctionBegin; 1559 ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1560 if (f) { 1561 ierr = (*f)(B,data);CHKERRQ(ierr); 1562 } 1563 PetscFunctionReturn(0); 1564 } 1565 1566 EXTERN_C_BEGIN 1567 #undef __FUNCT__ 1568 #define __FUNCT__ "MatSeqDensePreallocation_SeqDense" 1569 int MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 1570 { 1571 Mat_SeqDense *b; 1572 int ierr; 1573 1574 PetscFunctionBegin; 1575 B->preallocated = PETSC_TRUE; 1576 b = (Mat_SeqDense*)B->data; 1577 if (!data) { 1578 ierr = PetscMalloc((B->m*B->n+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 1579 ierr = PetscMemzero(b->v,B->m*B->n*sizeof(PetscScalar));CHKERRQ(ierr); 1580 b->user_alloc = PETSC_FALSE; 1581 PetscLogObjectMemory(B,B->n*B->m*sizeof(PetscScalar)); 1582 } else { /* user-allocated storage */ 1583 b->v = data; 1584 b->user_alloc = PETSC_TRUE; 1585 } 1586 PetscFunctionReturn(0); 1587 } 1588 EXTERN_C_END 1589 1590 #undef __FUNCT__ 1591 #define __FUNCT__ "MatSeqDenseSetLDA" 1592 /*@C 1593 MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 1594 1595 Input parameter: 1596 + A - the matrix 1597 - lda - the leading dimension 1598 1599 Notes: 1600 This routine is to be used in conjunction with MatSeqDenseSetPreallocation; 1601 it asserts that the preallocation has a leading dimension (the LDA parameter 1602 of Blas and Lapack fame) larger than M, the first dimension of the matrix. 1603 1604 Level: intermediate 1605 1606 .keywords: dense, matrix, LAPACK, BLAS 1607 1608 .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation() 1609 @*/ 1610 int MatSeqDenseSetLDA(Mat B,int lda) 1611 { 1612 Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1613 PetscFunctionBegin; 1614 if (lda<B->m) SETERRQ(1,"LDA must be at least matrix i dimension"); 1615 b->lda = lda; 1616 PetscFunctionReturn(0); 1617 } 1618 1619 /*MC 1620 MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 1621 1622 Options Database Keys: 1623 . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 1624 1625 Level: beginner 1626 1627 .seealso: MatCreateSeqDense 1628 M*/ 1629 1630 EXTERN_C_BEGIN 1631 #undef __FUNCT__ 1632 #define __FUNCT__ "MatCreate_SeqDense" 1633 int MatCreate_SeqDense(Mat B) 1634 { 1635 Mat_SeqDense *b; 1636 int ierr,size; 1637 1638 PetscFunctionBegin; 1639 ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1640 if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1641 1642 B->m = B->M = PetscMax(B->m,B->M); 1643 B->n = B->N = PetscMax(B->n,B->N); 1644 1645 ierr = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr); 1646 ierr = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr); 1647 ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1648 B->factor = 0; 1649 B->mapping = 0; 1650 PetscLogObjectMemory(B,sizeof(struct _p_Mat)); 1651 B->data = (void*)b; 1652 1653 ierr = PetscMapCreateMPI(B->comm,B->m,B->m,&B->rmap);CHKERRQ(ierr); 1654 ierr = PetscMapCreateMPI(B->comm,B->n,B->n,&B->cmap);CHKERRQ(ierr); 1655 1656 b->pivots = 0; 1657 b->roworiented = PETSC_TRUE; 1658 b->v = 0; 1659 b->lda = B->m; 1660 1661 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 1662 "MatSeqDenseSetPreallocation_SeqDense", 1663 MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 1664 PetscFunctionReturn(0); 1665 } 1666 EXTERN_C_END 1667