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