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