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