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