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