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 issocket,iascii,isbinary,isdraw; 997 998 PetscFunctionBegin; 999 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_SOCKET,&issocket);CHKERRQ(ierr); 1000 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1001 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 1002 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 1003 1004 if (iascii) { 1005 ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 1006 #if defined(PETSC_USE_SOCKET_VIEWER) 1007 } else if (issocket) { 1008 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1009 if (a->lda>A->rmap.n) SETERRQ(PETSC_ERR_SUP,"Case can not handle LDA"); 1010 ierr = PetscViewerSocketPutScalar(viewer,A->rmap.n,A->cmap.n,a->v);CHKERRQ(ierr); 1011 #endif 1012 } else if (isbinary) { 1013 ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1014 } else if (isdraw) { 1015 ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 1016 } else { 1017 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1018 } 1019 PetscFunctionReturn(0); 1020 } 1021 1022 #undef __FUNCT__ 1023 #define __FUNCT__ "MatDestroy_SeqDense" 1024 PetscErrorCode MatDestroy_SeqDense(Mat mat) 1025 { 1026 Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1027 PetscErrorCode ierr; 1028 1029 PetscFunctionBegin; 1030 #if defined(PETSC_USE_LOG) 1031 PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap.n,mat->cmap.n); 1032 #endif 1033 ierr = PetscFree(l->pivots);CHKERRQ(ierr); 1034 if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1035 ierr = PetscFree(l);CHKERRQ(ierr); 1036 ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 1037 ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 1038 ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 1039 ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 1040 PetscFunctionReturn(0); 1041 } 1042 1043 #undef __FUNCT__ 1044 #define __FUNCT__ "MatTranspose_SeqDense" 1045 PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout) 1046 { 1047 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1048 PetscErrorCode ierr; 1049 PetscInt k,j,m,n,M; 1050 PetscScalar *v,tmp; 1051 1052 PetscFunctionBegin; 1053 v = mat->v; m = A->rmap.n; M = mat->lda; n = A->cmap.n; 1054 if (!matout) { /* in place transpose */ 1055 if (m != n) { 1056 SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1057 } else { 1058 for (j=0; j<m; j++) { 1059 for (k=0; k<j; k++) { 1060 tmp = v[j + k*M]; 1061 v[j + k*M] = v[k + j*M]; 1062 v[k + j*M] = tmp; 1063 } 1064 } 1065 } 1066 } else { /* out-of-place transpose */ 1067 Mat tmat; 1068 Mat_SeqDense *tmatd; 1069 PetscScalar *v2; 1070 1071 ierr = MatCreate(A->comm,&tmat);CHKERRQ(ierr); 1072 ierr = MatSetSizes(tmat,A->cmap.n,A->rmap.n,A->cmap.n,A->rmap.n);CHKERRQ(ierr); 1073 ierr = MatSetType(tmat,A->type_name);CHKERRQ(ierr); 1074 ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1075 tmatd = (Mat_SeqDense*)tmat->data; 1076 v = mat->v; v2 = tmatd->v; 1077 for (j=0; j<n; j++) { 1078 for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1079 } 1080 ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1081 ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1082 *matout = tmat; 1083 } 1084 PetscFunctionReturn(0); 1085 } 1086 1087 #undef __FUNCT__ 1088 #define __FUNCT__ "MatEqual_SeqDense" 1089 PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1090 { 1091 Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1092 Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 1093 PetscInt i,j; 1094 PetscScalar *v1 = mat1->v,*v2 = mat2->v; 1095 1096 PetscFunctionBegin; 1097 if (A1->rmap.n != A2->rmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1098 if (A1->cmap.n != A2->cmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1099 for (i=0; i<A1->rmap.n; i++) { 1100 v1 = mat1->v+i; v2 = mat2->v+i; 1101 for (j=0; j<A1->cmap.n; j++) { 1102 if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1103 v1 += mat1->lda; v2 += mat2->lda; 1104 } 1105 } 1106 *flg = PETSC_TRUE; 1107 PetscFunctionReturn(0); 1108 } 1109 1110 #undef __FUNCT__ 1111 #define __FUNCT__ "MatGetDiagonal_SeqDense" 1112 PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1113 { 1114 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1115 PetscErrorCode ierr; 1116 PetscInt i,n,len; 1117 PetscScalar *x,zero = 0.0; 1118 1119 PetscFunctionBegin; 1120 ierr = VecSet(v,zero);CHKERRQ(ierr); 1121 ierr = VecGetSize(v,&n);CHKERRQ(ierr); 1122 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1123 len = PetscMin(A->rmap.n,A->cmap.n); 1124 if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 1125 for (i=0; i<len; i++) { 1126 x[i] = mat->v[i*mat->lda + i]; 1127 } 1128 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1129 PetscFunctionReturn(0); 1130 } 1131 1132 #undef __FUNCT__ 1133 #define __FUNCT__ "MatDiagonalScale_SeqDense" 1134 PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1135 { 1136 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1137 PetscScalar *l,*r,x,*v; 1138 PetscErrorCode ierr; 1139 PetscInt i,j,m = A->rmap.n,n = A->cmap.n; 1140 1141 PetscFunctionBegin; 1142 if (ll) { 1143 ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 1144 ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1145 if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1146 for (i=0; i<m; i++) { 1147 x = l[i]; 1148 v = mat->v + i; 1149 for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1150 } 1151 ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1152 ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1153 } 1154 if (rr) { 1155 ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 1156 ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1157 if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1158 for (i=0; i<n; i++) { 1159 x = r[i]; 1160 v = mat->v + i*m; 1161 for (j=0; j<m; j++) { (*v++) *= x;} 1162 } 1163 ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1164 ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1165 } 1166 PetscFunctionReturn(0); 1167 } 1168 1169 #undef __FUNCT__ 1170 #define __FUNCT__ "MatNorm_SeqDense" 1171 PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1172 { 1173 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1174 PetscScalar *v = mat->v; 1175 PetscReal sum = 0.0; 1176 PetscInt lda=mat->lda,m=A->rmap.n,i,j; 1177 PetscErrorCode ierr; 1178 1179 PetscFunctionBegin; 1180 if (type == NORM_FROBENIUS) { 1181 if (lda>m) { 1182 for (j=0; j<A->cmap.n; j++) { 1183 v = mat->v+j*lda; 1184 for (i=0; i<m; i++) { 1185 #if defined(PETSC_USE_COMPLEX) 1186 sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1187 #else 1188 sum += (*v)*(*v); v++; 1189 #endif 1190 } 1191 } 1192 } else { 1193 for (i=0; i<A->cmap.n*A->rmap.n; i++) { 1194 #if defined(PETSC_USE_COMPLEX) 1195 sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1196 #else 1197 sum += (*v)*(*v); v++; 1198 #endif 1199 } 1200 } 1201 *nrm = sqrt(sum); 1202 ierr = PetscLogFlops(2*A->cmap.n*A->rmap.n);CHKERRQ(ierr); 1203 } else if (type == NORM_1) { 1204 *nrm = 0.0; 1205 for (j=0; j<A->cmap.n; j++) { 1206 v = mat->v + j*mat->lda; 1207 sum = 0.0; 1208 for (i=0; i<A->rmap.n; i++) { 1209 sum += PetscAbsScalar(*v); v++; 1210 } 1211 if (sum > *nrm) *nrm = sum; 1212 } 1213 ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr); 1214 } else if (type == NORM_INFINITY) { 1215 *nrm = 0.0; 1216 for (j=0; j<A->rmap.n; j++) { 1217 v = mat->v + j; 1218 sum = 0.0; 1219 for (i=0; i<A->cmap.n; i++) { 1220 sum += PetscAbsScalar(*v); v += mat->lda; 1221 } 1222 if (sum > *nrm) *nrm = sum; 1223 } 1224 ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr); 1225 } else { 1226 SETERRQ(PETSC_ERR_SUP,"No two norm"); 1227 } 1228 PetscFunctionReturn(0); 1229 } 1230 1231 #undef __FUNCT__ 1232 #define __FUNCT__ "MatSetOption_SeqDense" 1233 PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op) 1234 { 1235 Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 1236 PetscErrorCode ierr; 1237 1238 PetscFunctionBegin; 1239 switch (op) { 1240 case MAT_ROW_ORIENTED: 1241 aij->roworiented = PETSC_TRUE; 1242 break; 1243 case MAT_COLUMN_ORIENTED: 1244 aij->roworiented = PETSC_FALSE; 1245 break; 1246 case MAT_ROWS_SORTED: 1247 case MAT_ROWS_UNSORTED: 1248 case MAT_COLUMNS_SORTED: 1249 case MAT_COLUMNS_UNSORTED: 1250 case MAT_NO_NEW_NONZERO_LOCATIONS: 1251 case MAT_YES_NEW_NONZERO_LOCATIONS: 1252 case MAT_NEW_NONZERO_LOCATION_ERR: 1253 case MAT_NO_NEW_DIAGONALS: 1254 case MAT_YES_NEW_DIAGONALS: 1255 case MAT_IGNORE_OFF_PROC_ENTRIES: 1256 case MAT_USE_HASH_TABLE: 1257 ierr = PetscInfo1(A,"Option %d ignored\n",op);CHKERRQ(ierr); 1258 break; 1259 case MAT_SYMMETRIC: 1260 case MAT_STRUCTURALLY_SYMMETRIC: 1261 case MAT_NOT_SYMMETRIC: 1262 case MAT_NOT_STRUCTURALLY_SYMMETRIC: 1263 case MAT_HERMITIAN: 1264 case MAT_NOT_HERMITIAN: 1265 case MAT_SYMMETRY_ETERNAL: 1266 case MAT_NOT_SYMMETRY_ETERNAL: 1267 break; 1268 default: 1269 SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op); 1270 } 1271 PetscFunctionReturn(0); 1272 } 1273 1274 #undef __FUNCT__ 1275 #define __FUNCT__ "MatZeroEntries_SeqDense" 1276 PetscErrorCode MatZeroEntries_SeqDense(Mat A) 1277 { 1278 Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1279 PetscErrorCode ierr; 1280 PetscInt lda=l->lda,m=A->rmap.n,j; 1281 1282 PetscFunctionBegin; 1283 if (lda>m) { 1284 for (j=0; j<A->cmap.n; j++) { 1285 ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1286 } 1287 } else { 1288 ierr = PetscMemzero(l->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 1289 } 1290 PetscFunctionReturn(0); 1291 } 1292 1293 #undef __FUNCT__ 1294 #define __FUNCT__ "MatZeroRows_SeqDense" 1295 PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag) 1296 { 1297 Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1298 PetscInt n = A->cmap.n,i,j; 1299 PetscScalar *slot; 1300 1301 PetscFunctionBegin; 1302 for (i=0; i<N; i++) { 1303 slot = l->v + rows[i]; 1304 for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 1305 } 1306 if (diag != 0.0) { 1307 for (i=0; i<N; i++) { 1308 slot = l->v + (n+1)*rows[i]; 1309 *slot = diag; 1310 } 1311 } 1312 PetscFunctionReturn(0); 1313 } 1314 1315 #undef __FUNCT__ 1316 #define __FUNCT__ "MatGetArray_SeqDense" 1317 PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 1318 { 1319 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1320 1321 PetscFunctionBegin; 1322 if (mat->lda != A->rmap.n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows"); 1323 *array = mat->v; 1324 PetscFunctionReturn(0); 1325 } 1326 1327 #undef __FUNCT__ 1328 #define __FUNCT__ "MatRestoreArray_SeqDense" 1329 PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1330 { 1331 PetscFunctionBegin; 1332 *array = 0; /* user cannot accidently use the array later */ 1333 PetscFunctionReturn(0); 1334 } 1335 1336 #undef __FUNCT__ 1337 #define __FUNCT__ "MatGetSubMatrix_SeqDense" 1338 static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 1339 { 1340 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1341 PetscErrorCode ierr; 1342 PetscInt i,j,*irow,*icol,nrows,ncols; 1343 PetscScalar *av,*bv,*v = mat->v; 1344 Mat newmat; 1345 1346 PetscFunctionBegin; 1347 ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 1348 ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1349 ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1350 ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 1351 1352 /* Check submatrixcall */ 1353 if (scall == MAT_REUSE_MATRIX) { 1354 PetscInt n_cols,n_rows; 1355 ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 1356 if (n_rows != nrows || n_cols != ncols) { 1357 /* resize the result result matrix to match number of requested rows/columns */ 1358 ierr = MatSetSizes(*B,nrows,nrows,nrows,nrows);CHKERRQ(ierr); 1359 } 1360 newmat = *B; 1361 } else { 1362 /* Create and fill new matrix */ 1363 ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr); 1364 ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 1365 ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 1366 ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1367 } 1368 1369 /* Now extract the data pointers and do the copy,column at a time */ 1370 bv = ((Mat_SeqDense*)newmat->data)->v; 1371 1372 for (i=0; i<ncols; i++) { 1373 av = v + mat->lda*icol[i]; 1374 for (j=0; j<nrows; j++) { 1375 *bv++ = av[irow[j]]; 1376 } 1377 } 1378 1379 /* Assemble the matrices so that the correct flags are set */ 1380 ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1381 ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1382 1383 /* Free work space */ 1384 ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 1385 ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1386 *B = newmat; 1387 PetscFunctionReturn(0); 1388 } 1389 1390 #undef __FUNCT__ 1391 #define __FUNCT__ "MatGetSubMatrices_SeqDense" 1392 PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1393 { 1394 PetscErrorCode ierr; 1395 PetscInt i; 1396 1397 PetscFunctionBegin; 1398 if (scall == MAT_INITIAL_MATRIX) { 1399 ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1400 } 1401 1402 for (i=0; i<n; i++) { 1403 ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1404 } 1405 PetscFunctionReturn(0); 1406 } 1407 1408 #undef __FUNCT__ 1409 #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1410 PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1411 { 1412 PetscFunctionBegin; 1413 PetscFunctionReturn(0); 1414 } 1415 1416 #undef __FUNCT__ 1417 #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1418 PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1419 { 1420 PetscFunctionBegin; 1421 PetscFunctionReturn(0); 1422 } 1423 1424 #undef __FUNCT__ 1425 #define __FUNCT__ "MatCopy_SeqDense" 1426 PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 1427 { 1428 Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 1429 PetscErrorCode ierr; 1430 PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap.n,n=A->cmap.n, j; 1431 1432 PetscFunctionBegin; 1433 /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 1434 if (A->ops->copy != B->ops->copy) { 1435 ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1436 PetscFunctionReturn(0); 1437 } 1438 if (m != B->rmap.n || n != B->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1439 if (lda1>m || lda2>m) { 1440 for (j=0; j<n; j++) { 1441 ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1442 } 1443 } else { 1444 ierr = PetscMemcpy(b->v,a->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 1445 } 1446 PetscFunctionReturn(0); 1447 } 1448 1449 #undef __FUNCT__ 1450 #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1451 PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A) 1452 { 1453 PetscErrorCode ierr; 1454 1455 PetscFunctionBegin; 1456 ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 1457 PetscFunctionReturn(0); 1458 } 1459 1460 #undef __FUNCT__ 1461 #define __FUNCT__ "MatSetSizes_SeqDense" 1462 PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N) 1463 { 1464 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1465 PetscErrorCode ierr; 1466 PetscFunctionBegin; 1467 /* this will not be called before lda, Mmax, and Nmax have been set */ 1468 m = PetscMax(m,M); 1469 n = PetscMax(n,N); 1470 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); 1471 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); 1472 A->rmap.n = A->rmap.n = m; 1473 A->cmap.n = A->cmap.N = n; 1474 if (a->changelda) a->lda = m; 1475 ierr = PetscMemzero(a->v,m*n*sizeof(PetscScalar));CHKERRQ(ierr); 1476 PetscFunctionReturn(0); 1477 } 1478 1479 /* ----------------------------------------------------------------*/ 1480 #undef __FUNCT__ 1481 #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1482 PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1483 { 1484 PetscErrorCode ierr; 1485 1486 PetscFunctionBegin; 1487 if (scall == MAT_INITIAL_MATRIX){ 1488 ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1489 } 1490 ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1491 PetscFunctionReturn(0); 1492 } 1493 1494 #undef __FUNCT__ 1495 #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1496 PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1497 { 1498 PetscErrorCode ierr; 1499 PetscInt m=A->rmap.n,n=B->cmap.n; 1500 Mat Cmat; 1501 1502 PetscFunctionBegin; 1503 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); 1504 ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1505 ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1506 ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1507 ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1508 Cmat->assembled = PETSC_TRUE; 1509 *C = Cmat; 1510 PetscFunctionReturn(0); 1511 } 1512 1513 #undef __FUNCT__ 1514 #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1515 PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1516 { 1517 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1518 Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1519 Mat_SeqDense *c = (Mat_SeqDense*)C->data; 1520 PetscBLASInt m=(PetscBLASInt)A->rmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->cmap.n; 1521 PetscScalar _DOne=1.0,_DZero=0.0; 1522 1523 PetscFunctionBegin; 1524 BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1525 PetscFunctionReturn(0); 1526 } 1527 1528 #undef __FUNCT__ 1529 #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense" 1530 PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1531 { 1532 PetscErrorCode ierr; 1533 1534 PetscFunctionBegin; 1535 if (scall == MAT_INITIAL_MATRIX){ 1536 ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1537 } 1538 ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1539 PetscFunctionReturn(0); 1540 } 1541 1542 #undef __FUNCT__ 1543 #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense" 1544 PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1545 { 1546 PetscErrorCode ierr; 1547 PetscInt m=A->cmap.n,n=B->cmap.n; 1548 Mat Cmat; 1549 1550 PetscFunctionBegin; 1551 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); 1552 ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1553 ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1554 ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1555 ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1556 Cmat->assembled = PETSC_TRUE; 1557 *C = Cmat; 1558 PetscFunctionReturn(0); 1559 } 1560 1561 #undef __FUNCT__ 1562 #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense" 1563 PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1564 { 1565 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1566 Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1567 Mat_SeqDense *c = (Mat_SeqDense*)C->data; 1568 PetscBLASInt m=(PetscBLASInt)A->cmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->rmap.n; 1569 PetscScalar _DOne=1.0,_DZero=0.0; 1570 1571 PetscFunctionBegin; 1572 BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1573 PetscFunctionReturn(0); 1574 } 1575 /* -------------------------------------------------------------------*/ 1576 static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1577 MatGetRow_SeqDense, 1578 MatRestoreRow_SeqDense, 1579 MatMult_SeqDense, 1580 /* 4*/ MatMultAdd_SeqDense, 1581 MatMultTranspose_SeqDense, 1582 MatMultTransposeAdd_SeqDense, 1583 MatSolve_SeqDense, 1584 MatSolveAdd_SeqDense, 1585 MatSolveTranspose_SeqDense, 1586 /*10*/ MatSolveTransposeAdd_SeqDense, 1587 MatLUFactor_SeqDense, 1588 MatCholeskyFactor_SeqDense, 1589 MatRelax_SeqDense, 1590 MatTranspose_SeqDense, 1591 /*15*/ MatGetInfo_SeqDense, 1592 MatEqual_SeqDense, 1593 MatGetDiagonal_SeqDense, 1594 MatDiagonalScale_SeqDense, 1595 MatNorm_SeqDense, 1596 /*20*/ MatAssemblyBegin_SeqDense, 1597 MatAssemblyEnd_SeqDense, 1598 0, 1599 MatSetOption_SeqDense, 1600 MatZeroEntries_SeqDense, 1601 /*25*/ MatZeroRows_SeqDense, 1602 MatLUFactorSymbolic_SeqDense, 1603 MatLUFactorNumeric_SeqDense, 1604 MatCholeskyFactorSymbolic_SeqDense, 1605 MatCholeskyFactorNumeric_SeqDense, 1606 /*30*/ MatSetUpPreallocation_SeqDense, 1607 0, 1608 0, 1609 MatGetArray_SeqDense, 1610 MatRestoreArray_SeqDense, 1611 /*35*/ MatDuplicate_SeqDense, 1612 0, 1613 0, 1614 0, 1615 0, 1616 /*40*/ MatAXPY_SeqDense, 1617 MatGetSubMatrices_SeqDense, 1618 0, 1619 MatGetValues_SeqDense, 1620 MatCopy_SeqDense, 1621 /*45*/ 0, 1622 MatScale_SeqDense, 1623 0, 1624 0, 1625 0, 1626 /*50*/ 0, 1627 0, 1628 0, 1629 0, 1630 0, 1631 /*55*/ 0, 1632 0, 1633 0, 1634 0, 1635 0, 1636 /*60*/ 0, 1637 MatDestroy_SeqDense, 1638 MatView_SeqDense, 1639 0, 1640 0, 1641 /*65*/ 0, 1642 0, 1643 0, 1644 0, 1645 0, 1646 /*70*/ 0, 1647 0, 1648 0, 1649 0, 1650 0, 1651 /*75*/ 0, 1652 0, 1653 0, 1654 0, 1655 0, 1656 /*80*/ 0, 1657 0, 1658 0, 1659 0, 1660 /*84*/ MatLoad_SeqDense, 1661 0, 1662 0, 1663 0, 1664 0, 1665 0, 1666 /*90*/ MatMatMult_SeqDense_SeqDense, 1667 MatMatMultSymbolic_SeqDense_SeqDense, 1668 MatMatMultNumeric_SeqDense_SeqDense, 1669 0, 1670 0, 1671 /*95*/ 0, 1672 MatMatMultTranspose_SeqDense_SeqDense, 1673 MatMatMultTransposeSymbolic_SeqDense_SeqDense, 1674 MatMatMultTransposeNumeric_SeqDense_SeqDense, 1675 0, 1676 /*100*/0, 1677 0, 1678 0, 1679 0, 1680 MatSetSizes_SeqDense}; 1681 1682 #undef __FUNCT__ 1683 #define __FUNCT__ "MatCreateSeqDense" 1684 /*@C 1685 MatCreateSeqDense - Creates a sequential dense matrix that 1686 is stored in column major order (the usual Fortran 77 manner). Many 1687 of the matrix operations use the BLAS and LAPACK routines. 1688 1689 Collective on MPI_Comm 1690 1691 Input Parameters: 1692 + comm - MPI communicator, set to PETSC_COMM_SELF 1693 . m - number of rows 1694 . n - number of columns 1695 - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1696 to control all matrix memory allocation. 1697 1698 Output Parameter: 1699 . A - the matrix 1700 1701 Notes: 1702 The data input variable is intended primarily for Fortran programmers 1703 who wish to allocate their own matrix memory space. Most users should 1704 set data=PETSC_NULL. 1705 1706 Level: intermediate 1707 1708 .keywords: dense, matrix, LAPACK, BLAS 1709 1710 .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1711 @*/ 1712 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 1713 { 1714 PetscErrorCode ierr; 1715 1716 PetscFunctionBegin; 1717 ierr = MatCreate(comm,A);CHKERRQ(ierr); 1718 ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1719 ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1720 ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1721 PetscFunctionReturn(0); 1722 } 1723 1724 #undef __FUNCT__ 1725 #define __FUNCT__ "MatSeqDensePreallocation" 1726 /*@C 1727 MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1728 1729 Collective on MPI_Comm 1730 1731 Input Parameters: 1732 + A - the matrix 1733 - data - the array (or PETSC_NULL) 1734 1735 Notes: 1736 The data input variable is intended primarily for Fortran programmers 1737 who wish to allocate their own matrix memory space. Most users should 1738 need not call this routine. 1739 1740 Level: intermediate 1741 1742 .keywords: dense, matrix, LAPACK, BLAS 1743 1744 .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1745 @*/ 1746 PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 1747 { 1748 PetscErrorCode ierr,(*f)(Mat,PetscScalar[]); 1749 1750 PetscFunctionBegin; 1751 ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1752 if (f) { 1753 ierr = (*f)(B,data);CHKERRQ(ierr); 1754 } 1755 PetscFunctionReturn(0); 1756 } 1757 1758 EXTERN_C_BEGIN 1759 #undef __FUNCT__ 1760 #define __FUNCT__ "MatSeqDensePreallocation_SeqDense" 1761 PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 1762 { 1763 Mat_SeqDense *b; 1764 PetscErrorCode ierr; 1765 1766 PetscFunctionBegin; 1767 B->preallocated = PETSC_TRUE; 1768 b = (Mat_SeqDense*)B->data; 1769 if (!data) { 1770 ierr = PetscMalloc((b->lda*b->Nmax+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 1771 ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 1772 b->user_alloc = PETSC_FALSE; 1773 ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 1774 } else { /* user-allocated storage */ 1775 b->v = data; 1776 b->user_alloc = PETSC_TRUE; 1777 } 1778 PetscFunctionReturn(0); 1779 } 1780 EXTERN_C_END 1781 1782 #undef __FUNCT__ 1783 #define __FUNCT__ "MatSeqDenseSetLDA" 1784 /*@C 1785 MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 1786 1787 Input parameter: 1788 + A - the matrix 1789 - lda - the leading dimension 1790 1791 Notes: 1792 This routine is to be used in conjunction with MatSeqDenseSetPreallocation; 1793 it asserts that the preallocation has a leading dimension (the LDA parameter 1794 of Blas and Lapack fame) larger than M, the first dimension of the matrix. 1795 1796 Level: intermediate 1797 1798 .keywords: dense, matrix, LAPACK, BLAS 1799 1800 .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 1801 @*/ 1802 PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda) 1803 { 1804 Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1805 1806 PetscFunctionBegin; 1807 if (lda < B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap.n); 1808 b->lda = lda; 1809 b->changelda = PETSC_FALSE; 1810 b->Mmax = PetscMax(b->Mmax,lda); 1811 PetscFunctionReturn(0); 1812 } 1813 1814 /*MC 1815 MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 1816 1817 Options Database Keys: 1818 . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 1819 1820 Level: beginner 1821 1822 .seealso: MatCreateSeqDense() 1823 1824 M*/ 1825 1826 EXTERN_C_BEGIN 1827 #undef __FUNCT__ 1828 #define __FUNCT__ "MatCreate_SeqDense" 1829 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B) 1830 { 1831 Mat_SeqDense *b; 1832 PetscErrorCode ierr; 1833 PetscMPIInt size; 1834 1835 PetscFunctionBegin; 1836 ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 1837 if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 1838 1839 B->rmap.bs = B->cmap.bs = 1; 1840 ierr = PetscMapInitialize(B->comm,&B->rmap);CHKERRQ(ierr); 1841 ierr = PetscMapInitialize(B->comm,&B->cmap);CHKERRQ(ierr); 1842 1843 ierr = PetscNew(Mat_SeqDense,&b);CHKERRQ(ierr); 1844 ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1845 B->factor = 0; 1846 B->mapping = 0; 1847 ierr = PetscLogObjectMemory(B,sizeof(struct _p_Mat));CHKERRQ(ierr); 1848 B->data = (void*)b; 1849 1850 1851 b->pivots = 0; 1852 b->roworiented = PETSC_TRUE; 1853 b->v = 0; 1854 b->lda = B->rmap.n; 1855 b->changelda = PETSC_FALSE; 1856 b->Mmax = B->rmap.n; 1857 b->Nmax = B->cmap.n; 1858 1859 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 1860 "MatSeqDenseSetPreallocation_SeqDense", 1861 MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 1862 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C", 1863 "MatMatMult_SeqAIJ_SeqDense", 1864 MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 1865 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C", 1866 "MatMatMultSymbolic_SeqAIJ_SeqDense", 1867 MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 1868 ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C", 1869 "MatMatMultNumeric_SeqAIJ_SeqDense", 1870 MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 1871 PetscFunctionReturn(0); 1872 } 1873 1874 1875 EXTERN_C_END 1876