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 #if defined(PETSC_USE_COMPLEX) 754 PetscTruth allreal = PETSC_TRUE; 755 #endif 756 757 PetscFunctionBegin; 758 ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 759 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 760 if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 761 PetscFunctionReturn(0); /* do nothing for now */ 762 } else if (format == PETSC_VIEWER_ASCII_COMMON) { 763 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 764 for (i=0; i<A->rmap.n; i++) { 765 v = a->v + i; 766 ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 767 for (j=0; j<A->cmap.n; j++) { 768 #if defined(PETSC_USE_COMPLEX) 769 if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 770 ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 771 } else if (PetscRealPart(*v)) { 772 ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr); 773 } 774 #else 775 if (*v) { 776 ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr); 777 } 778 #endif 779 v += a->lda; 780 } 781 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 782 } 783 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 784 } else { 785 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 786 #if defined(PETSC_USE_COMPLEX) 787 /* determine if matrix has all real values */ 788 v = a->v; 789 for (i=0; i<A->rmap.n*A->cmap.n; i++) { 790 if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 791 } 792 #endif 793 if (format == PETSC_VIEWER_ASCII_MATLAB) { 794 ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 795 ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap.n,A->cmap.n);CHKERRQ(ierr); 796 ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap.n,A->cmap.n);CHKERRQ(ierr); 797 ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 798 } 799 800 for (i=0; i<A->rmap.n; i++) { 801 v = a->v + i; 802 for (j=0; j<A->cmap.n; j++) { 803 #if defined(PETSC_USE_COMPLEX) 804 if (allreal) { 805 ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr); 806 } else { 807 ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 808 } 809 #else 810 ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr); 811 #endif 812 v += a->lda; 813 } 814 ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 815 } 816 if (format == PETSC_VIEWER_ASCII_MATLAB) { 817 ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 818 } 819 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 820 } 821 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 822 PetscFunctionReturn(0); 823 } 824 825 #undef __FUNCT__ 826 #define __FUNCT__ "MatView_SeqDense_Binary" 827 static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 828 { 829 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 830 PetscErrorCode ierr; 831 int fd; 832 PetscInt ict,j,n = A->cmap.n,m = A->rmap.n,i,*col_lens,nz = m*n; 833 PetscScalar *v,*anonz,*vals; 834 PetscViewerFormat format; 835 836 PetscFunctionBegin; 837 ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 838 839 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 840 if (format == PETSC_VIEWER_BINARY_NATIVE) { 841 /* store the matrix as a dense matrix */ 842 ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 843 col_lens[0] = MAT_FILE_COOKIE; 844 col_lens[1] = m; 845 col_lens[2] = n; 846 col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 847 ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 848 ierr = PetscFree(col_lens);CHKERRQ(ierr); 849 850 /* write out matrix, by rows */ 851 ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 852 v = a->v; 853 for (i=0; i<m; i++) { 854 for (j=0; j<n; j++) { 855 vals[i + j*m] = *v++; 856 } 857 } 858 ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 859 ierr = PetscFree(vals);CHKERRQ(ierr); 860 } else { 861 ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 862 col_lens[0] = MAT_FILE_COOKIE; 863 col_lens[1] = m; 864 col_lens[2] = n; 865 col_lens[3] = nz; 866 867 /* store lengths of each row and write (including header) to file */ 868 for (i=0; i<m; i++) col_lens[4+i] = n; 869 ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 870 871 /* Possibly should write in smaller increments, not whole matrix at once? */ 872 /* store column indices (zero start index) */ 873 ict = 0; 874 for (i=0; i<m; i++) { 875 for (j=0; j<n; j++) col_lens[ict++] = j; 876 } 877 ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 878 ierr = PetscFree(col_lens);CHKERRQ(ierr); 879 880 /* store nonzero values */ 881 ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 882 ict = 0; 883 for (i=0; i<m; i++) { 884 v = a->v + i; 885 for (j=0; j<n; j++) { 886 anonz[ict++] = *v; v += a->lda; 887 } 888 } 889 ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 890 ierr = PetscFree(anonz);CHKERRQ(ierr); 891 } 892 PetscFunctionReturn(0); 893 } 894 895 #undef __FUNCT__ 896 #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 897 PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 898 { 899 Mat A = (Mat) Aa; 900 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 901 PetscErrorCode ierr; 902 PetscInt m = A->rmap.n,n = A->cmap.n,color,i,j; 903 PetscScalar *v = a->v; 904 PetscViewer viewer; 905 PetscDraw popup; 906 PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 907 PetscViewerFormat format; 908 909 PetscFunctionBegin; 910 911 ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 912 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 913 ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 914 915 /* Loop over matrix elements drawing boxes */ 916 if (format != PETSC_VIEWER_DRAW_CONTOUR) { 917 /* Blue for negative and Red for positive */ 918 color = PETSC_DRAW_BLUE; 919 for(j = 0; j < n; j++) { 920 x_l = j; 921 x_r = x_l + 1.0; 922 for(i = 0; i < m; i++) { 923 y_l = m - i - 1.0; 924 y_r = y_l + 1.0; 925 #if defined(PETSC_USE_COMPLEX) 926 if (PetscRealPart(v[j*m+i]) > 0.) { 927 color = PETSC_DRAW_RED; 928 } else if (PetscRealPart(v[j*m+i]) < 0.) { 929 color = PETSC_DRAW_BLUE; 930 } else { 931 continue; 932 } 933 #else 934 if (v[j*m+i] > 0.) { 935 color = PETSC_DRAW_RED; 936 } else if (v[j*m+i] < 0.) { 937 color = PETSC_DRAW_BLUE; 938 } else { 939 continue; 940 } 941 #endif 942 ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 943 } 944 } 945 } else { 946 /* use contour shading to indicate magnitude of values */ 947 /* first determine max of all nonzero values */ 948 for(i = 0; i < m*n; i++) { 949 if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 950 } 951 scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 952 ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 953 if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 954 for(j = 0; j < n; j++) { 955 x_l = j; 956 x_r = x_l + 1.0; 957 for(i = 0; i < m; i++) { 958 y_l = m - i - 1.0; 959 y_r = y_l + 1.0; 960 color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 961 ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 962 } 963 } 964 } 965 PetscFunctionReturn(0); 966 } 967 968 #undef __FUNCT__ 969 #define __FUNCT__ "MatView_SeqDense_Draw" 970 PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 971 { 972 PetscDraw draw; 973 PetscTruth isnull; 974 PetscReal xr,yr,xl,yl,h,w; 975 PetscErrorCode ierr; 976 977 PetscFunctionBegin; 978 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 979 ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 980 if (isnull) PetscFunctionReturn(0); 981 982 ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 983 xr = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0; 984 xr += w; yr += h; xl = -w; yl = -h; 985 ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 986 ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 987 ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 988 PetscFunctionReturn(0); 989 } 990 991 #undef __FUNCT__ 992 #define __FUNCT__ "MatView_SeqDense" 993 PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 994 { 995 PetscErrorCode ierr; 996 PetscTruth iascii,isbinary,isdraw; 997 998 PetscFunctionBegin; 999 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1000 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 1001 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 1002 1003 if (iascii) { 1004 ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 1005 } else if (isbinary) { 1006 ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1007 } else if (isdraw) { 1008 ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 1009 } else { 1010 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1011 } 1012 PetscFunctionReturn(0); 1013 } 1014 1015 #undef __FUNCT__ 1016 #define __FUNCT__ "MatDestroy_SeqDense" 1017 PetscErrorCode MatDestroy_SeqDense(Mat mat) 1018 { 1019 Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1020 PetscErrorCode ierr; 1021 1022 PetscFunctionBegin; 1023 #if defined(PETSC_USE_LOG) 1024 PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap.n,mat->cmap.n); 1025 #endif 1026 ierr = PetscFree(l->pivots);CHKERRQ(ierr); 1027 if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1028 ierr = PetscFree(l);CHKERRQ(ierr); 1029 ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 1030 ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 1031 ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 1032 ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 1033 PetscFunctionReturn(0); 1034 } 1035 1036 #undef __FUNCT__ 1037 #define __FUNCT__ "MatTranspose_SeqDense" 1038 PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout) 1039 { 1040 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1041 PetscErrorCode ierr; 1042 PetscInt k,j,m,n,M; 1043 PetscScalar *v,tmp; 1044 1045 PetscFunctionBegin; 1046 v = mat->v; m = A->rmap.n; M = mat->lda; n = A->cmap.n; 1047 if (!matout) { /* in place transpose */ 1048 if (m != n) { 1049 SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1050 } else { 1051 for (j=0; j<m; j++) { 1052 for (k=0; k<j; k++) { 1053 tmp = v[j + k*M]; 1054 v[j + k*M] = v[k + j*M]; 1055 v[k + j*M] = tmp; 1056 } 1057 } 1058 } 1059 } else { /* out-of-place transpose */ 1060 Mat tmat; 1061 Mat_SeqDense *tmatd; 1062 PetscScalar *v2; 1063 1064 ierr = MatCreate(A->comm,&tmat);CHKERRQ(ierr); 1065 ierr = MatSetSizes(tmat,A->cmap.n,A->rmap.n,A->cmap.n,A->rmap.n);CHKERRQ(ierr); 1066 ierr = MatSetType(tmat,A->type_name);CHKERRQ(ierr); 1067 ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1068 tmatd = (Mat_SeqDense*)tmat->data; 1069 v = mat->v; v2 = tmatd->v; 1070 for (j=0; j<n; j++) { 1071 for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1072 } 1073 ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1074 ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1075 *matout = tmat; 1076 } 1077 PetscFunctionReturn(0); 1078 } 1079 1080 #undef __FUNCT__ 1081 #define __FUNCT__ "MatEqual_SeqDense" 1082 PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1083 { 1084 Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1085 Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 1086 PetscInt i,j; 1087 PetscScalar *v1 = mat1->v,*v2 = mat2->v; 1088 1089 PetscFunctionBegin; 1090 if (A1->rmap.n != A2->rmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1091 if (A1->cmap.n != A2->cmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1092 for (i=0; i<A1->rmap.n; i++) { 1093 v1 = mat1->v+i; v2 = mat2->v+i; 1094 for (j=0; j<A1->cmap.n; j++) { 1095 if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1096 v1 += mat1->lda; v2 += mat2->lda; 1097 } 1098 } 1099 *flg = PETSC_TRUE; 1100 PetscFunctionReturn(0); 1101 } 1102 1103 #undef __FUNCT__ 1104 #define __FUNCT__ "MatGetDiagonal_SeqDense" 1105 PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1106 { 1107 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1108 PetscErrorCode ierr; 1109 PetscInt i,n,len; 1110 PetscScalar *x,zero = 0.0; 1111 1112 PetscFunctionBegin; 1113 ierr = VecSet(v,zero);CHKERRQ(ierr); 1114 ierr = VecGetSize(v,&n);CHKERRQ(ierr); 1115 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1116 len = PetscMin(A->rmap.n,A->cmap.n); 1117 if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 1118 for (i=0; i<len; i++) { 1119 x[i] = mat->v[i*mat->lda + i]; 1120 } 1121 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1122 PetscFunctionReturn(0); 1123 } 1124 1125 #undef __FUNCT__ 1126 #define __FUNCT__ "MatDiagonalScale_SeqDense" 1127 PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1128 { 1129 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1130 PetscScalar *l,*r,x,*v; 1131 PetscErrorCode ierr; 1132 PetscInt i,j,m = A->rmap.n,n = A->cmap.n; 1133 1134 PetscFunctionBegin; 1135 if (ll) { 1136 ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 1137 ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1138 if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1139 for (i=0; i<m; i++) { 1140 x = l[i]; 1141 v = mat->v + i; 1142 for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1143 } 1144 ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1145 ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1146 } 1147 if (rr) { 1148 ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 1149 ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1150 if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1151 for (i=0; i<n; i++) { 1152 x = r[i]; 1153 v = mat->v + i*m; 1154 for (j=0; j<m; j++) { (*v++) *= x;} 1155 } 1156 ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1157 ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1158 } 1159 PetscFunctionReturn(0); 1160 } 1161 1162 #undef __FUNCT__ 1163 #define __FUNCT__ "MatNorm_SeqDense" 1164 PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1165 { 1166 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1167 PetscScalar *v = mat->v; 1168 PetscReal sum = 0.0; 1169 PetscInt lda=mat->lda,m=A->rmap.n,i,j; 1170 PetscErrorCode ierr; 1171 1172 PetscFunctionBegin; 1173 if (type == NORM_FROBENIUS) { 1174 if (lda>m) { 1175 for (j=0; j<A->cmap.n; j++) { 1176 v = mat->v+j*lda; 1177 for (i=0; i<m; i++) { 1178 #if defined(PETSC_USE_COMPLEX) 1179 sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1180 #else 1181 sum += (*v)*(*v); v++; 1182 #endif 1183 } 1184 } 1185 } else { 1186 for (i=0; i<A->cmap.n*A->rmap.n; i++) { 1187 #if defined(PETSC_USE_COMPLEX) 1188 sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1189 #else 1190 sum += (*v)*(*v); v++; 1191 #endif 1192 } 1193 } 1194 *nrm = sqrt(sum); 1195 ierr = PetscLogFlops(2*A->cmap.n*A->rmap.n);CHKERRQ(ierr); 1196 } else if (type == NORM_1) { 1197 *nrm = 0.0; 1198 for (j=0; j<A->cmap.n; j++) { 1199 v = mat->v + j*mat->lda; 1200 sum = 0.0; 1201 for (i=0; i<A->rmap.n; i++) { 1202 sum += PetscAbsScalar(*v); v++; 1203 } 1204 if (sum > *nrm) *nrm = sum; 1205 } 1206 ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr); 1207 } else if (type == NORM_INFINITY) { 1208 *nrm = 0.0; 1209 for (j=0; j<A->rmap.n; j++) { 1210 v = mat->v + j; 1211 sum = 0.0; 1212 for (i=0; i<A->cmap.n; i++) { 1213 sum += PetscAbsScalar(*v); v += mat->lda; 1214 } 1215 if (sum > *nrm) *nrm = sum; 1216 } 1217 ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr); 1218 } else { 1219 SETERRQ(PETSC_ERR_SUP,"No two norm"); 1220 } 1221 PetscFunctionReturn(0); 1222 } 1223 1224 #undef __FUNCT__ 1225 #define __FUNCT__ "MatSetOption_SeqDense" 1226 PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op) 1227 { 1228 Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 1229 PetscErrorCode ierr; 1230 1231 PetscFunctionBegin; 1232 switch (op) { 1233 case MAT_ROW_ORIENTED: 1234 aij->roworiented = PETSC_TRUE; 1235 break; 1236 case MAT_COLUMN_ORIENTED: 1237 aij->roworiented = PETSC_FALSE; 1238 break; 1239 case MAT_ROWS_SORTED: 1240 case MAT_ROWS_UNSORTED: 1241 case MAT_COLUMNS_SORTED: 1242 case MAT_COLUMNS_UNSORTED: 1243 case MAT_NO_NEW_NONZERO_LOCATIONS: 1244 case MAT_YES_NEW_NONZERO_LOCATIONS: 1245 case MAT_NEW_NONZERO_LOCATION_ERR: 1246 case MAT_NO_NEW_DIAGONALS: 1247 case MAT_YES_NEW_DIAGONALS: 1248 case MAT_IGNORE_OFF_PROC_ENTRIES: 1249 case MAT_USE_HASH_TABLE: 1250 case MAT_SYMMETRIC: 1251 case MAT_STRUCTURALLY_SYMMETRIC: 1252 case MAT_NOT_SYMMETRIC: 1253 case MAT_NOT_STRUCTURALLY_SYMMETRIC: 1254 case MAT_HERMITIAN: 1255 case MAT_NOT_HERMITIAN: 1256 case MAT_SYMMETRY_ETERNAL: 1257 case MAT_NOT_SYMMETRY_ETERNAL: 1258 ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 1259 break; 1260 default: 1261 SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op); 1262 } 1263 PetscFunctionReturn(0); 1264 } 1265 1266 #undef __FUNCT__ 1267 #define __FUNCT__ "MatZeroEntries_SeqDense" 1268 PetscErrorCode MatZeroEntries_SeqDense(Mat A) 1269 { 1270 Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1271 PetscErrorCode ierr; 1272 PetscInt lda=l->lda,m=A->rmap.n,j; 1273 1274 PetscFunctionBegin; 1275 if (lda>m) { 1276 for (j=0; j<A->cmap.n; j++) { 1277 ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1278 } 1279 } else { 1280 ierr = PetscMemzero(l->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 1281 } 1282 PetscFunctionReturn(0); 1283 } 1284 1285 #undef __FUNCT__ 1286 #define __FUNCT__ "MatZeroRows_SeqDense" 1287 PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag) 1288 { 1289 Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1290 PetscInt n = A->cmap.n,i,j; 1291 PetscScalar *slot; 1292 1293 PetscFunctionBegin; 1294 for (i=0; i<N; i++) { 1295 slot = l->v + rows[i]; 1296 for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 1297 } 1298 if (diag != 0.0) { 1299 for (i=0; i<N; i++) { 1300 slot = l->v + (n+1)*rows[i]; 1301 *slot = diag; 1302 } 1303 } 1304 PetscFunctionReturn(0); 1305 } 1306 1307 #undef __FUNCT__ 1308 #define __FUNCT__ "MatGetArray_SeqDense" 1309 PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 1310 { 1311 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1312 1313 PetscFunctionBegin; 1314 if (mat->lda != A->rmap.n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows"); 1315 *array = mat->v; 1316 PetscFunctionReturn(0); 1317 } 1318 1319 #undef __FUNCT__ 1320 #define __FUNCT__ "MatRestoreArray_SeqDense" 1321 PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1322 { 1323 PetscFunctionBegin; 1324 *array = 0; /* user cannot accidently use the array later */ 1325 PetscFunctionReturn(0); 1326 } 1327 1328 #undef __FUNCT__ 1329 #define __FUNCT__ "MatGetSubMatrix_SeqDense" 1330 static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 1331 { 1332 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1333 PetscErrorCode ierr; 1334 PetscInt i,j,*irow,*icol,nrows,ncols; 1335 PetscScalar *av,*bv,*v = mat->v; 1336 Mat newmat; 1337 1338 PetscFunctionBegin; 1339 ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 1340 ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1341 ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1342 ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 1343 1344 /* Check submatrixcall */ 1345 if (scall == MAT_REUSE_MATRIX) { 1346 PetscInt n_cols,n_rows; 1347 ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 1348 if (n_rows != nrows || n_cols != ncols) { 1349 /* resize the result result matrix to match number of requested rows/columns */ 1350 ierr = MatSetSizes(*B,nrows,nrows,nrows,nrows);CHKERRQ(ierr); 1351 } 1352 newmat = *B; 1353 } else { 1354 /* Create and fill new matrix */ 1355 ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr); 1356 ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 1357 ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 1358 ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1359 } 1360 1361 /* Now extract the data pointers and do the copy,column at a time */ 1362 bv = ((Mat_SeqDense*)newmat->data)->v; 1363 1364 for (i=0; i<ncols; i++) { 1365 av = v + mat->lda*icol[i]; 1366 for (j=0; j<nrows; j++) { 1367 *bv++ = av[irow[j]]; 1368 } 1369 } 1370 1371 /* Assemble the matrices so that the correct flags are set */ 1372 ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1373 ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1374 1375 /* Free work space */ 1376 ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 1377 ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1378 *B = newmat; 1379 PetscFunctionReturn(0); 1380 } 1381 1382 #undef __FUNCT__ 1383 #define __FUNCT__ "MatGetSubMatrices_SeqDense" 1384 PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1385 { 1386 PetscErrorCode ierr; 1387 PetscInt i; 1388 1389 PetscFunctionBegin; 1390 if (scall == MAT_INITIAL_MATRIX) { 1391 ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1392 } 1393 1394 for (i=0; i<n; i++) { 1395 ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1396 } 1397 PetscFunctionReturn(0); 1398 } 1399 1400 #undef __FUNCT__ 1401 #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1402 PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1403 { 1404 PetscFunctionBegin; 1405 PetscFunctionReturn(0); 1406 } 1407 1408 #undef __FUNCT__ 1409 #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1410 PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1411 { 1412 PetscFunctionBegin; 1413 PetscFunctionReturn(0); 1414 } 1415 1416 #undef __FUNCT__ 1417 #define __FUNCT__ "MatCopy_SeqDense" 1418 PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 1419 { 1420 Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 1421 PetscErrorCode ierr; 1422 PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap.n,n=A->cmap.n, j; 1423 1424 PetscFunctionBegin; 1425 /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 1426 if (A->ops->copy != B->ops->copy) { 1427 ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1428 PetscFunctionReturn(0); 1429 } 1430 if (m != B->rmap.n || n != B->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1431 if (lda1>m || lda2>m) { 1432 for (j=0; j<n; j++) { 1433 ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1434 } 1435 } else { 1436 ierr = PetscMemcpy(b->v,a->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 1437 } 1438 PetscFunctionReturn(0); 1439 } 1440 1441 #undef __FUNCT__ 1442 #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1443 PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A) 1444 { 1445 PetscErrorCode ierr; 1446 1447 PetscFunctionBegin; 1448 ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 1449 PetscFunctionReturn(0); 1450 } 1451 1452 #undef __FUNCT__ 1453 #define __FUNCT__ "MatSetSizes_SeqDense" 1454 PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N) 1455 { 1456 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1457 PetscErrorCode ierr; 1458 PetscFunctionBegin; 1459 /* this will not be called before lda, Mmax, and Nmax have been set */ 1460 m = PetscMax(m,M); 1461 n = PetscMax(n,N); 1462 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); 1463 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); 1464 A->rmap.n = A->rmap.n = m; 1465 A->cmap.n = A->cmap.N = n; 1466 if (a->changelda) a->lda = m; 1467 ierr = PetscMemzero(a->v,m*n*sizeof(PetscScalar));CHKERRQ(ierr); 1468 PetscFunctionReturn(0); 1469 } 1470 1471 /* ----------------------------------------------------------------*/ 1472 #undef __FUNCT__ 1473 #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1474 PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1475 { 1476 PetscErrorCode ierr; 1477 1478 PetscFunctionBegin; 1479 if (scall == MAT_INITIAL_MATRIX){ 1480 ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1481 } 1482 ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1483 PetscFunctionReturn(0); 1484 } 1485 1486 #undef __FUNCT__ 1487 #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1488 PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1489 { 1490 PetscErrorCode ierr; 1491 PetscInt m=A->rmap.n,n=B->cmap.n; 1492 Mat Cmat; 1493 1494 PetscFunctionBegin; 1495 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); 1496 ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1497 ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1498 ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1499 ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1500 Cmat->assembled = PETSC_TRUE; 1501 *C = Cmat; 1502 PetscFunctionReturn(0); 1503 } 1504 1505 #undef __FUNCT__ 1506 #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1507 PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1508 { 1509 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1510 Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1511 Mat_SeqDense *c = (Mat_SeqDense*)C->data; 1512 PetscBLASInt m=(PetscBLASInt)A->rmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->cmap.n; 1513 PetscScalar _DOne=1.0,_DZero=0.0; 1514 1515 PetscFunctionBegin; 1516 BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1517 PetscFunctionReturn(0); 1518 } 1519 1520 #undef __FUNCT__ 1521 #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense" 1522 PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1523 { 1524 PetscErrorCode ierr; 1525 1526 PetscFunctionBegin; 1527 if (scall == MAT_INITIAL_MATRIX){ 1528 ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1529 } 1530 ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1531 PetscFunctionReturn(0); 1532 } 1533 1534 #undef __FUNCT__ 1535 #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense" 1536 PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1537 { 1538 PetscErrorCode ierr; 1539 PetscInt m=A->cmap.n,n=B->cmap.n; 1540 Mat Cmat; 1541 1542 PetscFunctionBegin; 1543 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); 1544 ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1545 ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1546 ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1547 ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1548 Cmat->assembled = PETSC_TRUE; 1549 *C = Cmat; 1550 PetscFunctionReturn(0); 1551 } 1552 1553 #undef __FUNCT__ 1554 #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense" 1555 PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1556 { 1557 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1558 Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1559 Mat_SeqDense *c = (Mat_SeqDense*)C->data; 1560 PetscBLASInt m=(PetscBLASInt)A->cmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->rmap.n; 1561 PetscScalar _DOne=1.0,_DZero=0.0; 1562 1563 PetscFunctionBegin; 1564 BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1565 PetscFunctionReturn(0); 1566 } 1567 /* -------------------------------------------------------------------*/ 1568 static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1569 MatGetRow_SeqDense, 1570 MatRestoreRow_SeqDense, 1571 MatMult_SeqDense, 1572 /* 4*/ MatMultAdd_SeqDense, 1573 MatMultTranspose_SeqDense, 1574 MatMultTransposeAdd_SeqDense, 1575 MatSolve_SeqDense, 1576 MatSolveAdd_SeqDense, 1577 MatSolveTranspose_SeqDense, 1578 /*10*/ MatSolveTransposeAdd_SeqDense, 1579 MatLUFactor_SeqDense, 1580 MatCholeskyFactor_SeqDense, 1581 MatRelax_SeqDense, 1582 MatTranspose_SeqDense, 1583 /*15*/ MatGetInfo_SeqDense, 1584 MatEqual_SeqDense, 1585 MatGetDiagonal_SeqDense, 1586 MatDiagonalScale_SeqDense, 1587 MatNorm_SeqDense, 1588 /*20*/ MatAssemblyBegin_SeqDense, 1589 MatAssemblyEnd_SeqDense, 1590 0, 1591 MatSetOption_SeqDense, 1592 MatZeroEntries_SeqDense, 1593 /*25*/ MatZeroRows_SeqDense, 1594 MatLUFactorSymbolic_SeqDense, 1595 MatLUFactorNumeric_SeqDense, 1596 MatCholeskyFactorSymbolic_SeqDense, 1597 MatCholeskyFactorNumeric_SeqDense, 1598 /*30*/ MatSetUpPreallocation_SeqDense, 1599 0, 1600 0, 1601 MatGetArray_SeqDense, 1602 MatRestoreArray_SeqDense, 1603 /*35*/ MatDuplicate_SeqDense, 1604 0, 1605 0, 1606 0, 1607 0, 1608 /*40*/ MatAXPY_SeqDense, 1609 MatGetSubMatrices_SeqDense, 1610 0, 1611 MatGetValues_SeqDense, 1612 MatCopy_SeqDense, 1613 /*45*/ 0, 1614 MatScale_SeqDense, 1615 0, 1616 0, 1617 0, 1618 /*50*/ 0, 1619 0, 1620 0, 1621 0, 1622 0, 1623 /*55*/ 0, 1624 0, 1625 0, 1626 0, 1627 0, 1628 /*60*/ 0, 1629 MatDestroy_SeqDense, 1630 MatView_SeqDense, 1631 0, 1632 0, 1633 /*65*/ 0, 1634 0, 1635 0, 1636 0, 1637 0, 1638 /*70*/ 0, 1639 0, 1640 0, 1641 0, 1642 0, 1643 /*75*/ 0, 1644 0, 1645 0, 1646 0, 1647 0, 1648 /*80*/ 0, 1649 0, 1650 0, 1651 0, 1652 /*84*/ MatLoad_SeqDense, 1653 0, 1654 0, 1655 0, 1656 0, 1657 0, 1658 /*90*/ MatMatMult_SeqDense_SeqDense, 1659 MatMatMultSymbolic_SeqDense_SeqDense, 1660 MatMatMultNumeric_SeqDense_SeqDense, 1661 0, 1662 0, 1663 /*95*/ 0, 1664 MatMatMultTranspose_SeqDense_SeqDense, 1665 MatMatMultTransposeSymbolic_SeqDense_SeqDense, 1666 MatMatMultTransposeNumeric_SeqDense_SeqDense, 1667 0, 1668 /*100*/0, 1669 0, 1670 0, 1671 0, 1672 MatSetSizes_SeqDense}; 1673 1674 #undef __FUNCT__ 1675 #define __FUNCT__ "MatCreateSeqDense" 1676 /*@C 1677 MatCreateSeqDense - Creates a sequential dense matrix that 1678 is stored in column major order (the usual Fortran 77 manner). Many 1679 of the matrix operations use the BLAS and LAPACK routines. 1680 1681 Collective on MPI_Comm 1682 1683 Input Parameters: 1684 + comm - MPI communicator, set to PETSC_COMM_SELF 1685 . m - number of rows 1686 . n - number of columns 1687 - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1688 to control all matrix memory allocation. 1689 1690 Output Parameter: 1691 . A - the matrix 1692 1693 Notes: 1694 The data input variable is intended primarily for Fortran programmers 1695 who wish to allocate their own matrix memory space. Most users should 1696 set data=PETSC_NULL. 1697 1698 Level: intermediate 1699 1700 .keywords: dense, matrix, LAPACK, BLAS 1701 1702 .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1703 @*/ 1704 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 1705 { 1706 PetscErrorCode ierr; 1707 1708 PetscFunctionBegin; 1709 ierr = MatCreate(comm,A);CHKERRQ(ierr); 1710 ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1711 ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1712 ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1713 PetscFunctionReturn(0); 1714 } 1715 1716 #undef __FUNCT__ 1717 #define __FUNCT__ "MatSeqDensePreallocation" 1718 /*@C 1719 MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1720 1721 Collective on MPI_Comm 1722 1723 Input Parameters: 1724 + A - the matrix 1725 - data - the array (or PETSC_NULL) 1726 1727 Notes: 1728 The data input variable is intended primarily for Fortran programmers 1729 who wish to allocate their own matrix memory space. Most users should 1730 need not call this routine. 1731 1732 Level: intermediate 1733 1734 .keywords: dense, matrix, LAPACK, BLAS 1735 1736 .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1737 @*/ 1738 PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 1739 { 1740 PetscErrorCode ierr,(*f)(Mat,PetscScalar[]); 1741 1742 PetscFunctionBegin; 1743 ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1744 if (f) { 1745 ierr = (*f)(B,data);CHKERRQ(ierr); 1746 } 1747 PetscFunctionReturn(0); 1748 } 1749 1750 EXTERN_C_BEGIN 1751 #undef __FUNCT__ 1752 #define __FUNCT__ "MatSeqDensePreallocation_SeqDense" 1753 PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 1754 { 1755 Mat_SeqDense *b; 1756 PetscErrorCode ierr; 1757 1758 PetscFunctionBegin; 1759 B->preallocated = PETSC_TRUE; 1760 b = (Mat_SeqDense*)B->data; 1761 if (!data) { 1762 ierr = PetscMalloc((b->lda*b->Nmax+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 1763 ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 1764 b->user_alloc = PETSC_FALSE; 1765 ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 1766 } else { /* user-allocated storage */ 1767 b->v = data; 1768 b->user_alloc = PETSC_TRUE; 1769 } 1770 PetscFunctionReturn(0); 1771 } 1772 EXTERN_C_END 1773 1774 #undef __FUNCT__ 1775 #define __FUNCT__ "MatSeqDenseSetLDA" 1776 /*@C 1777 MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 1778 1779 Input parameter: 1780 + A - the matrix 1781 - lda - the leading dimension 1782 1783 Notes: 1784 This routine is to be used in conjunction with MatSeqDenseSetPreallocation; 1785 it asserts that the preallocation has a leading dimension (the LDA parameter 1786 of Blas and Lapack fame) larger than M, the first dimension of the matrix. 1787 1788 Level: intermediate 1789 1790 .keywords: dense, matrix, LAPACK, BLAS 1791 1792 .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 1793 @*/ 1794 PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda) 1795 { 1796 Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1797 1798 PetscFunctionBegin; 1799 if (lda < B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap.n); 1800 b->lda = lda; 1801 b->changelda = PETSC_FALSE; 1802 b->Mmax = PetscMax(b->Mmax,lda); 1803 PetscFunctionReturn(0); 1804 } 1805 1806 /*MC 1807 MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 1808 1809 Options Database Keys: 1810 . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 1811 1812 Level: beginner 1813 1814 .seealso: MatCreateSeqDense() 1815 1816 M*/ 1817 1818 EXTERN_C_BEGIN 1819 #undef __FUNCT__ 1820 #define __FUNCT__ "MatCreate_SeqDense" 1821 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B) 1822 { 1823 Mat_SeqDense *b; 1824 PetscErrorCode ierr; 1825 PetscMPIInt size; 1826 1827 PetscFunctionBegin; 1828 ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1829 if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1830 1831 B->rmap.bs = B->cmap.bs = 1; 1832 ierr = PetscMapInitialize(B->comm,&B->rmap);CHKERRQ(ierr); 1833 ierr = PetscMapInitialize(B->comm,&B->cmap);CHKERRQ(ierr); 1834 1835 ierr = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr); 1836 ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1837 B->factor = 0; 1838 B->mapping = 0; 1839 ierr = PetscLogObjectMemory(B,sizeof(struct _p_Mat));CHKERRQ(ierr); 1840 B->data = (void*)b; 1841 1842 1843 b->pivots = 0; 1844 b->roworiented = PETSC_TRUE; 1845 b->v = 0; 1846 b->lda = B->rmap.n; 1847 b->changelda = PETSC_FALSE; 1848 b->Mmax = B->rmap.n; 1849 b->Nmax = B->cmap.n; 1850 1851 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 1852 "MatSeqDenseSetPreallocation_SeqDense", 1853 MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 1854 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C", 1855 "MatMatMult_SeqAIJ_SeqDense", 1856 MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 1857 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C", 1858 "MatMatMultSymbolic_SeqAIJ_SeqDense", 1859 MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 1860 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C", 1861 "MatMatMultNumeric_SeqAIJ_SeqDense", 1862 MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 1863 ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr); 1864 PetscFunctionReturn(0); 1865 } 1866 1867 1868 EXTERN_C_END 1869