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