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