1 2 /* 3 Defines the basic matrix operations for sequential dense. 4 */ 5 6 #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 7 #include <petscblaslapack.h> 8 9 #include <../src/mat/impls/aij/seq/aij.h> 10 11 PetscErrorCode MatSeqDenseSymmetrize_Private(Mat A, PetscBool hermitian) 12 { 13 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 14 PetscInt j, k, n = A->rmap->n; 15 PetscScalar *v; 16 17 PetscFunctionBegin; 18 PetscCheck(A->rmap->n == A->cmap->n,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Cannot symmetrize a rectangular matrix"); 19 PetscCall(MatDenseGetArray(A,&v)); 20 if (!hermitian) { 21 for (k=0;k<n;k++) { 22 for (j=k;j<n;j++) { 23 v[j*mat->lda + k] = v[k*mat->lda + j]; 24 } 25 } 26 } else { 27 for (k=0;k<n;k++) { 28 for (j=k;j<n;j++) { 29 v[j*mat->lda + k] = PetscConj(v[k*mat->lda + j]); 30 } 31 } 32 } 33 PetscCall(MatDenseRestoreArray(A,&v)); 34 PetscFunctionReturn(0); 35 } 36 37 PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat A) 38 { 39 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 40 PetscBLASInt info,n; 41 42 PetscFunctionBegin; 43 if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 44 PetscCall(PetscBLASIntCast(A->cmap->n,&n)); 45 if (A->factortype == MAT_FACTOR_LU) { 46 PetscCheck(mat->pivots,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present"); 47 if (!mat->fwork) { 48 mat->lfwork = n; 49 PetscCall(PetscMalloc1(mat->lfwork,&mat->fwork)); 50 PetscCall(PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt))); 51 } 52 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 53 PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info)); 54 PetscCall(PetscFPTrapPop()); 55 PetscCall(PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0)); 56 } else if (A->factortype == MAT_FACTOR_CHOLESKY) { 57 if (A->spd) { 58 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 59 PetscStackCallBLAS("LAPACKpotri",LAPACKpotri_("L",&n,mat->v,&mat->lda,&info)); 60 PetscCall(PetscFPTrapPop()); 61 PetscCall(MatSeqDenseSymmetrize_Private(A,PETSC_TRUE)); 62 #if defined(PETSC_USE_COMPLEX) 63 } else if (A->hermitian) { 64 PetscCheck(mat->pivots,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present"); 65 PetscCheck(mat->fwork,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Fwork not present"); 66 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 67 PetscStackCallBLAS("LAPACKhetri",LAPACKhetri_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&info)); 68 PetscCall(PetscFPTrapPop()); 69 PetscCall(MatSeqDenseSymmetrize_Private(A,PETSC_TRUE)); 70 #endif 71 } else { /* symmetric case */ 72 PetscCheck(mat->pivots,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Pivots not present"); 73 PetscCheck(mat->fwork,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Fwork not present"); 74 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 75 PetscStackCallBLAS("LAPACKsytri",LAPACKsytri_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&info)); 76 PetscCall(PetscFPTrapPop()); 77 PetscCall(MatSeqDenseSymmetrize_Private(A,PETSC_FALSE)); 78 } 79 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad Inversion: zero pivot in row %" PetscInt_FMT,(PetscInt)info-1); 80 PetscCall(PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0)); 81 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 82 83 A->ops->solve = NULL; 84 A->ops->matsolve = NULL; 85 A->ops->solvetranspose = NULL; 86 A->ops->matsolvetranspose = NULL; 87 A->ops->solveadd = NULL; 88 A->ops->solvetransposeadd = NULL; 89 A->factortype = MAT_FACTOR_NONE; 90 PetscCall(PetscFree(A->solvertype)); 91 PetscFunctionReturn(0); 92 } 93 94 PetscErrorCode MatZeroRowsColumns_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 95 { 96 Mat_SeqDense *l = (Mat_SeqDense*)A->data; 97 PetscInt m = l->lda, n = A->cmap->n,r = A->rmap->n, i,j; 98 PetscScalar *slot,*bb,*v; 99 const PetscScalar *xx; 100 101 PetscFunctionBegin; 102 if (PetscDefined(USE_DEBUG)) { 103 for (i=0; i<N; i++) { 104 PetscCheck(rows[i] >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed"); 105 PetscCheck(rows[i] < A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %" PetscInt_FMT " requested to be zeroed greater than or equal number of rows %" PetscInt_FMT,rows[i],A->rmap->n); 106 PetscCheck(rows[i] < A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Col %" PetscInt_FMT " requested to be zeroed greater than or equal number of cols %" PetscInt_FMT,rows[i],A->cmap->n); 107 } 108 } 109 if (!N) PetscFunctionReturn(0); 110 111 /* fix right hand side if needed */ 112 if (x && b) { 113 Vec xt; 114 115 PetscCheck(A->rmap->n == A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices"); 116 PetscCall(VecDuplicate(x,&xt)); 117 PetscCall(VecCopy(x,xt)); 118 PetscCall(VecScale(xt,-1.0)); 119 PetscCall(MatMultAdd(A,xt,b,b)); 120 PetscCall(VecDestroy(&xt)); 121 PetscCall(VecGetArrayRead(x,&xx)); 122 PetscCall(VecGetArray(b,&bb)); 123 for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]]; 124 PetscCall(VecRestoreArrayRead(x,&xx)); 125 PetscCall(VecRestoreArray(b,&bb)); 126 } 127 128 PetscCall(MatDenseGetArray(A,&v)); 129 for (i=0; i<N; i++) { 130 slot = v + rows[i]*m; 131 PetscCall(PetscArrayzero(slot,r)); 132 } 133 for (i=0; i<N; i++) { 134 slot = v + rows[i]; 135 for (j=0; j<n; j++) { *slot = 0.0; slot += m;} 136 } 137 if (diag != 0.0) { 138 PetscCheck(A->rmap->n == A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices"); 139 for (i=0; i<N; i++) { 140 slot = v + (m+1)*rows[i]; 141 *slot = diag; 142 } 143 } 144 PetscCall(MatDenseRestoreArray(A,&v)); 145 PetscFunctionReturn(0); 146 } 147 148 PetscErrorCode MatPtAPNumeric_SeqDense_SeqDense(Mat A,Mat P,Mat C) 149 { 150 Mat_SeqDense *c = (Mat_SeqDense*)(C->data); 151 152 PetscFunctionBegin; 153 if (c->ptapwork) { 154 PetscCall((*C->ops->matmultnumeric)(A,P,c->ptapwork)); 155 PetscCall((*C->ops->transposematmultnumeric)(P,c->ptapwork,C)); 156 } else SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Must call MatPtAPSymbolic_SeqDense_SeqDense() first"); 157 PetscFunctionReturn(0); 158 } 159 160 PetscErrorCode MatPtAPSymbolic_SeqDense_SeqDense(Mat A,Mat P,PetscReal fill,Mat C) 161 { 162 Mat_SeqDense *c; 163 PetscBool cisdense; 164 165 PetscFunctionBegin; 166 PetscCall(MatSetSizes(C,P->cmap->n,P->cmap->n,P->cmap->N,P->cmap->N)); 167 PetscCall(PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATSEQDENSE,MATSEQDENSECUDA,"")); 168 if (!cisdense) { 169 PetscBool flg; 170 171 PetscCall(PetscObjectTypeCompare((PetscObject)P,((PetscObject)A)->type_name,&flg)); 172 PetscCall(MatSetType(C,flg ? ((PetscObject)A)->type_name : MATDENSE)); 173 } 174 PetscCall(MatSetUp(C)); 175 c = (Mat_SeqDense*)C->data; 176 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&c->ptapwork)); 177 PetscCall(MatSetSizes(c->ptapwork,A->rmap->n,P->cmap->n,A->rmap->N,P->cmap->N)); 178 PetscCall(MatSetType(c->ptapwork,((PetscObject)C)->type_name)); 179 PetscCall(MatSetUp(c->ptapwork)); 180 PetscFunctionReturn(0); 181 } 182 183 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqDense(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 184 { 185 Mat B = NULL; 186 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 187 Mat_SeqDense *b; 188 PetscInt *ai=a->i,*aj=a->j,m=A->rmap->N,n=A->cmap->N,i; 189 const MatScalar *av; 190 PetscBool isseqdense; 191 192 PetscFunctionBegin; 193 if (reuse == MAT_REUSE_MATRIX) { 194 PetscCall(PetscObjectTypeCompare((PetscObject)*newmat,MATSEQDENSE,&isseqdense)); 195 PetscCheck(isseqdense,PetscObjectComm((PetscObject)*newmat),PETSC_ERR_USER,"Cannot reuse matrix of type %s",((PetscObject)(*newmat))->type_name); 196 } 197 if (reuse != MAT_REUSE_MATRIX) { 198 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B)); 199 PetscCall(MatSetSizes(B,m,n,m,n)); 200 PetscCall(MatSetType(B,MATSEQDENSE)); 201 PetscCall(MatSeqDenseSetPreallocation(B,NULL)); 202 b = (Mat_SeqDense*)(B->data); 203 } else { 204 b = (Mat_SeqDense*)((*newmat)->data); 205 PetscCall(PetscArrayzero(b->v,m*n)); 206 } 207 PetscCall(MatSeqAIJGetArrayRead(A,&av)); 208 for (i=0; i<m; i++) { 209 PetscInt j; 210 for (j=0;j<ai[1]-ai[0];j++) { 211 b->v[*aj*m+i] = *av; 212 aj++; 213 av++; 214 } 215 ai++; 216 } 217 PetscCall(MatSeqAIJRestoreArrayRead(A,&av)); 218 219 if (reuse == MAT_INPLACE_MATRIX) { 220 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 221 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 222 PetscCall(MatHeaderReplace(A,&B)); 223 } else { 224 if (B) *newmat = B; 225 PetscCall(MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY)); 226 PetscCall(MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY)); 227 } 228 PetscFunctionReturn(0); 229 } 230 231 PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 232 { 233 Mat B = NULL; 234 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 235 PetscInt i, j; 236 PetscInt *rows, *nnz; 237 MatScalar *aa = a->v, *vals; 238 239 PetscFunctionBegin; 240 PetscCall(PetscCalloc3(A->rmap->n,&rows,A->rmap->n,&nnz,A->rmap->n,&vals)); 241 if (reuse != MAT_REUSE_MATRIX) { 242 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&B)); 243 PetscCall(MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N)); 244 PetscCall(MatSetType(B,MATSEQAIJ)); 245 for (j=0; j<A->cmap->n; j++) { 246 for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || (i == j && A->cmap->n == A->rmap->n)) ++nnz[i]; 247 aa += a->lda; 248 } 249 PetscCall(MatSeqAIJSetPreallocation(B,PETSC_DETERMINE,nnz)); 250 } else B = *newmat; 251 aa = a->v; 252 for (j=0; j<A->cmap->n; j++) { 253 PetscInt numRows = 0; 254 for (i=0; i<A->rmap->n; i++) if (aa[i] != 0.0 || (i == j && A->cmap->n == A->rmap->n)) {rows[numRows] = i; vals[numRows++] = aa[i];} 255 PetscCall(MatSetValues(B,numRows,rows,1,&j,vals,INSERT_VALUES)); 256 aa += a->lda; 257 } 258 PetscCall(PetscFree3(rows,nnz,vals)); 259 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 260 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 261 262 if (reuse == MAT_INPLACE_MATRIX) { 263 PetscCall(MatHeaderReplace(A,&B)); 264 } else if (reuse != MAT_REUSE_MATRIX) *newmat = B; 265 PetscFunctionReturn(0); 266 } 267 268 PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str) 269 { 270 Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data; 271 const PetscScalar *xv; 272 PetscScalar *yv; 273 PetscBLASInt N,m,ldax = 0,lday = 0,one = 1; 274 275 PetscFunctionBegin; 276 PetscCall(MatDenseGetArrayRead(X,&xv)); 277 PetscCall(MatDenseGetArray(Y,&yv)); 278 PetscCall(PetscBLASIntCast(X->rmap->n*X->cmap->n,&N)); 279 PetscCall(PetscBLASIntCast(X->rmap->n,&m)); 280 PetscCall(PetscBLASIntCast(x->lda,&ldax)); 281 PetscCall(PetscBLASIntCast(y->lda,&lday)); 282 if (ldax>m || lday>m) { 283 PetscInt j; 284 285 for (j=0; j<X->cmap->n; j++) { 286 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&m,&alpha,xv+j*ldax,&one,yv+j*lday,&one)); 287 } 288 } else { 289 PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&N,&alpha,xv,&one,yv,&one)); 290 } 291 PetscCall(MatDenseRestoreArrayRead(X,&xv)); 292 PetscCall(MatDenseRestoreArray(Y,&yv)); 293 PetscCall(PetscLogFlops(PetscMax(2.0*N-1,0))); 294 PetscFunctionReturn(0); 295 } 296 297 static PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 298 { 299 PetscLogDouble N = A->rmap->n*A->cmap->n; 300 301 PetscFunctionBegin; 302 info->block_size = 1.0; 303 info->nz_allocated = N; 304 info->nz_used = N; 305 info->nz_unneeded = 0; 306 info->assemblies = A->num_ass; 307 info->mallocs = 0; 308 info->memory = ((PetscObject)A)->mem; 309 info->fill_ratio_given = 0; 310 info->fill_ratio_needed = 0; 311 info->factor_mallocs = 0; 312 PetscFunctionReturn(0); 313 } 314 315 PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha) 316 { 317 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 318 PetscScalar *v; 319 PetscBLASInt one = 1,j,nz,lda = 0; 320 321 PetscFunctionBegin; 322 PetscCall(MatDenseGetArray(A,&v)); 323 PetscCall(PetscBLASIntCast(a->lda,&lda)); 324 if (lda>A->rmap->n) { 325 PetscCall(PetscBLASIntCast(A->rmap->n,&nz)); 326 for (j=0; j<A->cmap->n; j++) { 327 PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&alpha,v+j*lda,&one)); 328 } 329 } else { 330 PetscCall(PetscBLASIntCast(A->rmap->n*A->cmap->n,&nz)); 331 PetscStackCallBLAS("BLASscal",BLASscal_(&nz,&alpha,v,&one)); 332 } 333 PetscCall(PetscLogFlops(A->rmap->n*A->cmap->n)); 334 PetscCall(MatDenseRestoreArray(A,&v)); 335 PetscFunctionReturn(0); 336 } 337 338 PetscErrorCode MatShift_SeqDense(Mat A,PetscScalar alpha) 339 { 340 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 341 PetscScalar *v; 342 PetscInt j,k; 343 344 PetscFunctionBegin; 345 PetscCall(MatDenseGetArray(A,&v)); 346 k = PetscMin(A->rmap->n,A->cmap->n); 347 for (j=0; j<k; j++) v[j+j*a->lda] += alpha; 348 PetscCall(PetscLogFlops(k)); 349 PetscCall(MatDenseRestoreArray(A,&v)); 350 PetscFunctionReturn(0); 351 } 352 353 static PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool *fl) 354 { 355 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 356 PetscInt i,j,m = A->rmap->n,N = a->lda; 357 const PetscScalar *v; 358 359 PetscFunctionBegin; 360 *fl = PETSC_FALSE; 361 if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0); 362 PetscCall(MatDenseGetArrayRead(A,&v)); 363 for (i=0; i<m; i++) { 364 for (j=i; j<m; j++) { 365 if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) { 366 goto restore; 367 } 368 } 369 } 370 *fl = PETSC_TRUE; 371 restore: 372 PetscCall(MatDenseRestoreArrayRead(A,&v)); 373 PetscFunctionReturn(0); 374 } 375 376 static PetscErrorCode MatIsSymmetric_SeqDense(Mat A,PetscReal rtol,PetscBool *fl) 377 { 378 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 379 PetscInt i,j,m = A->rmap->n,N = a->lda; 380 const PetscScalar *v; 381 382 PetscFunctionBegin; 383 *fl = PETSC_FALSE; 384 if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0); 385 PetscCall(MatDenseGetArrayRead(A,&v)); 386 for (i=0; i<m; i++) { 387 for (j=i; j<m; j++) { 388 if (PetscAbsScalar(v[i+j*N] - v[j+i*N]) > rtol) { 389 goto restore; 390 } 391 } 392 } 393 *fl = PETSC_TRUE; 394 restore: 395 PetscCall(MatDenseRestoreArrayRead(A,&v)); 396 PetscFunctionReturn(0); 397 } 398 399 PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues) 400 { 401 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 402 PetscInt lda = (PetscInt)mat->lda,j,m,nlda = lda; 403 PetscBool isdensecpu; 404 405 PetscFunctionBegin; 406 PetscCall(PetscLayoutReference(A->rmap,&newi->rmap)); 407 PetscCall(PetscLayoutReference(A->cmap,&newi->cmap)); 408 if (cpvalues == MAT_SHARE_NONZERO_PATTERN) { /* propagate LDA */ 409 PetscCall(MatDenseSetLDA(newi,lda)); 410 } 411 PetscCall(PetscObjectTypeCompare((PetscObject)newi,MATSEQDENSE,&isdensecpu)); 412 if (isdensecpu) PetscCall(MatSeqDenseSetPreallocation(newi,NULL)); 413 if (cpvalues == MAT_COPY_VALUES) { 414 const PetscScalar *av; 415 PetscScalar *v; 416 417 PetscCall(MatDenseGetArrayRead(A,&av)); 418 PetscCall(MatDenseGetArrayWrite(newi,&v)); 419 PetscCall(MatDenseGetLDA(newi,&nlda)); 420 m = A->rmap->n; 421 if (lda>m || nlda>m) { 422 for (j=0; j<A->cmap->n; j++) { 423 PetscCall(PetscArraycpy(v+j*nlda,av+j*lda,m)); 424 } 425 } else { 426 PetscCall(PetscArraycpy(v,av,A->rmap->n*A->cmap->n)); 427 } 428 PetscCall(MatDenseRestoreArrayWrite(newi,&v)); 429 PetscCall(MatDenseRestoreArrayRead(A,&av)); 430 } 431 PetscFunctionReturn(0); 432 } 433 434 PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 435 { 436 PetscFunctionBegin; 437 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),newmat)); 438 PetscCall(MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n)); 439 PetscCall(MatSetType(*newmat,((PetscObject)A)->type_name)); 440 PetscCall(MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues)); 441 PetscFunctionReturn(0); 442 } 443 444 static PetscErrorCode MatSolve_SeqDense_Internal_LU(Mat A, PetscScalar *x, PetscBLASInt ldx, PetscBLASInt m, PetscBLASInt nrhs, PetscBLASInt k, PetscBool T) 445 { 446 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 447 PetscBLASInt info; 448 449 PetscFunctionBegin; 450 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 451 PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_(T ? "T" : "N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info)); 452 PetscCall(PetscFPTrapPop()); 453 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 454 PetscCall(PetscLogFlops(nrhs*(2.0*m*m - m))); 455 PetscFunctionReturn(0); 456 } 457 458 static PetscErrorCode MatConjugate_SeqDense(Mat); 459 460 static PetscErrorCode MatSolve_SeqDense_Internal_Cholesky(Mat A, PetscScalar *x, PetscBLASInt ldx, PetscBLASInt m, PetscBLASInt nrhs, PetscBLASInt k, PetscBool T) 461 { 462 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 463 PetscBLASInt info; 464 465 PetscFunctionBegin; 466 if (A->spd) { 467 if (PetscDefined(USE_COMPLEX) && T) PetscCall(MatConjugate_SeqDense(A)); 468 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 469 PetscStackCallBLAS("LAPACKpotrs",LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info)); 470 PetscCall(PetscFPTrapPop()); 471 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve"); 472 if (PetscDefined(USE_COMPLEX) && T) PetscCall(MatConjugate_SeqDense(A)); 473 #if defined(PETSC_USE_COMPLEX) 474 } else if (A->hermitian) { 475 if (T) PetscCall(MatConjugate_SeqDense(A)); 476 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 477 PetscStackCallBLAS("LAPACKhetrs",LAPACKhetrs_("L",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info)); 478 PetscCall(PetscFPTrapPop()); 479 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"HETRS Bad solve"); 480 if (T) PetscCall(MatConjugate_SeqDense(A)); 481 #endif 482 } else { /* symmetric case */ 483 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 484 PetscStackCallBLAS("LAPACKsytrs",LAPACKsytrs_("L",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info)); 485 PetscCall(PetscFPTrapPop()); 486 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"SYTRS Bad solve"); 487 } 488 PetscCall(PetscLogFlops(nrhs*(2.0*m*m - m))); 489 PetscFunctionReturn(0); 490 } 491 492 static PetscErrorCode MatSolve_SeqDense_Internal_QR(Mat A, PetscScalar *x, PetscBLASInt ldx, PetscBLASInt m, PetscBLASInt nrhs, PetscBLASInt k) 493 { 494 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 495 PetscBLASInt info; 496 char trans; 497 498 PetscFunctionBegin; 499 if (PetscDefined(USE_COMPLEX)) { 500 trans = 'C'; 501 } else { 502 trans = 'T'; 503 } 504 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 505 PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("L", &trans, &m,&nrhs,&mat->rank,mat->v,&mat->lda,mat->tau,x,&ldx,mat->fwork,&mat->lfwork,&info)); 506 PetscCall(PetscFPTrapPop()); 507 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"ORMQR - Bad orthogonal transform"); 508 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 509 PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U", "N", "N", &mat->rank,&nrhs,mat->v,&mat->lda,x,&ldx,&info)); 510 PetscCall(PetscFPTrapPop()); 511 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"TRTRS - Bad triangular solve"); 512 for (PetscInt j = 0; j < nrhs; j++) { 513 for (PetscInt i = mat->rank; i < k; i++) { 514 x[j*ldx + i] = 0.; 515 } 516 } 517 PetscCall(PetscLogFlops(nrhs*(4.0*m*mat->rank - PetscSqr(mat->rank)))); 518 PetscFunctionReturn(0); 519 } 520 521 static PetscErrorCode MatSolveTranspose_SeqDense_Internal_QR(Mat A, PetscScalar *x, PetscBLASInt ldx, PetscBLASInt m, PetscBLASInt nrhs, PetscBLASInt k) 522 { 523 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 524 PetscBLASInt info; 525 526 PetscFunctionBegin; 527 if (A->rmap->n == A->cmap->n && mat->rank == A->rmap->n) { 528 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 529 PetscStackCallBLAS("LAPACKtrtrs",LAPACKtrtrs_("U", "T", "N", &m,&nrhs,mat->v,&mat->lda,x,&ldx,&info)); 530 PetscCall(PetscFPTrapPop()); 531 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"TRTRS - Bad triangular solve"); 532 if (PetscDefined(USE_COMPLEX)) PetscCall(MatConjugate_SeqDense(A)); 533 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 534 PetscStackCallBLAS("LAPACKormqr",LAPACKormqr_("L", "N", &m,&nrhs,&mat->rank,mat->v,&mat->lda,mat->tau,x,&ldx,mat->fwork,&mat->lfwork,&info)); 535 PetscCall(PetscFPTrapPop()); 536 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"ORMQR - Bad orthogonal transform"); 537 if (PetscDefined(USE_COMPLEX)) PetscCall(MatConjugate_SeqDense(A)); 538 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"QR factored matrix cannot be used for transpose solve"); 539 PetscCall(PetscLogFlops(nrhs*(4.0*m*mat->rank - PetscSqr(mat->rank)))); 540 PetscFunctionReturn(0); 541 } 542 543 static PetscErrorCode MatSolve_SeqDense_SetUp(Mat A, Vec xx, Vec yy, PetscScalar **_y, PetscBLASInt *_m, PetscBLASInt *_k) 544 { 545 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 546 PetscScalar *y; 547 PetscBLASInt m=0, k=0; 548 549 PetscFunctionBegin; 550 PetscCall(PetscBLASIntCast(A->rmap->n,&m)); 551 PetscCall(PetscBLASIntCast(A->cmap->n,&k)); 552 if (k < m) { 553 PetscCall(VecCopy(xx, mat->qrrhs)); 554 PetscCall(VecGetArray(mat->qrrhs,&y)); 555 } else { 556 PetscCall(VecCopy(xx, yy)); 557 PetscCall(VecGetArray(yy,&y)); 558 } 559 *_y = y; 560 *_k = k; 561 *_m = m; 562 PetscFunctionReturn(0); 563 } 564 565 static PetscErrorCode MatSolve_SeqDense_TearDown(Mat A, Vec xx, Vec yy, PetscScalar **_y, PetscBLASInt *_m, PetscBLASInt *_k) 566 { 567 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 568 PetscScalar *y = NULL; 569 PetscBLASInt m, k; 570 571 PetscFunctionBegin; 572 y = *_y; 573 *_y = NULL; 574 k = *_k; 575 m = *_m; 576 if (k < m) { 577 PetscScalar *yv; 578 PetscCall(VecGetArray(yy,&yv)); 579 PetscCall(PetscArraycpy(yv, y, k)); 580 PetscCall(VecRestoreArray(yy,&yv)); 581 PetscCall(VecRestoreArray(mat->qrrhs, &y)); 582 } else { 583 PetscCall(VecRestoreArray(yy,&y)); 584 } 585 PetscFunctionReturn(0); 586 } 587 588 static PetscErrorCode MatSolve_SeqDense_LU(Mat A, Vec xx, Vec yy) 589 { 590 PetscScalar *y = NULL; 591 PetscBLASInt m = 0, k = 0; 592 593 PetscFunctionBegin; 594 PetscCall(MatSolve_SeqDense_SetUp(A, xx, yy, &y, &m, &k)); 595 PetscCall(MatSolve_SeqDense_Internal_LU(A, y, m, m, 1, k, PETSC_FALSE)); 596 PetscCall(MatSolve_SeqDense_TearDown(A, xx, yy, &y, &m, &k)); 597 PetscFunctionReturn(0); 598 } 599 600 static PetscErrorCode MatSolveTranspose_SeqDense_LU(Mat A, Vec xx, Vec yy) 601 { 602 PetscScalar *y = NULL; 603 PetscBLASInt m = 0, k = 0; 604 605 PetscFunctionBegin; 606 PetscCall(MatSolve_SeqDense_SetUp(A, xx, yy, &y, &m, &k)); 607 PetscCall(MatSolve_SeqDense_Internal_LU(A, y, m, m, 1, k, PETSC_TRUE)); 608 PetscCall(MatSolve_SeqDense_TearDown(A, xx, yy, &y, &m, &k)); 609 PetscFunctionReturn(0); 610 } 611 612 static PetscErrorCode MatSolve_SeqDense_Cholesky(Mat A, Vec xx, Vec yy) 613 { 614 PetscScalar *y = NULL; 615 PetscBLASInt m = 0, k = 0; 616 617 PetscFunctionBegin; 618 PetscCall(MatSolve_SeqDense_SetUp(A, xx, yy, &y, &m, &k)); 619 PetscCall(MatSolve_SeqDense_Internal_Cholesky(A, y, m, m, 1, k, PETSC_FALSE)); 620 PetscCall(MatSolve_SeqDense_TearDown(A, xx, yy, &y, &m, &k)); 621 PetscFunctionReturn(0); 622 } 623 624 static PetscErrorCode MatSolveTranspose_SeqDense_Cholesky(Mat A, Vec xx, Vec yy) 625 { 626 PetscScalar *y = NULL; 627 PetscBLASInt m = 0, k = 0; 628 629 PetscFunctionBegin; 630 PetscCall(MatSolve_SeqDense_SetUp(A, xx, yy, &y, &m, &k)); 631 PetscCall(MatSolve_SeqDense_Internal_Cholesky(A, y, m, m, 1, k, PETSC_TRUE)); 632 PetscCall(MatSolve_SeqDense_TearDown(A, xx, yy, &y, &m, &k)); 633 PetscFunctionReturn(0); 634 } 635 636 static PetscErrorCode MatSolve_SeqDense_QR(Mat A, Vec xx, Vec yy) 637 { 638 PetscScalar *y = NULL; 639 PetscBLASInt m = 0, k = 0; 640 641 PetscFunctionBegin; 642 PetscCall(MatSolve_SeqDense_SetUp(A, xx, yy, &y, &m, &k)); 643 PetscCall(MatSolve_SeqDense_Internal_QR(A, y, PetscMax(m,k), m, 1, k)); 644 PetscCall(MatSolve_SeqDense_TearDown(A, xx, yy, &y, &m, &k)); 645 PetscFunctionReturn(0); 646 } 647 648 static PetscErrorCode MatSolveTranspose_SeqDense_QR(Mat A, Vec xx, Vec yy) 649 { 650 PetscScalar *y = NULL; 651 PetscBLASInt m = 0, k = 0; 652 653 PetscFunctionBegin; 654 PetscCall(MatSolve_SeqDense_SetUp(A, xx, yy, &y, &m, &k)); 655 PetscCall(MatSolveTranspose_SeqDense_Internal_QR(A, y, PetscMax(m,k), m, 1, k)); 656 PetscCall(MatSolve_SeqDense_TearDown(A, xx, yy, &y, &m, &k)); 657 PetscFunctionReturn(0); 658 } 659 660 static PetscErrorCode MatMatSolve_SeqDense_SetUp(Mat A, Mat B, Mat X, PetscScalar **_y, PetscBLASInt *_ldy, PetscBLASInt *_m, PetscBLASInt *_nrhs, PetscBLASInt *_k) 661 { 662 const PetscScalar *b; 663 PetscScalar *y; 664 PetscInt n, _ldb, _ldx; 665 PetscBLASInt nrhs=0,m=0,k=0,ldb=0,ldx=0,ldy=0; 666 667 PetscFunctionBegin; 668 *_ldy=0; *_m=0; *_nrhs=0; *_k=0; *_y = NULL; 669 PetscCall(PetscBLASIntCast(A->rmap->n,&m)); 670 PetscCall(PetscBLASIntCast(A->cmap->n,&k)); 671 PetscCall(MatGetSize(B,NULL,&n)); 672 PetscCall(PetscBLASIntCast(n,&nrhs)); 673 PetscCall(MatDenseGetLDA(B,&_ldb)); 674 PetscCall(PetscBLASIntCast(_ldb, &ldb)); 675 PetscCall(MatDenseGetLDA(X,&_ldx)); 676 PetscCall(PetscBLASIntCast(_ldx, &ldx)); 677 if (ldx < m) { 678 PetscCall(MatDenseGetArrayRead(B,&b)); 679 PetscCall(PetscMalloc1(nrhs * m, &y)); 680 if (ldb == m) { 681 PetscCall(PetscArraycpy(y,b,ldb*nrhs)); 682 } else { 683 for (PetscInt j = 0; j < nrhs; j++) { 684 PetscCall(PetscArraycpy(&y[j*m],&b[j*ldb],m)); 685 } 686 } 687 ldy = m; 688 PetscCall(MatDenseRestoreArrayRead(B,&b)); 689 } else { 690 if (ldb == ldx) { 691 PetscCall(MatCopy(B, X, SAME_NONZERO_PATTERN)); 692 PetscCall(MatDenseGetArray(X,&y)); 693 } else { 694 PetscCall(MatDenseGetArray(X,&y)); 695 PetscCall(MatDenseGetArrayRead(B,&b)); 696 for (PetscInt j = 0; j < nrhs; j++) { 697 PetscCall(PetscArraycpy(&y[j*ldx],&b[j*ldb],m)); 698 } 699 PetscCall(MatDenseRestoreArrayRead(B,&b)); 700 } 701 ldy = ldx; 702 } 703 *_y = y; 704 *_ldy = ldy; 705 *_k = k; 706 *_m = m; 707 *_nrhs = nrhs; 708 PetscFunctionReturn(0); 709 } 710 711 static PetscErrorCode MatMatSolve_SeqDense_TearDown(Mat A, Mat B, Mat X, PetscScalar **_y, PetscBLASInt *_ldy, PetscBLASInt *_m, PetscBLASInt *_nrhs, PetscBLASInt *_k) 712 { 713 PetscScalar *y; 714 PetscInt _ldx; 715 PetscBLASInt k,ldy,nrhs,ldx=0; 716 717 PetscFunctionBegin; 718 y = *_y; 719 *_y = NULL; 720 k = *_k; 721 ldy = *_ldy; 722 nrhs = *_nrhs; 723 PetscCall(MatDenseGetLDA(X,&_ldx)); 724 PetscCall(PetscBLASIntCast(_ldx, &ldx)); 725 if (ldx != ldy) { 726 PetscScalar *xv; 727 PetscCall(MatDenseGetArray(X,&xv)); 728 for (PetscInt j = 0; j < nrhs; j++) { 729 PetscCall(PetscArraycpy(&xv[j*ldx],&y[j*ldy],k)); 730 } 731 PetscCall(MatDenseRestoreArray(X,&xv)); 732 PetscCall(PetscFree(y)); 733 } else { 734 PetscCall(MatDenseRestoreArray(X,&y)); 735 } 736 PetscFunctionReturn(0); 737 } 738 739 static PetscErrorCode MatMatSolve_SeqDense_LU(Mat A, Mat B, Mat X) 740 { 741 PetscScalar *y; 742 PetscBLASInt m, k, ldy, nrhs; 743 744 PetscFunctionBegin; 745 PetscCall(MatMatSolve_SeqDense_SetUp(A, B, X, &y, &ldy, &m, &nrhs, &k)); 746 PetscCall(MatSolve_SeqDense_Internal_LU(A, y, ldy, m, nrhs, k, PETSC_FALSE)); 747 PetscCall(MatMatSolve_SeqDense_TearDown(A, B, X, &y, &ldy, &m, &nrhs, &k)); 748 PetscFunctionReturn(0); 749 } 750 751 static PetscErrorCode MatMatSolveTranspose_SeqDense_LU(Mat A, Mat B, Mat X) 752 { 753 PetscScalar *y; 754 PetscBLASInt m, k, ldy, nrhs; 755 756 PetscFunctionBegin; 757 PetscCall(MatMatSolve_SeqDense_SetUp(A, B, X, &y, &ldy, &m, &nrhs, &k)); 758 PetscCall(MatSolve_SeqDense_Internal_LU(A, y, ldy, m, nrhs, k, PETSC_TRUE)); 759 PetscCall(MatMatSolve_SeqDense_TearDown(A, B, X, &y, &ldy, &m, &nrhs, &k)); 760 PetscFunctionReturn(0); 761 } 762 763 static PetscErrorCode MatMatSolve_SeqDense_Cholesky(Mat A, Mat B, Mat X) 764 { 765 PetscScalar *y; 766 PetscBLASInt m, k, ldy, nrhs; 767 768 PetscFunctionBegin; 769 PetscCall(MatMatSolve_SeqDense_SetUp(A, B, X, &y, &ldy, &m, &nrhs, &k)); 770 PetscCall(MatSolve_SeqDense_Internal_Cholesky(A, y, ldy, m, nrhs, k, PETSC_FALSE)); 771 PetscCall(MatMatSolve_SeqDense_TearDown(A, B, X, &y, &ldy, &m, &nrhs, &k)); 772 PetscFunctionReturn(0); 773 } 774 775 static PetscErrorCode MatMatSolveTranspose_SeqDense_Cholesky(Mat A, Mat B, Mat X) 776 { 777 PetscScalar *y; 778 PetscBLASInt m, k, ldy, nrhs; 779 780 PetscFunctionBegin; 781 PetscCall(MatMatSolve_SeqDense_SetUp(A, B, X, &y, &ldy, &m, &nrhs, &k)); 782 PetscCall(MatSolve_SeqDense_Internal_Cholesky(A, y, ldy, m, nrhs, k, PETSC_TRUE)); 783 PetscCall(MatMatSolve_SeqDense_TearDown(A, B, X, &y, &ldy, &m, &nrhs, &k)); 784 PetscFunctionReturn(0); 785 } 786 787 static PetscErrorCode MatMatSolve_SeqDense_QR(Mat A, Mat B, Mat X) 788 { 789 PetscScalar *y; 790 PetscBLASInt m, k, ldy, nrhs; 791 792 PetscFunctionBegin; 793 PetscCall(MatMatSolve_SeqDense_SetUp(A, B, X, &y, &ldy, &m, &nrhs, &k)); 794 PetscCall(MatSolve_SeqDense_Internal_QR(A, y, ldy, m, nrhs, k)); 795 PetscCall(MatMatSolve_SeqDense_TearDown(A, B, X, &y, &ldy, &m, &nrhs, &k)); 796 PetscFunctionReturn(0); 797 } 798 799 static PetscErrorCode MatMatSolveTranspose_SeqDense_QR(Mat A, Mat B, Mat X) 800 { 801 PetscScalar *y; 802 PetscBLASInt m, k, ldy, nrhs; 803 804 PetscFunctionBegin; 805 PetscCall(MatMatSolve_SeqDense_SetUp(A, B, X, &y, &ldy, &m, &nrhs, &k)); 806 PetscCall(MatSolveTranspose_SeqDense_Internal_QR(A, y, ldy, m, nrhs, k)); 807 PetscCall(MatMatSolve_SeqDense_TearDown(A, B, X, &y, &ldy, &m, &nrhs, &k)); 808 PetscFunctionReturn(0); 809 } 810 811 static PetscErrorCode MatConjugate_SeqDense(Mat); 812 813 /* ---------------------------------------------------------------*/ 814 /* COMMENT: I have chosen to hide row permutation in the pivots, 815 rather than put it in the Mat->row slot.*/ 816 PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo) 817 { 818 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 819 PetscBLASInt n,m,info; 820 821 PetscFunctionBegin; 822 PetscCall(PetscBLASIntCast(A->cmap->n,&n)); 823 PetscCall(PetscBLASIntCast(A->rmap->n,&m)); 824 if (!mat->pivots) { 825 PetscCall(PetscMalloc1(A->rmap->n,&mat->pivots)); 826 PetscCall(PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt))); 827 } 828 if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 829 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 830 PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info)); 831 PetscCall(PetscFPTrapPop()); 832 833 PetscCheck(info>=0,PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization"); 834 PetscCheck(info<=0,PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 835 836 A->ops->solve = MatSolve_SeqDense_LU; 837 A->ops->matsolve = MatMatSolve_SeqDense_LU; 838 A->ops->solvetranspose = MatSolveTranspose_SeqDense_LU; 839 A->ops->matsolvetranspose = MatMatSolveTranspose_SeqDense_LU; 840 A->factortype = MAT_FACTOR_LU; 841 842 PetscCall(PetscFree(A->solvertype)); 843 PetscCall(PetscStrallocpy(MATSOLVERPETSC,&A->solvertype)); 844 845 PetscCall(PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3)); 846 PetscFunctionReturn(0); 847 } 848 849 static PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 850 { 851 MatFactorInfo info; 852 853 PetscFunctionBegin; 854 PetscCall(MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES)); 855 PetscCall((*fact->ops->lufactor)(fact,NULL,NULL,&info)); 856 PetscFunctionReturn(0); 857 } 858 859 PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info) 860 { 861 PetscFunctionBegin; 862 fact->preallocated = PETSC_TRUE; 863 fact->assembled = PETSC_TRUE; 864 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense; 865 PetscFunctionReturn(0); 866 } 867 868 /* Cholesky as L*L^T or L*D*L^T and the symmetric/hermitian complex variants */ 869 PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo) 870 { 871 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 872 PetscBLASInt info,n; 873 874 PetscFunctionBegin; 875 PetscCall(PetscBLASIntCast(A->cmap->n,&n)); 876 if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 877 if (A->spd) { 878 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 879 PetscStackCallBLAS("LAPACKpotrf",LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info)); 880 PetscCall(PetscFPTrapPop()); 881 #if defined(PETSC_USE_COMPLEX) 882 } else if (A->hermitian) { 883 if (!mat->pivots) { 884 PetscCall(PetscMalloc1(A->rmap->n,&mat->pivots)); 885 PetscCall(PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt))); 886 } 887 if (!mat->fwork) { 888 PetscScalar dummy; 889 890 mat->lfwork = -1; 891 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 892 PetscStackCallBLAS("LAPACKhetrf",LAPACKhetrf_("L",&n,mat->v,&mat->lda,mat->pivots,&dummy,&mat->lfwork,&info)); 893 PetscCall(PetscFPTrapPop()); 894 mat->lfwork = (PetscInt)PetscRealPart(dummy); 895 PetscCall(PetscMalloc1(mat->lfwork,&mat->fwork)); 896 PetscCall(PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt))); 897 } 898 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 899 PetscStackCallBLAS("LAPACKhetrf",LAPACKhetrf_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info)); 900 PetscCall(PetscFPTrapPop()); 901 #endif 902 } else { /* symmetric case */ 903 if (!mat->pivots) { 904 PetscCall(PetscMalloc1(A->rmap->n,&mat->pivots)); 905 PetscCall(PetscLogObjectMemory((PetscObject)A,A->rmap->n*sizeof(PetscBLASInt))); 906 } 907 if (!mat->fwork) { 908 PetscScalar dummy; 909 910 mat->lfwork = -1; 911 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 912 PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&n,mat->v,&mat->lda,mat->pivots,&dummy,&mat->lfwork,&info)); 913 PetscCall(PetscFPTrapPop()); 914 mat->lfwork = (PetscInt)PetscRealPart(dummy); 915 PetscCall(PetscMalloc1(mat->lfwork,&mat->fwork)); 916 PetscCall(PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt))); 917 } 918 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 919 PetscStackCallBLAS("LAPACKsytrf",LAPACKsytrf_("L",&n,mat->v,&mat->lda,mat->pivots,mat->fwork,&mat->lfwork,&info)); 920 PetscCall(PetscFPTrapPop()); 921 } 922 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %" PetscInt_FMT,(PetscInt)info-1); 923 924 A->ops->solve = MatSolve_SeqDense_Cholesky; 925 A->ops->matsolve = MatMatSolve_SeqDense_Cholesky; 926 A->ops->solvetranspose = MatSolveTranspose_SeqDense_Cholesky; 927 A->ops->matsolvetranspose = MatMatSolveTranspose_SeqDense_Cholesky; 928 A->factortype = MAT_FACTOR_CHOLESKY; 929 930 PetscCall(PetscFree(A->solvertype)); 931 PetscCall(PetscStrallocpy(MATSOLVERPETSC,&A->solvertype)); 932 933 PetscCall(PetscLogFlops((1.0*A->cmap->n*A->cmap->n*A->cmap->n)/3.0)); 934 PetscFunctionReturn(0); 935 } 936 937 static PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 938 { 939 MatFactorInfo info; 940 941 PetscFunctionBegin; 942 info.fill = 1.0; 943 944 PetscCall(MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES)); 945 PetscCall((*fact->ops->choleskyfactor)(fact,NULL,&info)); 946 PetscFunctionReturn(0); 947 } 948 949 PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info) 950 { 951 PetscFunctionBegin; 952 fact->assembled = PETSC_TRUE; 953 fact->preallocated = PETSC_TRUE; 954 fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense; 955 PetscFunctionReturn(0); 956 } 957 958 PetscErrorCode MatQRFactor_SeqDense(Mat A,IS col,const MatFactorInfo *minfo) 959 { 960 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 961 PetscBLASInt n,m,info, min, max; 962 963 PetscFunctionBegin; 964 PetscCall(PetscBLASIntCast(A->cmap->n,&n)); 965 PetscCall(PetscBLASIntCast(A->rmap->n,&m)); 966 max = PetscMax(m, n); 967 min = PetscMin(m, n); 968 if (!mat->tau) { 969 PetscCall(PetscMalloc1(min,&mat->tau)); 970 PetscCall(PetscLogObjectMemory((PetscObject)A,min*sizeof(PetscScalar))); 971 } 972 if (!mat->pivots) { 973 PetscCall(PetscMalloc1(n,&mat->pivots)); 974 PetscCall(PetscLogObjectMemory((PetscObject)A,n*sizeof(PetscScalar))); 975 } 976 if (!mat->qrrhs) { 977 PetscCall(MatCreateVecs(A, NULL, &(mat->qrrhs))); 978 } 979 if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 980 if (!mat->fwork) { 981 PetscScalar dummy; 982 983 mat->lfwork = -1; 984 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 985 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&m,&n,mat->v,&mat->lda,mat->tau,&dummy,&mat->lfwork,&info)); 986 PetscCall(PetscFPTrapPop()); 987 mat->lfwork = (PetscInt)PetscRealPart(dummy); 988 PetscCall(PetscMalloc1(mat->lfwork,&mat->fwork)); 989 PetscCall(PetscLogObjectMemory((PetscObject)A,mat->lfwork*sizeof(PetscBLASInt))); 990 } 991 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 992 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&m,&n,mat->v,&mat->lda,mat->tau,mat->fwork,&mat->lfwork,&info)); 993 PetscCall(PetscFPTrapPop()); 994 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to QR factorization"); 995 // TODO: try to estimate rank or test for and use geqp3 for rank revealing QR. For now just say rank is min of m and n 996 mat->rank = min; 997 998 A->ops->solve = MatSolve_SeqDense_QR; 999 A->ops->matsolve = MatMatSolve_SeqDense_QR; 1000 A->factortype = MAT_FACTOR_QR; 1001 if (m == n) { 1002 A->ops->solvetranspose = MatSolveTranspose_SeqDense_QR; 1003 A->ops->matsolvetranspose = MatMatSolveTranspose_SeqDense_QR; 1004 } 1005 1006 PetscCall(PetscFree(A->solvertype)); 1007 PetscCall(PetscStrallocpy(MATSOLVERPETSC,&A->solvertype)); 1008 1009 PetscCall(PetscLogFlops(2.0*min*min*(max-min/3.0))); 1010 PetscFunctionReturn(0); 1011 } 1012 1013 static PetscErrorCode MatQRFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 1014 { 1015 MatFactorInfo info; 1016 1017 PetscFunctionBegin; 1018 info.fill = 1.0; 1019 1020 PetscCall(MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES)); 1021 PetscUseMethod(fact,"MatQRFactor_C",(Mat,IS,const MatFactorInfo *),(fact,NULL,&info)); 1022 PetscFunctionReturn(0); 1023 } 1024 1025 PetscErrorCode MatQRFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info) 1026 { 1027 PetscFunctionBegin; 1028 fact->assembled = PETSC_TRUE; 1029 fact->preallocated = PETSC_TRUE; 1030 PetscCall(PetscObjectComposeFunction((PetscObject)fact,"MatQRFactorNumeric_C",MatQRFactorNumeric_SeqDense)); 1031 PetscFunctionReturn(0); 1032 } 1033 1034 /* uses LAPACK */ 1035 PETSC_INTERN PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact) 1036 { 1037 PetscFunctionBegin; 1038 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),fact)); 1039 PetscCall(MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n)); 1040 PetscCall(MatSetType(*fact,MATDENSE)); 1041 (*fact)->trivialsymbolic = PETSC_TRUE; 1042 if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU) { 1043 (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense; 1044 (*fact)->ops->ilufactorsymbolic = MatLUFactorSymbolic_SeqDense; 1045 } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 1046 (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense; 1047 } else if (ftype == MAT_FACTOR_QR) { 1048 PetscCall(PetscObjectComposeFunction((PetscObject)(*fact),"MatQRFactorSymbolic_C",MatQRFactorSymbolic_SeqDense)); 1049 } 1050 (*fact)->factortype = ftype; 1051 1052 PetscCall(PetscFree((*fact)->solvertype)); 1053 PetscCall(PetscStrallocpy(MATSOLVERPETSC,&(*fact)->solvertype)); 1054 PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&(*fact)->preferredordering[MAT_FACTOR_LU])); 1055 PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&(*fact)->preferredordering[MAT_FACTOR_ILU])); 1056 PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&(*fact)->preferredordering[MAT_FACTOR_CHOLESKY])); 1057 PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL,(char**)&(*fact)->preferredordering[MAT_FACTOR_ICC])); 1058 PetscFunctionReturn(0); 1059 } 1060 1061 /* ------------------------------------------------------------------*/ 1062 static PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 1063 { 1064 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1065 PetscScalar *x,*v = mat->v,zero = 0.0,xt; 1066 const PetscScalar *b; 1067 PetscInt m = A->rmap->n,i; 1068 PetscBLASInt o = 1,bm = 0; 1069 1070 PetscFunctionBegin; 1071 #if defined(PETSC_HAVE_CUDA) 1072 PetscCheck(A->offloadmask != PETSC_OFFLOAD_GPU,PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented"); 1073 #endif 1074 if (shift == -1) shift = 0.0; /* negative shift indicates do not error on zero diagonal; this code never zeros on zero diagonal */ 1075 PetscCall(PetscBLASIntCast(m,&bm)); 1076 if (flag & SOR_ZERO_INITIAL_GUESS) { 1077 /* this is a hack fix, should have another version without the second BLASdotu */ 1078 PetscCall(VecSet(xx,zero)); 1079 } 1080 PetscCall(VecGetArray(xx,&x)); 1081 PetscCall(VecGetArrayRead(bb,&b)); 1082 its = its*lits; 1083 PetscCheck(its > 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive",its,lits); 1084 while (its--) { 1085 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 1086 for (i=0; i<m; i++) { 1087 PetscStackCallBLAS("BLASdotu",xt = b[i] - BLASdotu_(&bm,v+i,&bm,x,&o)); 1088 x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 1089 } 1090 } 1091 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 1092 for (i=m-1; i>=0; i--) { 1093 PetscStackCallBLAS("BLASdotu",xt = b[i] - BLASdotu_(&bm,v+i,&bm,x,&o)); 1094 x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 1095 } 1096 } 1097 } 1098 PetscCall(VecRestoreArrayRead(bb,&b)); 1099 PetscCall(VecRestoreArray(xx,&x)); 1100 PetscFunctionReturn(0); 1101 } 1102 1103 /* -----------------------------------------------------------------*/ 1104 PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 1105 { 1106 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1107 const PetscScalar *v = mat->v,*x; 1108 PetscScalar *y; 1109 PetscBLASInt m, n,_One=1; 1110 PetscScalar _DOne=1.0,_DZero=0.0; 1111 1112 PetscFunctionBegin; 1113 PetscCall(PetscBLASIntCast(A->rmap->n,&m)); 1114 PetscCall(PetscBLASIntCast(A->cmap->n,&n)); 1115 PetscCall(VecGetArrayRead(xx,&x)); 1116 PetscCall(VecGetArrayWrite(yy,&y)); 1117 if (!A->rmap->n || !A->cmap->n) { 1118 PetscBLASInt i; 1119 for (i=0; i<n; i++) y[i] = 0.0; 1120 } else { 1121 PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One)); 1122 PetscCall(PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n)); 1123 } 1124 PetscCall(VecRestoreArrayRead(xx,&x)); 1125 PetscCall(VecRestoreArrayWrite(yy,&y)); 1126 PetscFunctionReturn(0); 1127 } 1128 1129 PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 1130 { 1131 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1132 PetscScalar *y,_DOne=1.0,_DZero=0.0; 1133 PetscBLASInt m, n, _One=1; 1134 const PetscScalar *v = mat->v,*x; 1135 1136 PetscFunctionBegin; 1137 PetscCall(PetscBLASIntCast(A->rmap->n,&m)); 1138 PetscCall(PetscBLASIntCast(A->cmap->n,&n)); 1139 PetscCall(VecGetArrayRead(xx,&x)); 1140 PetscCall(VecGetArrayWrite(yy,&y)); 1141 if (!A->rmap->n || !A->cmap->n) { 1142 PetscBLASInt i; 1143 for (i=0; i<m; i++) y[i] = 0.0; 1144 } else { 1145 PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One)); 1146 PetscCall(PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n)); 1147 } 1148 PetscCall(VecRestoreArrayRead(xx,&x)); 1149 PetscCall(VecRestoreArrayWrite(yy,&y)); 1150 PetscFunctionReturn(0); 1151 } 1152 1153 PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 1154 { 1155 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1156 const PetscScalar *v = mat->v,*x; 1157 PetscScalar *y,_DOne=1.0; 1158 PetscBLASInt m, n, _One=1; 1159 1160 PetscFunctionBegin; 1161 PetscCall(PetscBLASIntCast(A->rmap->n,&m)); 1162 PetscCall(PetscBLASIntCast(A->cmap->n,&n)); 1163 PetscCall(VecCopy(zz,yy)); 1164 if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 1165 PetscCall(VecGetArrayRead(xx,&x)); 1166 PetscCall(VecGetArray(yy,&y)); 1167 PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One)); 1168 PetscCall(VecRestoreArrayRead(xx,&x)); 1169 PetscCall(VecRestoreArray(yy,&y)); 1170 PetscCall(PetscLogFlops(2.0*A->rmap->n*A->cmap->n)); 1171 PetscFunctionReturn(0); 1172 } 1173 1174 PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 1175 { 1176 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1177 const PetscScalar *v = mat->v,*x; 1178 PetscScalar *y; 1179 PetscBLASInt m, n, _One=1; 1180 PetscScalar _DOne=1.0; 1181 1182 PetscFunctionBegin; 1183 PetscCall(PetscBLASIntCast(A->rmap->n,&m)); 1184 PetscCall(PetscBLASIntCast(A->cmap->n,&n)); 1185 PetscCall(VecCopy(zz,yy)); 1186 if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 1187 PetscCall(VecGetArrayRead(xx,&x)); 1188 PetscCall(VecGetArray(yy,&y)); 1189 PetscStackCallBLAS("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One)); 1190 PetscCall(VecRestoreArrayRead(xx,&x)); 1191 PetscCall(VecRestoreArray(yy,&y)); 1192 PetscCall(PetscLogFlops(2.0*A->rmap->n*A->cmap->n)); 1193 PetscFunctionReturn(0); 1194 } 1195 1196 /* -----------------------------------------------------------------*/ 1197 static PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 1198 { 1199 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1200 PetscInt i; 1201 1202 PetscFunctionBegin; 1203 *ncols = A->cmap->n; 1204 if (cols) { 1205 PetscCall(PetscMalloc1(A->cmap->n,cols)); 1206 for (i=0; i<A->cmap->n; i++) (*cols)[i] = i; 1207 } 1208 if (vals) { 1209 const PetscScalar *v; 1210 1211 PetscCall(MatDenseGetArrayRead(A,&v)); 1212 PetscCall(PetscMalloc1(A->cmap->n,vals)); 1213 v += row; 1214 for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;} 1215 PetscCall(MatDenseRestoreArrayRead(A,&v)); 1216 } 1217 PetscFunctionReturn(0); 1218 } 1219 1220 static PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 1221 { 1222 PetscFunctionBegin; 1223 if (ncols) *ncols = 0; 1224 if (cols) PetscCall(PetscFree(*cols)); 1225 if (vals) PetscCall(PetscFree(*vals)); 1226 PetscFunctionReturn(0); 1227 } 1228 /* ----------------------------------------------------------------*/ 1229 static PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 1230 { 1231 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1232 PetscScalar *av; 1233 PetscInt i,j,idx=0; 1234 #if defined(PETSC_HAVE_CUDA) 1235 PetscOffloadMask oldf; 1236 #endif 1237 1238 PetscFunctionBegin; 1239 PetscCall(MatDenseGetArray(A,&av)); 1240 if (!mat->roworiented) { 1241 if (addv == INSERT_VALUES) { 1242 for (j=0; j<n; j++) { 1243 if (indexn[j] < 0) {idx += m; continue;} 1244 PetscCheck(indexn[j] < A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT,indexn[j],A->cmap->n-1); 1245 for (i=0; i<m; i++) { 1246 if (indexm[i] < 0) {idx++; continue;} 1247 PetscCheck(indexm[i] < A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT,indexm[i],A->rmap->n-1); 1248 av[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 1249 } 1250 } 1251 } else { 1252 for (j=0; j<n; j++) { 1253 if (indexn[j] < 0) {idx += m; continue;} 1254 PetscCheck(indexn[j] < A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT,indexn[j],A->cmap->n-1); 1255 for (i=0; i<m; i++) { 1256 if (indexm[i] < 0) {idx++; continue;} 1257 PetscCheck(indexm[i] < A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT,indexm[i],A->rmap->n-1); 1258 av[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 1259 } 1260 } 1261 } 1262 } else { 1263 if (addv == INSERT_VALUES) { 1264 for (i=0; i<m; i++) { 1265 if (indexm[i] < 0) { idx += n; continue;} 1266 PetscCheck(indexm[i] < A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT,indexm[i],A->rmap->n-1); 1267 for (j=0; j<n; j++) { 1268 if (indexn[j] < 0) { idx++; continue;} 1269 PetscCheck(indexn[j] < A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT,indexn[j],A->cmap->n-1); 1270 av[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 1271 } 1272 } 1273 } else { 1274 for (i=0; i<m; i++) { 1275 if (indexm[i] < 0) { idx += n; continue;} 1276 PetscCheck(indexm[i] < A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT,indexm[i],A->rmap->n-1); 1277 for (j=0; j<n; j++) { 1278 if (indexn[j] < 0) { idx++; continue;} 1279 PetscCheck(indexn[j] < A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %" PetscInt_FMT " max %" PetscInt_FMT,indexn[j],A->cmap->n-1); 1280 av[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 1281 } 1282 } 1283 } 1284 } 1285 /* hack to prevent unneeded copy to the GPU while returning the array */ 1286 #if defined(PETSC_HAVE_CUDA) 1287 oldf = A->offloadmask; 1288 A->offloadmask = PETSC_OFFLOAD_GPU; 1289 #endif 1290 PetscCall(MatDenseRestoreArray(A,&av)); 1291 #if defined(PETSC_HAVE_CUDA) 1292 A->offloadmask = (oldf == PETSC_OFFLOAD_UNALLOCATED ? PETSC_OFFLOAD_UNALLOCATED : PETSC_OFFLOAD_CPU); 1293 #endif 1294 PetscFunctionReturn(0); 1295 } 1296 1297 static PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 1298 { 1299 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1300 const PetscScalar *vv; 1301 PetscInt i,j; 1302 1303 PetscFunctionBegin; 1304 PetscCall(MatDenseGetArrayRead(A,&vv)); 1305 /* row-oriented output */ 1306 for (i=0; i<m; i++) { 1307 if (indexm[i] < 0) {v += n;continue;} 1308 PetscCheck(indexm[i] < A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %" PetscInt_FMT " requested larger than number rows %" PetscInt_FMT,indexm[i],A->rmap->n); 1309 for (j=0; j<n; j++) { 1310 if (indexn[j] < 0) {v++; continue;} 1311 PetscCheck(indexn[j] < A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column %" PetscInt_FMT " requested larger than number columns %" PetscInt_FMT,indexn[j],A->cmap->n); 1312 *v++ = vv[indexn[j]*mat->lda + indexm[i]]; 1313 } 1314 } 1315 PetscCall(MatDenseRestoreArrayRead(A,&vv)); 1316 PetscFunctionReturn(0); 1317 } 1318 1319 /* -----------------------------------------------------------------*/ 1320 1321 PetscErrorCode MatView_Dense_Binary(Mat mat,PetscViewer viewer) 1322 { 1323 PetscBool skipHeader; 1324 PetscViewerFormat format; 1325 PetscInt header[4],M,N,m,lda,i,j,k; 1326 const PetscScalar *v; 1327 PetscScalar *vwork; 1328 1329 PetscFunctionBegin; 1330 PetscCall(PetscViewerSetUp(viewer)); 1331 PetscCall(PetscViewerBinaryGetSkipHeader(viewer,&skipHeader)); 1332 PetscCall(PetscViewerGetFormat(viewer,&format)); 1333 if (skipHeader) format = PETSC_VIEWER_NATIVE; 1334 1335 PetscCall(MatGetSize(mat,&M,&N)); 1336 1337 /* write matrix header */ 1338 header[0] = MAT_FILE_CLASSID; header[1] = M; header[2] = N; 1339 header[3] = (format == PETSC_VIEWER_NATIVE) ? MATRIX_BINARY_FORMAT_DENSE : M*N; 1340 if (!skipHeader) PetscCall(PetscViewerBinaryWrite(viewer,header,4,PETSC_INT)); 1341 1342 PetscCall(MatGetLocalSize(mat,&m,NULL)); 1343 if (format != PETSC_VIEWER_NATIVE) { 1344 PetscInt nnz = m*N, *iwork; 1345 /* store row lengths for each row */ 1346 PetscCall(PetscMalloc1(nnz,&iwork)); 1347 for (i=0; i<m; i++) iwork[i] = N; 1348 PetscCall(PetscViewerBinaryWriteAll(viewer,iwork,m,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT)); 1349 /* store column indices (zero start index) */ 1350 for (k=0, i=0; i<m; i++) 1351 for (j=0; j<N; j++, k++) 1352 iwork[k] = j; 1353 PetscCall(PetscViewerBinaryWriteAll(viewer,iwork,nnz,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT)); 1354 PetscCall(PetscFree(iwork)); 1355 } 1356 /* store matrix values as a dense matrix in row major order */ 1357 PetscCall(PetscMalloc1(m*N,&vwork)); 1358 PetscCall(MatDenseGetArrayRead(mat,&v)); 1359 PetscCall(MatDenseGetLDA(mat,&lda)); 1360 for (k=0, i=0; i<m; i++) 1361 for (j=0; j<N; j++, k++) 1362 vwork[k] = v[i+lda*j]; 1363 PetscCall(MatDenseRestoreArrayRead(mat,&v)); 1364 PetscCall(PetscViewerBinaryWriteAll(viewer,vwork,m*N,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_SCALAR)); 1365 PetscCall(PetscFree(vwork)); 1366 PetscFunctionReturn(0); 1367 } 1368 1369 PetscErrorCode MatLoad_Dense_Binary(Mat mat,PetscViewer viewer) 1370 { 1371 PetscBool skipHeader; 1372 PetscInt header[4],M,N,m,nz,lda,i,j,k; 1373 PetscInt rows,cols; 1374 PetscScalar *v,*vwork; 1375 1376 PetscFunctionBegin; 1377 PetscCall(PetscViewerSetUp(viewer)); 1378 PetscCall(PetscViewerBinaryGetSkipHeader(viewer,&skipHeader)); 1379 1380 if (!skipHeader) { 1381 PetscCall(PetscViewerBinaryRead(viewer,header,4,NULL,PETSC_INT)); 1382 PetscCheck(header[0] == MAT_FILE_CLASSID,PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Not a matrix object in file"); 1383 M = header[1]; N = header[2]; 1384 PetscCheck(M >= 0,PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Matrix row size (%" PetscInt_FMT ") in file is negative",M); 1385 PetscCheck(N >= 0,PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Matrix column size (%" PetscInt_FMT ") in file is negative",N); 1386 nz = header[3]; 1387 PetscCheck(nz == MATRIX_BINARY_FORMAT_DENSE || nz >= 0,PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED,"Unknown matrix format %" PetscInt_FMT " in file",nz); 1388 } else { 1389 PetscCall(MatGetSize(mat,&M,&N)); 1390 PetscCheck(M >= 0 && N >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"Matrix binary file header was skipped, thus the user must specify the global sizes of input matrix"); 1391 nz = MATRIX_BINARY_FORMAT_DENSE; 1392 } 1393 1394 /* setup global sizes if not set */ 1395 if (mat->rmap->N < 0) mat->rmap->N = M; 1396 if (mat->cmap->N < 0) mat->cmap->N = N; 1397 PetscCall(MatSetUp(mat)); 1398 /* check if global sizes are correct */ 1399 PetscCall(MatGetSize(mat,&rows,&cols)); 1400 PetscCheck(M == rows && N == cols,PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different sizes (%" PetscInt_FMT ", %" PetscInt_FMT ") than the input matrix (%" PetscInt_FMT ", %" PetscInt_FMT ")",M,N,rows,cols); 1401 1402 PetscCall(MatGetSize(mat,NULL,&N)); 1403 PetscCall(MatGetLocalSize(mat,&m,NULL)); 1404 PetscCall(MatDenseGetArray(mat,&v)); 1405 PetscCall(MatDenseGetLDA(mat,&lda)); 1406 if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense format */ 1407 PetscInt nnz = m*N; 1408 /* read in matrix values */ 1409 PetscCall(PetscMalloc1(nnz,&vwork)); 1410 PetscCall(PetscViewerBinaryReadAll(viewer,vwork,nnz,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_SCALAR)); 1411 /* store values in column major order */ 1412 for (j=0; j<N; j++) 1413 for (i=0; i<m; i++) 1414 v[i+lda*j] = vwork[i*N+j]; 1415 PetscCall(PetscFree(vwork)); 1416 } else { /* matrix in file is sparse format */ 1417 PetscInt nnz = 0, *rlens, *icols; 1418 /* read in row lengths */ 1419 PetscCall(PetscMalloc1(m,&rlens)); 1420 PetscCall(PetscViewerBinaryReadAll(viewer,rlens,m,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT)); 1421 for (i=0; i<m; i++) nnz += rlens[i]; 1422 /* read in column indices and values */ 1423 PetscCall(PetscMalloc2(nnz,&icols,nnz,&vwork)); 1424 PetscCall(PetscViewerBinaryReadAll(viewer,icols,nnz,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_INT)); 1425 PetscCall(PetscViewerBinaryReadAll(viewer,vwork,nnz,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_SCALAR)); 1426 /* store values in column major order */ 1427 for (k=0, i=0; i<m; i++) 1428 for (j=0; j<rlens[i]; j++, k++) 1429 v[i+lda*icols[k]] = vwork[k]; 1430 PetscCall(PetscFree(rlens)); 1431 PetscCall(PetscFree2(icols,vwork)); 1432 } 1433 PetscCall(MatDenseRestoreArray(mat,&v)); 1434 PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 1435 PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 1436 PetscFunctionReturn(0); 1437 } 1438 1439 PetscErrorCode MatLoad_SeqDense(Mat newMat, PetscViewer viewer) 1440 { 1441 PetscBool isbinary, ishdf5; 1442 1443 PetscFunctionBegin; 1444 PetscValidHeaderSpecific(newMat,MAT_CLASSID,1); 1445 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 1446 /* force binary viewer to load .info file if it has not yet done so */ 1447 PetscCall(PetscViewerSetUp(viewer)); 1448 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary)); 1449 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5, &ishdf5)); 1450 if (isbinary) { 1451 PetscCall(MatLoad_Dense_Binary(newMat,viewer)); 1452 } else if (ishdf5) { 1453 #if defined(PETSC_HAVE_HDF5) 1454 PetscCall(MatLoad_Dense_HDF5(newMat,viewer)); 1455 #else 1456 SETERRQ(PetscObjectComm((PetscObject)newMat),PETSC_ERR_SUP,"HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5"); 1457 #endif 1458 } else { 1459 SETERRQ(PetscObjectComm((PetscObject)newMat),PETSC_ERR_SUP,"Viewer type %s not yet supported for reading %s matrices",((PetscObject)viewer)->type_name,((PetscObject)newMat)->type_name); 1460 } 1461 PetscFunctionReturn(0); 1462 } 1463 1464 static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 1465 { 1466 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1467 PetscInt i,j; 1468 const char *name; 1469 PetscScalar *v,*av; 1470 PetscViewerFormat format; 1471 #if defined(PETSC_USE_COMPLEX) 1472 PetscBool allreal = PETSC_TRUE; 1473 #endif 1474 1475 PetscFunctionBegin; 1476 PetscCall(MatDenseGetArrayRead(A,(const PetscScalar**)&av)); 1477 PetscCall(PetscViewerGetFormat(viewer,&format)); 1478 if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 1479 PetscFunctionReturn(0); /* do nothing for now */ 1480 } else if (format == PETSC_VIEWER_ASCII_COMMON) { 1481 PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); 1482 for (i=0; i<A->rmap->n; i++) { 1483 v = av + i; 1484 PetscCall(PetscViewerASCIIPrintf(viewer,"row %" PetscInt_FMT ":",i)); 1485 for (j=0; j<A->cmap->n; j++) { 1486 #if defined(PETSC_USE_COMPLEX) 1487 if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 1488 PetscCall(PetscViewerASCIIPrintf(viewer," (%" PetscInt_FMT ", %g + %g i) ",j,(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v))); 1489 } else if (PetscRealPart(*v)) { 1490 PetscCall(PetscViewerASCIIPrintf(viewer," (%" PetscInt_FMT ", %g) ",j,(double)PetscRealPart(*v))); 1491 } 1492 #else 1493 if (*v) { 1494 PetscCall(PetscViewerASCIIPrintf(viewer," (%" PetscInt_FMT ", %g) ",j,(double)*v)); 1495 } 1496 #endif 1497 v += a->lda; 1498 } 1499 PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); 1500 } 1501 PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); 1502 } else { 1503 PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_FALSE)); 1504 #if defined(PETSC_USE_COMPLEX) 1505 /* determine if matrix has all real values */ 1506 v = av; 1507 for (i=0; i<A->rmap->n*A->cmap->n; i++) { 1508 if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break;} 1509 } 1510 #endif 1511 if (format == PETSC_VIEWER_ASCII_MATLAB) { 1512 PetscCall(PetscObjectGetName((PetscObject)A,&name)); 1513 PetscCall(PetscViewerASCIIPrintf(viewer,"%% Size = %" PetscInt_FMT " %" PetscInt_FMT " \n",A->rmap->n,A->cmap->n)); 1514 PetscCall(PetscViewerASCIIPrintf(viewer,"%s = zeros(%" PetscInt_FMT ",%" PetscInt_FMT ");\n",name,A->rmap->n,A->cmap->n)); 1515 PetscCall(PetscViewerASCIIPrintf(viewer,"%s = [\n",name)); 1516 } 1517 1518 for (i=0; i<A->rmap->n; i++) { 1519 v = av + i; 1520 for (j=0; j<A->cmap->n; j++) { 1521 #if defined(PETSC_USE_COMPLEX) 1522 if (allreal) { 1523 PetscCall(PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v))); 1524 } else { 1525 PetscCall(PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v))); 1526 } 1527 #else 1528 PetscCall(PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v)); 1529 #endif 1530 v += a->lda; 1531 } 1532 PetscCall(PetscViewerASCIIPrintf(viewer,"\n")); 1533 } 1534 if (format == PETSC_VIEWER_ASCII_MATLAB) { 1535 PetscCall(PetscViewerASCIIPrintf(viewer,"];\n")); 1536 } 1537 PetscCall(PetscViewerASCIIUseTabs(viewer,PETSC_TRUE)); 1538 } 1539 PetscCall(MatDenseRestoreArrayRead(A,(const PetscScalar**)&av)); 1540 PetscCall(PetscViewerFlush(viewer)); 1541 PetscFunctionReturn(0); 1542 } 1543 1544 #include <petscdraw.h> 1545 static PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 1546 { 1547 Mat A = (Mat) Aa; 1548 PetscInt m = A->rmap->n,n = A->cmap->n,i,j; 1549 int color = PETSC_DRAW_WHITE; 1550 const PetscScalar *v; 1551 PetscViewer viewer; 1552 PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r; 1553 PetscViewerFormat format; 1554 1555 PetscFunctionBegin; 1556 PetscCall(PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer)); 1557 PetscCall(PetscViewerGetFormat(viewer,&format)); 1558 PetscCall(PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr)); 1559 1560 /* Loop over matrix elements drawing boxes */ 1561 PetscCall(MatDenseGetArrayRead(A,&v)); 1562 if (format != PETSC_VIEWER_DRAW_CONTOUR) { 1563 PetscDrawCollectiveBegin(draw); 1564 /* Blue for negative and Red for positive */ 1565 for (j = 0; j < n; j++) { 1566 x_l = j; x_r = x_l + 1.0; 1567 for (i = 0; i < m; i++) { 1568 y_l = m - i - 1.0; 1569 y_r = y_l + 1.0; 1570 if (PetscRealPart(v[j*m+i]) > 0.) color = PETSC_DRAW_RED; 1571 else if (PetscRealPart(v[j*m+i]) < 0.) color = PETSC_DRAW_BLUE; 1572 else continue; 1573 PetscCall(PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color)); 1574 } 1575 } 1576 PetscDrawCollectiveEnd(draw); 1577 } else { 1578 /* use contour shading to indicate magnitude of values */ 1579 /* first determine max of all nonzero values */ 1580 PetscReal minv = 0.0, maxv = 0.0; 1581 PetscDraw popup; 1582 1583 for (i=0; i < m*n; i++) { 1584 if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 1585 } 1586 if (minv >= maxv) maxv = minv + PETSC_SMALL; 1587 PetscCall(PetscDrawGetPopup(draw,&popup)); 1588 PetscCall(PetscDrawScalePopup(popup,minv,maxv)); 1589 1590 PetscDrawCollectiveBegin(draw); 1591 for (j=0; j<n; j++) { 1592 x_l = j; 1593 x_r = x_l + 1.0; 1594 for (i=0; i<m; i++) { 1595 y_l = m - i - 1.0; 1596 y_r = y_l + 1.0; 1597 color = PetscDrawRealToColor(PetscAbsScalar(v[j*m+i]),minv,maxv); 1598 PetscCall(PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color)); 1599 } 1600 } 1601 PetscDrawCollectiveEnd(draw); 1602 } 1603 PetscCall(MatDenseRestoreArrayRead(A,&v)); 1604 PetscFunctionReturn(0); 1605 } 1606 1607 static PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 1608 { 1609 PetscDraw draw; 1610 PetscBool isnull; 1611 PetscReal xr,yr,xl,yl,h,w; 1612 1613 PetscFunctionBegin; 1614 PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw)); 1615 PetscCall(PetscDrawIsNull(draw,&isnull)); 1616 if (isnull) PetscFunctionReturn(0); 1617 1618 xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0; 1619 xr += w; yr += h; xl = -w; yl = -h; 1620 PetscCall(PetscDrawSetCoordinates(draw,xl,yl,xr,yr)); 1621 PetscCall(PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer)); 1622 PetscCall(PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A)); 1623 PetscCall(PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL)); 1624 PetscCall(PetscDrawSave(draw)); 1625 PetscFunctionReturn(0); 1626 } 1627 1628 PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 1629 { 1630 PetscBool iascii,isbinary,isdraw; 1631 1632 PetscFunctionBegin; 1633 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 1634 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary)); 1635 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw)); 1636 if (iascii) { 1637 PetscCall(MatView_SeqDense_ASCII(A,viewer)); 1638 } else if (isbinary) { 1639 PetscCall(MatView_Dense_Binary(A,viewer)); 1640 } else if (isdraw) { 1641 PetscCall(MatView_SeqDense_Draw(A,viewer)); 1642 } 1643 PetscFunctionReturn(0); 1644 } 1645 1646 static PetscErrorCode MatDensePlaceArray_SeqDense(Mat A,const PetscScalar *array) 1647 { 1648 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1649 1650 PetscFunctionBegin; 1651 PetscCheck(!a->vecinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first"); 1652 PetscCheck(!a->matinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first"); 1653 PetscCheck(!a->unplacedarray,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreArray() first"); 1654 a->unplacedarray = a->v; 1655 a->unplaced_user_alloc = a->user_alloc; 1656 a->v = (PetscScalar*) array; 1657 a->user_alloc = PETSC_TRUE; 1658 #if defined(PETSC_HAVE_CUDA) 1659 A->offloadmask = PETSC_OFFLOAD_CPU; 1660 #endif 1661 PetscFunctionReturn(0); 1662 } 1663 1664 static PetscErrorCode MatDenseResetArray_SeqDense(Mat A) 1665 { 1666 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1667 1668 PetscFunctionBegin; 1669 PetscCheck(!a->vecinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first"); 1670 PetscCheck(!a->matinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first"); 1671 a->v = a->unplacedarray; 1672 a->user_alloc = a->unplaced_user_alloc; 1673 a->unplacedarray = NULL; 1674 #if defined(PETSC_HAVE_CUDA) 1675 A->offloadmask = PETSC_OFFLOAD_CPU; 1676 #endif 1677 PetscFunctionReturn(0); 1678 } 1679 1680 static PetscErrorCode MatDenseReplaceArray_SeqDense(Mat A,const PetscScalar *array) 1681 { 1682 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1683 1684 PetscFunctionBegin; 1685 PetscCheck(!a->vecinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first"); 1686 PetscCheck(!a->matinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first"); 1687 if (!a->user_alloc) PetscCall(PetscFree(a->v)); 1688 a->v = (PetscScalar*) array; 1689 a->user_alloc = PETSC_FALSE; 1690 #if defined(PETSC_HAVE_CUDA) 1691 A->offloadmask = PETSC_OFFLOAD_CPU; 1692 #endif 1693 PetscFunctionReturn(0); 1694 } 1695 1696 PetscErrorCode MatDestroy_SeqDense(Mat mat) 1697 { 1698 Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1699 1700 PetscFunctionBegin; 1701 #if defined(PETSC_USE_LOG) 1702 PetscLogObjectState((PetscObject)mat,"Rows %" PetscInt_FMT " Cols %" PetscInt_FMT,mat->rmap->n,mat->cmap->n); 1703 #endif 1704 PetscCall(VecDestroy(&(l->qrrhs))); 1705 PetscCall(PetscFree(l->tau)); 1706 PetscCall(PetscFree(l->pivots)); 1707 PetscCall(PetscFree(l->fwork)); 1708 PetscCall(MatDestroy(&l->ptapwork)); 1709 if (!l->user_alloc) PetscCall(PetscFree(l->v)); 1710 if (!l->unplaced_user_alloc) PetscCall(PetscFree(l->unplacedarray)); 1711 PetscCheck(!l->vecinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first"); 1712 PetscCheck(!l->matinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first"); 1713 PetscCall(VecDestroy(&l->cvec)); 1714 PetscCall(MatDestroy(&l->cmat)); 1715 PetscCall(PetscFree(mat->data)); 1716 1717 PetscCall(PetscObjectChangeTypeName((PetscObject)mat,NULL)); 1718 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatQRFactor_C",NULL)); 1719 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetLDA_C",NULL)); 1720 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseSetLDA_C",NULL)); 1721 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C",NULL)); 1722 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C",NULL)); 1723 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDensePlaceArray_C",NULL)); 1724 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseResetArray_C",NULL)); 1725 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseReplaceArray_C",NULL)); 1726 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayRead_C",NULL)); 1727 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayRead_C",NULL)); 1728 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArrayWrite_C",NULL)); 1729 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArrayWrite_C",NULL)); 1730 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_seqaij_C",NULL)); 1731 #if defined(PETSC_HAVE_ELEMENTAL) 1732 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_elemental_C",NULL)); 1733 #endif 1734 #if defined(PETSC_HAVE_SCALAPACK) 1735 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_scalapack_C",NULL)); 1736 #endif 1737 #if defined(PETSC_HAVE_CUDA) 1738 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatConvert_seqdense_seqdensecuda_C",NULL)); 1739 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqdensecuda_seqdensecuda_C",NULL)); 1740 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqdensecuda_seqdense_C",NULL)); 1741 #endif 1742 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatSeqDenseSetPreallocation_C",NULL)); 1743 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqaij_seqdense_C",NULL)); 1744 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqdense_seqdense_C",NULL)); 1745 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqbaij_seqdense_C",NULL)); 1746 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatProductSetFromOptions_seqsbaij_seqdense_C",NULL)); 1747 1748 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumn_C",NULL)); 1749 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumn_C",NULL)); 1750 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVec_C",NULL)); 1751 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVec_C",NULL)); 1752 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVecRead_C",NULL)); 1753 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVecRead_C",NULL)); 1754 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetColumnVecWrite_C",NULL)); 1755 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreColumnVecWrite_C",NULL)); 1756 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetSubMatrix_C",NULL)); 1757 PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreSubMatrix_C",NULL)); 1758 PetscFunctionReturn(0); 1759 } 1760 1761 static PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout) 1762 { 1763 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1764 PetscInt k,j,m = A->rmap->n, M = mat->lda, n = A->cmap->n; 1765 PetscScalar *v,tmp; 1766 1767 PetscFunctionBegin; 1768 if (reuse == MAT_INPLACE_MATRIX) { 1769 if (m == n) { /* in place transpose */ 1770 PetscCall(MatDenseGetArray(A,&v)); 1771 for (j=0; j<m; j++) { 1772 for (k=0; k<j; k++) { 1773 tmp = v[j + k*M]; 1774 v[j + k*M] = v[k + j*M]; 1775 v[k + j*M] = tmp; 1776 } 1777 } 1778 PetscCall(MatDenseRestoreArray(A,&v)); 1779 } else { /* reuse memory, temporary allocates new memory */ 1780 PetscScalar *v2; 1781 PetscLayout tmplayout; 1782 1783 PetscCall(PetscMalloc1((size_t)m*n,&v2)); 1784 PetscCall(MatDenseGetArray(A,&v)); 1785 for (j=0; j<n; j++) { 1786 for (k=0; k<m; k++) v2[j + (size_t)k*n] = v[k + (size_t)j*M]; 1787 } 1788 PetscCall(PetscArraycpy(v,v2,(size_t)m*n)); 1789 PetscCall(PetscFree(v2)); 1790 PetscCall(MatDenseRestoreArray(A,&v)); 1791 /* cleanup size dependent quantities */ 1792 PetscCall(VecDestroy(&mat->cvec)); 1793 PetscCall(MatDestroy(&mat->cmat)); 1794 PetscCall(PetscFree(mat->pivots)); 1795 PetscCall(PetscFree(mat->fwork)); 1796 PetscCall(MatDestroy(&mat->ptapwork)); 1797 /* swap row/col layouts */ 1798 mat->lda = n; 1799 tmplayout = A->rmap; 1800 A->rmap = A->cmap; 1801 A->cmap = tmplayout; 1802 } 1803 } else { /* out-of-place transpose */ 1804 Mat tmat; 1805 Mat_SeqDense *tmatd; 1806 PetscScalar *v2; 1807 PetscInt M2; 1808 1809 if (reuse == MAT_INITIAL_MATRIX) { 1810 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&tmat)); 1811 PetscCall(MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n)); 1812 PetscCall(MatSetType(tmat,((PetscObject)A)->type_name)); 1813 PetscCall(MatSeqDenseSetPreallocation(tmat,NULL)); 1814 } else tmat = *matout; 1815 1816 PetscCall(MatDenseGetArrayRead(A,(const PetscScalar**)&v)); 1817 PetscCall(MatDenseGetArray(tmat,&v2)); 1818 tmatd = (Mat_SeqDense*)tmat->data; 1819 M2 = tmatd->lda; 1820 for (j=0; j<n; j++) { 1821 for (k=0; k<m; k++) v2[j + k*M2] = v[k + j*M]; 1822 } 1823 PetscCall(MatDenseRestoreArray(tmat,&v2)); 1824 PetscCall(MatDenseRestoreArrayRead(A,(const PetscScalar**)&v)); 1825 PetscCall(MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY)); 1826 PetscCall(MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY)); 1827 *matout = tmat; 1828 } 1829 PetscFunctionReturn(0); 1830 } 1831 1832 static PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool *flg) 1833 { 1834 Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1835 Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 1836 PetscInt i; 1837 const PetscScalar *v1,*v2; 1838 1839 PetscFunctionBegin; 1840 if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1841 if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1842 PetscCall(MatDenseGetArrayRead(A1,&v1)); 1843 PetscCall(MatDenseGetArrayRead(A2,&v2)); 1844 for (i=0; i<A1->cmap->n; i++) { 1845 PetscCall(PetscArraycmp(v1,v2,A1->rmap->n,flg)); 1846 if (*flg == PETSC_FALSE) PetscFunctionReturn(0); 1847 v1 += mat1->lda; 1848 v2 += mat2->lda; 1849 } 1850 PetscCall(MatDenseRestoreArrayRead(A1,&v1)); 1851 PetscCall(MatDenseRestoreArrayRead(A2,&v2)); 1852 *flg = PETSC_TRUE; 1853 PetscFunctionReturn(0); 1854 } 1855 1856 static PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1857 { 1858 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1859 PetscInt i,n,len; 1860 PetscScalar *x; 1861 const PetscScalar *vv; 1862 1863 PetscFunctionBegin; 1864 PetscCall(VecGetSize(v,&n)); 1865 PetscCall(VecGetArray(v,&x)); 1866 len = PetscMin(A->rmap->n,A->cmap->n); 1867 PetscCall(MatDenseGetArrayRead(A,&vv)); 1868 PetscCheck(n == A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 1869 for (i=0; i<len; i++) { 1870 x[i] = vv[i*mat->lda + i]; 1871 } 1872 PetscCall(MatDenseRestoreArrayRead(A,&vv)); 1873 PetscCall(VecRestoreArray(v,&x)); 1874 PetscFunctionReturn(0); 1875 } 1876 1877 static PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1878 { 1879 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1880 const PetscScalar *l,*r; 1881 PetscScalar x,*v,*vv; 1882 PetscInt i,j,m = A->rmap->n,n = A->cmap->n; 1883 1884 PetscFunctionBegin; 1885 PetscCall(MatDenseGetArray(A,&vv)); 1886 if (ll) { 1887 PetscCall(VecGetSize(ll,&m)); 1888 PetscCall(VecGetArrayRead(ll,&l)); 1889 PetscCheck(m == A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1890 for (i=0; i<m; i++) { 1891 x = l[i]; 1892 v = vv + i; 1893 for (j=0; j<n; j++) { (*v) *= x; v+= mat->lda;} 1894 } 1895 PetscCall(VecRestoreArrayRead(ll,&l)); 1896 PetscCall(PetscLogFlops(1.0*n*m)); 1897 } 1898 if (rr) { 1899 PetscCall(VecGetSize(rr,&n)); 1900 PetscCall(VecGetArrayRead(rr,&r)); 1901 PetscCheck(n == A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1902 for (i=0; i<n; i++) { 1903 x = r[i]; 1904 v = vv + i*mat->lda; 1905 for (j=0; j<m; j++) (*v++) *= x; 1906 } 1907 PetscCall(VecRestoreArrayRead(rr,&r)); 1908 PetscCall(PetscLogFlops(1.0*n*m)); 1909 } 1910 PetscCall(MatDenseRestoreArray(A,&vv)); 1911 PetscFunctionReturn(0); 1912 } 1913 1914 PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1915 { 1916 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1917 PetscScalar *v,*vv; 1918 PetscReal sum = 0.0; 1919 PetscInt lda, m=A->rmap->n,i,j; 1920 1921 PetscFunctionBegin; 1922 PetscCall(MatDenseGetArrayRead(A,(const PetscScalar**)&vv)); 1923 PetscCall(MatDenseGetLDA(A,&lda)); 1924 v = vv; 1925 if (type == NORM_FROBENIUS) { 1926 if (lda>m) { 1927 for (j=0; j<A->cmap->n; j++) { 1928 v = vv+j*lda; 1929 for (i=0; i<m; i++) { 1930 sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1931 } 1932 } 1933 } else { 1934 #if defined(PETSC_USE_REAL___FP16) 1935 PetscBLASInt one = 1,cnt = A->cmap->n*A->rmap->n; 1936 PetscStackCallBLAS("BLASnrm2",*nrm = BLASnrm2_(&cnt,v,&one)); 1937 } 1938 #else 1939 for (i=0; i<A->cmap->n*A->rmap->n; i++) { 1940 sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1941 } 1942 } 1943 *nrm = PetscSqrtReal(sum); 1944 #endif 1945 PetscCall(PetscLogFlops(2.0*A->cmap->n*A->rmap->n)); 1946 } else if (type == NORM_1) { 1947 *nrm = 0.0; 1948 for (j=0; j<A->cmap->n; j++) { 1949 v = vv + j*mat->lda; 1950 sum = 0.0; 1951 for (i=0; i<A->rmap->n; i++) { 1952 sum += PetscAbsScalar(*v); v++; 1953 } 1954 if (sum > *nrm) *nrm = sum; 1955 } 1956 PetscCall(PetscLogFlops(1.0*A->cmap->n*A->rmap->n)); 1957 } else if (type == NORM_INFINITY) { 1958 *nrm = 0.0; 1959 for (j=0; j<A->rmap->n; j++) { 1960 v = vv + j; 1961 sum = 0.0; 1962 for (i=0; i<A->cmap->n; i++) { 1963 sum += PetscAbsScalar(*v); v += mat->lda; 1964 } 1965 if (sum > *nrm) *nrm = sum; 1966 } 1967 PetscCall(PetscLogFlops(1.0*A->cmap->n*A->rmap->n)); 1968 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm"); 1969 PetscCall(MatDenseRestoreArrayRead(A,(const PetscScalar**)&vv)); 1970 PetscFunctionReturn(0); 1971 } 1972 1973 static PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg) 1974 { 1975 Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 1976 1977 PetscFunctionBegin; 1978 switch (op) { 1979 case MAT_ROW_ORIENTED: 1980 aij->roworiented = flg; 1981 break; 1982 case MAT_NEW_NONZERO_LOCATIONS: 1983 case MAT_NEW_NONZERO_LOCATION_ERR: 1984 case MAT_NEW_NONZERO_ALLOCATION_ERR: 1985 case MAT_FORCE_DIAGONAL_ENTRIES: 1986 case MAT_KEEP_NONZERO_PATTERN: 1987 case MAT_IGNORE_OFF_PROC_ENTRIES: 1988 case MAT_USE_HASH_TABLE: 1989 case MAT_IGNORE_ZERO_ENTRIES: 1990 case MAT_IGNORE_LOWER_TRIANGULAR: 1991 case MAT_SORTED_FULL: 1992 PetscCall(PetscInfo(A,"Option %s ignored\n",MatOptions[op])); 1993 break; 1994 case MAT_SPD: 1995 case MAT_SYMMETRIC: 1996 case MAT_STRUCTURALLY_SYMMETRIC: 1997 case MAT_HERMITIAN: 1998 case MAT_SYMMETRY_ETERNAL: 1999 /* These options are handled directly by MatSetOption() */ 2000 break; 2001 default: 2002 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 2003 } 2004 PetscFunctionReturn(0); 2005 } 2006 2007 PetscErrorCode MatZeroEntries_SeqDense(Mat A) 2008 { 2009 Mat_SeqDense *l = (Mat_SeqDense*)A->data; 2010 PetscInt lda=l->lda,m=A->rmap->n,n=A->cmap->n,j; 2011 PetscScalar *v; 2012 2013 PetscFunctionBegin; 2014 PetscCall(MatDenseGetArrayWrite(A,&v)); 2015 if (lda>m) { 2016 for (j=0; j<n; j++) { 2017 PetscCall(PetscArrayzero(v+j*lda,m)); 2018 } 2019 } else { 2020 PetscCall(PetscArrayzero(v,PetscInt64Mult(m,n))); 2021 } 2022 PetscCall(MatDenseRestoreArrayWrite(A,&v)); 2023 PetscFunctionReturn(0); 2024 } 2025 2026 static PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 2027 { 2028 Mat_SeqDense *l = (Mat_SeqDense*)A->data; 2029 PetscInt m = l->lda, n = A->cmap->n, i,j; 2030 PetscScalar *slot,*bb,*v; 2031 const PetscScalar *xx; 2032 2033 PetscFunctionBegin; 2034 if (PetscDefined(USE_DEBUG)) { 2035 for (i=0; i<N; i++) { 2036 PetscCheck(rows[i] >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed"); 2037 PetscCheck(rows[i] < A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %" PetscInt_FMT " requested to be zeroed greater than or equal number of rows %" PetscInt_FMT,rows[i],A->rmap->n); 2038 } 2039 } 2040 if (!N) PetscFunctionReturn(0); 2041 2042 /* fix right hand side if needed */ 2043 if (x && b) { 2044 PetscCall(VecGetArrayRead(x,&xx)); 2045 PetscCall(VecGetArray(b,&bb)); 2046 for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]]; 2047 PetscCall(VecRestoreArrayRead(x,&xx)); 2048 PetscCall(VecRestoreArray(b,&bb)); 2049 } 2050 2051 PetscCall(MatDenseGetArray(A,&v)); 2052 for (i=0; i<N; i++) { 2053 slot = v + rows[i]; 2054 for (j=0; j<n; j++) { *slot = 0.0; slot += m;} 2055 } 2056 if (diag != 0.0) { 2057 PetscCheck(A->rmap->n == A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices"); 2058 for (i=0; i<N; i++) { 2059 slot = v + (m+1)*rows[i]; 2060 *slot = diag; 2061 } 2062 } 2063 PetscCall(MatDenseRestoreArray(A,&v)); 2064 PetscFunctionReturn(0); 2065 } 2066 2067 static PetscErrorCode MatDenseGetLDA_SeqDense(Mat A,PetscInt *lda) 2068 { 2069 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2070 2071 PetscFunctionBegin; 2072 *lda = mat->lda; 2073 PetscFunctionReturn(0); 2074 } 2075 2076 PetscErrorCode MatDenseGetArray_SeqDense(Mat A,PetscScalar **array) 2077 { 2078 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2079 2080 PetscFunctionBegin; 2081 PetscCheck(!mat->matinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first"); 2082 *array = mat->v; 2083 PetscFunctionReturn(0); 2084 } 2085 2086 PetscErrorCode MatDenseRestoreArray_SeqDense(Mat A,PetscScalar **array) 2087 { 2088 PetscFunctionBegin; 2089 if (array) *array = NULL; 2090 PetscFunctionReturn(0); 2091 } 2092 2093 /*@ 2094 MatDenseGetLDA - gets the leading dimension of the array returned from MatDenseGetArray() 2095 2096 Not collective 2097 2098 Input Parameter: 2099 . mat - a MATSEQDENSE or MATMPIDENSE matrix 2100 2101 Output Parameter: 2102 . lda - the leading dimension 2103 2104 Level: intermediate 2105 2106 .seealso: `MatDenseGetArray()`, `MatDenseRestoreArray()`, `MatDenseGetArrayRead()`, `MatDenseRestoreArrayRead()`, `MatDenseSetLDA()` 2107 @*/ 2108 PetscErrorCode MatDenseGetLDA(Mat A,PetscInt *lda) 2109 { 2110 PetscFunctionBegin; 2111 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2112 PetscValidIntPointer(lda,2); 2113 MatCheckPreallocated(A,1); 2114 PetscUseMethod(A,"MatDenseGetLDA_C",(Mat,PetscInt*),(A,lda)); 2115 PetscFunctionReturn(0); 2116 } 2117 2118 /*@ 2119 MatDenseSetLDA - Sets the leading dimension of the array used by the dense matrix 2120 2121 Not collective 2122 2123 Input Parameters: 2124 + mat - a MATSEQDENSE or MATMPIDENSE matrix 2125 - lda - the leading dimension 2126 2127 Level: intermediate 2128 2129 .seealso: `MatDenseGetArray()`, `MatDenseRestoreArray()`, `MatDenseGetArrayRead()`, `MatDenseRestoreArrayRead()`, `MatDenseGetLDA()` 2130 @*/ 2131 PetscErrorCode MatDenseSetLDA(Mat A,PetscInt lda) 2132 { 2133 PetscFunctionBegin; 2134 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2135 PetscTryMethod(A,"MatDenseSetLDA_C",(Mat,PetscInt),(A,lda)); 2136 PetscFunctionReturn(0); 2137 } 2138 2139 /*@C 2140 MatDenseGetArray - gives read-write access to the array where the data for a dense matrix is stored 2141 2142 Logically Collective on Mat 2143 2144 Input Parameter: 2145 . mat - a dense matrix 2146 2147 Output Parameter: 2148 . array - pointer to the data 2149 2150 Level: intermediate 2151 2152 .seealso: `MatDenseRestoreArray()`, `MatDenseGetArrayRead()`, `MatDenseRestoreArrayRead()`, `MatDenseGetArrayWrite()`, `MatDenseRestoreArrayWrite()` 2153 @*/ 2154 PetscErrorCode MatDenseGetArray(Mat A,PetscScalar **array) 2155 { 2156 PetscFunctionBegin; 2157 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2158 PetscValidPointer(array,2); 2159 PetscUseMethod(A,"MatDenseGetArray_C",(Mat,PetscScalar**),(A,array)); 2160 PetscFunctionReturn(0); 2161 } 2162 2163 /*@C 2164 MatDenseRestoreArray - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArray() 2165 2166 Logically Collective on Mat 2167 2168 Input Parameters: 2169 + mat - a dense matrix 2170 - array - pointer to the data 2171 2172 Level: intermediate 2173 2174 .seealso: `MatDenseGetArray()`, `MatDenseGetArrayRead()`, `MatDenseRestoreArrayRead()`, `MatDenseGetArrayWrite()`, `MatDenseRestoreArrayWrite()` 2175 @*/ 2176 PetscErrorCode MatDenseRestoreArray(Mat A,PetscScalar **array) 2177 { 2178 PetscFunctionBegin; 2179 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2180 PetscValidPointer(array,2); 2181 PetscUseMethod(A,"MatDenseRestoreArray_C",(Mat,PetscScalar**),(A,array)); 2182 PetscCall(PetscObjectStateIncrease((PetscObject)A)); 2183 #if defined(PETSC_HAVE_CUDA) 2184 A->offloadmask = PETSC_OFFLOAD_CPU; 2185 #endif 2186 PetscFunctionReturn(0); 2187 } 2188 2189 /*@C 2190 MatDenseGetArrayRead - gives read-only access to the array where the data for a dense matrix is stored 2191 2192 Not Collective 2193 2194 Input Parameter: 2195 . mat - a dense matrix 2196 2197 Output Parameter: 2198 . array - pointer to the data 2199 2200 Level: intermediate 2201 2202 .seealso: `MatDenseRestoreArrayRead()`, `MatDenseGetArray()`, `MatDenseRestoreArray()`, `MatDenseGetArrayWrite()`, `MatDenseRestoreArrayWrite()` 2203 @*/ 2204 PetscErrorCode MatDenseGetArrayRead(Mat A,const PetscScalar **array) 2205 { 2206 PetscFunctionBegin; 2207 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2208 PetscValidPointer(array,2); 2209 PetscUseMethod(A,"MatDenseGetArrayRead_C",(Mat,const PetscScalar**),(A,array)); 2210 PetscFunctionReturn(0); 2211 } 2212 2213 /*@C 2214 MatDenseRestoreArrayRead - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArrayRead() 2215 2216 Not Collective 2217 2218 Input Parameters: 2219 + mat - a dense matrix 2220 - array - pointer to the data 2221 2222 Level: intermediate 2223 2224 .seealso: `MatDenseGetArrayRead()`, `MatDenseGetArray()`, `MatDenseRestoreArray()`, `MatDenseGetArrayWrite()`, `MatDenseRestoreArrayWrite()` 2225 @*/ 2226 PetscErrorCode MatDenseRestoreArrayRead(Mat A,const PetscScalar **array) 2227 { 2228 PetscFunctionBegin; 2229 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2230 PetscValidPointer(array,2); 2231 PetscUseMethod(A,"MatDenseRestoreArrayRead_C",(Mat,const PetscScalar**),(A,array)); 2232 PetscFunctionReturn(0); 2233 } 2234 2235 /*@C 2236 MatDenseGetArrayWrite - gives write-only access to the array where the data for a dense matrix is stored 2237 2238 Not Collective 2239 2240 Input Parameter: 2241 . mat - a dense matrix 2242 2243 Output Parameter: 2244 . array - pointer to the data 2245 2246 Level: intermediate 2247 2248 .seealso: `MatDenseRestoreArrayWrite()`, `MatDenseGetArray()`, `MatDenseRestoreArray()`, `MatDenseGetArrayRead()`, `MatDenseRestoreArrayRead()` 2249 @*/ 2250 PetscErrorCode MatDenseGetArrayWrite(Mat A,PetscScalar **array) 2251 { 2252 PetscFunctionBegin; 2253 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2254 PetscValidPointer(array,2); 2255 PetscUseMethod(A,"MatDenseGetArrayWrite_C",(Mat,PetscScalar**),(A,array)); 2256 PetscFunctionReturn(0); 2257 } 2258 2259 /*@C 2260 MatDenseRestoreArrayWrite - returns access to the array where the data for a dense matrix is stored obtained by MatDenseGetArrayWrite() 2261 2262 Not Collective 2263 2264 Input Parameters: 2265 + mat - a dense matrix 2266 - array - pointer to the data 2267 2268 Level: intermediate 2269 2270 .seealso: `MatDenseGetArrayWrite()`, `MatDenseGetArray()`, `MatDenseRestoreArray()`, `MatDenseGetArrayRead()`, `MatDenseRestoreArrayRead()` 2271 @*/ 2272 PetscErrorCode MatDenseRestoreArrayWrite(Mat A,PetscScalar **array) 2273 { 2274 PetscFunctionBegin; 2275 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2276 PetscValidPointer(array,2); 2277 PetscUseMethod(A,"MatDenseRestoreArrayWrite_C",(Mat,PetscScalar**),(A,array)); 2278 PetscCall(PetscObjectStateIncrease((PetscObject)A)); 2279 #if defined(PETSC_HAVE_CUDA) 2280 A->offloadmask = PETSC_OFFLOAD_CPU; 2281 #endif 2282 PetscFunctionReturn(0); 2283 } 2284 2285 static PetscErrorCode MatCreateSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatReuse scall,Mat *B) 2286 { 2287 Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2288 PetscInt i,j,nrows,ncols,ldb; 2289 const PetscInt *irow,*icol; 2290 PetscScalar *av,*bv,*v = mat->v; 2291 Mat newmat; 2292 2293 PetscFunctionBegin; 2294 PetscCall(ISGetIndices(isrow,&irow)); 2295 PetscCall(ISGetIndices(iscol,&icol)); 2296 PetscCall(ISGetLocalSize(isrow,&nrows)); 2297 PetscCall(ISGetLocalSize(iscol,&ncols)); 2298 2299 /* Check submatrixcall */ 2300 if (scall == MAT_REUSE_MATRIX) { 2301 PetscInt n_cols,n_rows; 2302 PetscCall(MatGetSize(*B,&n_rows,&n_cols)); 2303 if (n_rows != nrows || n_cols != ncols) { 2304 /* resize the result matrix to match number of requested rows/columns */ 2305 PetscCall(MatSetSizes(*B,nrows,ncols,nrows,ncols)); 2306 } 2307 newmat = *B; 2308 } else { 2309 /* Create and fill new matrix */ 2310 PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&newmat)); 2311 PetscCall(MatSetSizes(newmat,nrows,ncols,nrows,ncols)); 2312 PetscCall(MatSetType(newmat,((PetscObject)A)->type_name)); 2313 PetscCall(MatSeqDenseSetPreallocation(newmat,NULL)); 2314 } 2315 2316 /* Now extract the data pointers and do the copy,column at a time */ 2317 PetscCall(MatDenseGetArray(newmat,&bv)); 2318 PetscCall(MatDenseGetLDA(newmat,&ldb)); 2319 for (i=0; i<ncols; i++) { 2320 av = v + mat->lda*icol[i]; 2321 for (j=0; j<nrows; j++) bv[j] = av[irow[j]]; 2322 bv += ldb; 2323 } 2324 PetscCall(MatDenseRestoreArray(newmat,&bv)); 2325 2326 /* Assemble the matrices so that the correct flags are set */ 2327 PetscCall(MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY)); 2328 PetscCall(MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY)); 2329 2330 /* Free work space */ 2331 PetscCall(ISRestoreIndices(isrow,&irow)); 2332 PetscCall(ISRestoreIndices(iscol,&icol)); 2333 *B = newmat; 2334 PetscFunctionReturn(0); 2335 } 2336 2337 static PetscErrorCode MatCreateSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 2338 { 2339 PetscInt i; 2340 2341 PetscFunctionBegin; 2342 if (scall == MAT_INITIAL_MATRIX) { 2343 PetscCall(PetscCalloc1(n,B)); 2344 } 2345 2346 for (i=0; i<n; i++) { 2347 PetscCall(MatCreateSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i])); 2348 } 2349 PetscFunctionReturn(0); 2350 } 2351 2352 static PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 2353 { 2354 PetscFunctionBegin; 2355 PetscFunctionReturn(0); 2356 } 2357 2358 static PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 2359 { 2360 PetscFunctionBegin; 2361 PetscFunctionReturn(0); 2362 } 2363 2364 PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 2365 { 2366 Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense*)B->data; 2367 const PetscScalar *va; 2368 PetscScalar *vb; 2369 PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j; 2370 2371 PetscFunctionBegin; 2372 /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 2373 if (A->ops->copy != B->ops->copy) { 2374 PetscCall(MatCopy_Basic(A,B,str)); 2375 PetscFunctionReturn(0); 2376 } 2377 PetscCheck(m == B->rmap->n && n == B->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 2378 PetscCall(MatDenseGetArrayRead(A,&va)); 2379 PetscCall(MatDenseGetArray(B,&vb)); 2380 if (lda1>m || lda2>m) { 2381 for (j=0; j<n; j++) { 2382 PetscCall(PetscArraycpy(vb+j*lda2,va+j*lda1,m)); 2383 } 2384 } else { 2385 PetscCall(PetscArraycpy(vb,va,A->rmap->n*A->cmap->n)); 2386 } 2387 PetscCall(MatDenseRestoreArray(B,&vb)); 2388 PetscCall(MatDenseRestoreArrayRead(A,&va)); 2389 PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 2390 PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 2391 PetscFunctionReturn(0); 2392 } 2393 2394 PetscErrorCode MatSetUp_SeqDense(Mat A) 2395 { 2396 PetscFunctionBegin; 2397 PetscCall(PetscLayoutSetUp(A->rmap)); 2398 PetscCall(PetscLayoutSetUp(A->cmap)); 2399 if (!A->preallocated) { 2400 PetscCall(MatSeqDenseSetPreallocation(A,NULL)); 2401 } 2402 PetscFunctionReturn(0); 2403 } 2404 2405 static PetscErrorCode MatConjugate_SeqDense(Mat A) 2406 { 2407 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 2408 PetscInt i,nz = A->rmap->n*A->cmap->n; 2409 PetscInt min = PetscMin(A->rmap->n,A->cmap->n); 2410 PetscScalar *aa; 2411 2412 PetscFunctionBegin; 2413 PetscCall(MatDenseGetArray(A,&aa)); 2414 for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]); 2415 PetscCall(MatDenseRestoreArray(A,&aa)); 2416 if (mat->tau) for (i = 0; i < min; i++) mat->tau[i] = PetscConj(mat->tau[i]); 2417 PetscFunctionReturn(0); 2418 } 2419 2420 static PetscErrorCode MatRealPart_SeqDense(Mat A) 2421 { 2422 PetscInt i,nz = A->rmap->n*A->cmap->n; 2423 PetscScalar *aa; 2424 2425 PetscFunctionBegin; 2426 PetscCall(MatDenseGetArray(A,&aa)); 2427 for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 2428 PetscCall(MatDenseRestoreArray(A,&aa)); 2429 PetscFunctionReturn(0); 2430 } 2431 2432 static PetscErrorCode MatImaginaryPart_SeqDense(Mat A) 2433 { 2434 PetscInt i,nz = A->rmap->n*A->cmap->n; 2435 PetscScalar *aa; 2436 2437 PetscFunctionBegin; 2438 PetscCall(MatDenseGetArray(A,&aa)); 2439 for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 2440 PetscCall(MatDenseRestoreArray(A,&aa)); 2441 PetscFunctionReturn(0); 2442 } 2443 2444 /* ----------------------------------------------------------------*/ 2445 PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat C) 2446 { 2447 PetscInt m=A->rmap->n,n=B->cmap->n; 2448 PetscBool cisdense; 2449 2450 PetscFunctionBegin; 2451 PetscCall(MatSetSizes(C,m,n,m,n)); 2452 PetscCall(PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATSEQDENSE,MATSEQDENSECUDA,"")); 2453 if (!cisdense) { 2454 PetscBool flg; 2455 2456 PetscCall(PetscObjectTypeCompare((PetscObject)B,((PetscObject)A)->type_name,&flg)); 2457 PetscCall(MatSetType(C,flg ? ((PetscObject)A)->type_name : MATDENSE)); 2458 } 2459 PetscCall(MatSetUp(C)); 2460 PetscFunctionReturn(0); 2461 } 2462 2463 PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 2464 { 2465 Mat_SeqDense *a=(Mat_SeqDense*)A->data,*b=(Mat_SeqDense*)B->data,*c=(Mat_SeqDense*)C->data; 2466 PetscBLASInt m,n,k; 2467 const PetscScalar *av,*bv; 2468 PetscScalar *cv; 2469 PetscScalar _DOne=1.0,_DZero=0.0; 2470 2471 PetscFunctionBegin; 2472 PetscCall(PetscBLASIntCast(C->rmap->n,&m)); 2473 PetscCall(PetscBLASIntCast(C->cmap->n,&n)); 2474 PetscCall(PetscBLASIntCast(A->cmap->n,&k)); 2475 if (!m || !n || !k) PetscFunctionReturn(0); 2476 PetscCall(MatDenseGetArrayRead(A,&av)); 2477 PetscCall(MatDenseGetArrayRead(B,&bv)); 2478 PetscCall(MatDenseGetArrayWrite(C,&cv)); 2479 PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&n,&k,&_DOne,av,&a->lda,bv,&b->lda,&_DZero,cv,&c->lda)); 2480 PetscCall(PetscLogFlops(1.0*m*n*k + 1.0*m*n*(k-1))); 2481 PetscCall(MatDenseRestoreArrayRead(A,&av)); 2482 PetscCall(MatDenseRestoreArrayRead(B,&bv)); 2483 PetscCall(MatDenseRestoreArrayWrite(C,&cv)); 2484 PetscFunctionReturn(0); 2485 } 2486 2487 PetscErrorCode MatMatTransposeMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat C) 2488 { 2489 PetscInt m=A->rmap->n,n=B->rmap->n; 2490 PetscBool cisdense; 2491 2492 PetscFunctionBegin; 2493 PetscCall(MatSetSizes(C,m,n,m,n)); 2494 PetscCall(PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATSEQDENSE,MATSEQDENSECUDA,"")); 2495 if (!cisdense) { 2496 PetscBool flg; 2497 2498 PetscCall(PetscObjectTypeCompare((PetscObject)B,((PetscObject)A)->type_name,&flg)); 2499 PetscCall(MatSetType(C,flg ? ((PetscObject)A)->type_name : MATDENSE)); 2500 } 2501 PetscCall(MatSetUp(C)); 2502 PetscFunctionReturn(0); 2503 } 2504 2505 PetscErrorCode MatMatTransposeMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 2506 { 2507 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 2508 Mat_SeqDense *b = (Mat_SeqDense*)B->data; 2509 Mat_SeqDense *c = (Mat_SeqDense*)C->data; 2510 const PetscScalar *av,*bv; 2511 PetscScalar *cv; 2512 PetscBLASInt m,n,k; 2513 PetscScalar _DOne=1.0,_DZero=0.0; 2514 2515 PetscFunctionBegin; 2516 PetscCall(PetscBLASIntCast(C->rmap->n,&m)); 2517 PetscCall(PetscBLASIntCast(C->cmap->n,&n)); 2518 PetscCall(PetscBLASIntCast(A->cmap->n,&k)); 2519 if (!m || !n || !k) PetscFunctionReturn(0); 2520 PetscCall(MatDenseGetArrayRead(A,&av)); 2521 PetscCall(MatDenseGetArrayRead(B,&bv)); 2522 PetscCall(MatDenseGetArrayWrite(C,&cv)); 2523 PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&m,&n,&k,&_DOne,av,&a->lda,bv,&b->lda,&_DZero,cv,&c->lda)); 2524 PetscCall(MatDenseRestoreArrayRead(A,&av)); 2525 PetscCall(MatDenseRestoreArrayRead(B,&bv)); 2526 PetscCall(MatDenseRestoreArrayWrite(C,&cv)); 2527 PetscCall(PetscLogFlops(1.0*m*n*k + 1.0*m*n*(k-1))); 2528 PetscFunctionReturn(0); 2529 } 2530 2531 PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat C) 2532 { 2533 PetscInt m=A->cmap->n,n=B->cmap->n; 2534 PetscBool cisdense; 2535 2536 PetscFunctionBegin; 2537 PetscCall(MatSetSizes(C,m,n,m,n)); 2538 PetscCall(PetscObjectTypeCompareAny((PetscObject)C,&cisdense,MATSEQDENSE,MATSEQDENSECUDA,"")); 2539 if (!cisdense) { 2540 PetscBool flg; 2541 2542 PetscCall(PetscObjectTypeCompare((PetscObject)B,((PetscObject)A)->type_name,&flg)); 2543 PetscCall(MatSetType(C,flg ? ((PetscObject)A)->type_name : MATDENSE)); 2544 } 2545 PetscCall(MatSetUp(C)); 2546 PetscFunctionReturn(0); 2547 } 2548 2549 PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 2550 { 2551 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 2552 Mat_SeqDense *b = (Mat_SeqDense*)B->data; 2553 Mat_SeqDense *c = (Mat_SeqDense*)C->data; 2554 const PetscScalar *av,*bv; 2555 PetscScalar *cv; 2556 PetscBLASInt m,n,k; 2557 PetscScalar _DOne=1.0,_DZero=0.0; 2558 2559 PetscFunctionBegin; 2560 PetscCall(PetscBLASIntCast(C->rmap->n,&m)); 2561 PetscCall(PetscBLASIntCast(C->cmap->n,&n)); 2562 PetscCall(PetscBLASIntCast(A->rmap->n,&k)); 2563 if (!m || !n || !k) PetscFunctionReturn(0); 2564 PetscCall(MatDenseGetArrayRead(A,&av)); 2565 PetscCall(MatDenseGetArrayRead(B,&bv)); 2566 PetscCall(MatDenseGetArrayWrite(C,&cv)); 2567 PetscStackCallBLAS("BLASgemm",BLASgemm_("T","N",&m,&n,&k,&_DOne,av,&a->lda,bv,&b->lda,&_DZero,cv,&c->lda)); 2568 PetscCall(MatDenseRestoreArrayRead(A,&av)); 2569 PetscCall(MatDenseRestoreArrayRead(B,&bv)); 2570 PetscCall(MatDenseRestoreArrayWrite(C,&cv)); 2571 PetscCall(PetscLogFlops(1.0*m*n*k + 1.0*m*n*(k-1))); 2572 PetscFunctionReturn(0); 2573 } 2574 2575 /* ----------------------------------------------- */ 2576 static PetscErrorCode MatProductSetFromOptions_SeqDense_AB(Mat C) 2577 { 2578 PetscFunctionBegin; 2579 C->ops->matmultsymbolic = MatMatMultSymbolic_SeqDense_SeqDense; 2580 C->ops->productsymbolic = MatProductSymbolic_AB; 2581 PetscFunctionReturn(0); 2582 } 2583 2584 static PetscErrorCode MatProductSetFromOptions_SeqDense_AtB(Mat C) 2585 { 2586 PetscFunctionBegin; 2587 C->ops->transposematmultsymbolic = MatTransposeMatMultSymbolic_SeqDense_SeqDense; 2588 C->ops->productsymbolic = MatProductSymbolic_AtB; 2589 PetscFunctionReturn(0); 2590 } 2591 2592 static PetscErrorCode MatProductSetFromOptions_SeqDense_ABt(Mat C) 2593 { 2594 PetscFunctionBegin; 2595 C->ops->mattransposemultsymbolic = MatMatTransposeMultSymbolic_SeqDense_SeqDense; 2596 C->ops->productsymbolic = MatProductSymbolic_ABt; 2597 PetscFunctionReturn(0); 2598 } 2599 2600 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense(Mat C) 2601 { 2602 Mat_Product *product = C->product; 2603 2604 PetscFunctionBegin; 2605 switch (product->type) { 2606 case MATPRODUCT_AB: 2607 PetscCall(MatProductSetFromOptions_SeqDense_AB(C)); 2608 break; 2609 case MATPRODUCT_AtB: 2610 PetscCall(MatProductSetFromOptions_SeqDense_AtB(C)); 2611 break; 2612 case MATPRODUCT_ABt: 2613 PetscCall(MatProductSetFromOptions_SeqDense_ABt(C)); 2614 break; 2615 default: 2616 break; 2617 } 2618 PetscFunctionReturn(0); 2619 } 2620 /* ----------------------------------------------- */ 2621 2622 static PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[]) 2623 { 2624 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 2625 PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 2626 PetscScalar *x; 2627 const PetscScalar *aa; 2628 2629 PetscFunctionBegin; 2630 PetscCheck(!A->factortype,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 2631 PetscCall(VecGetArray(v,&x)); 2632 PetscCall(VecGetLocalSize(v,&p)); 2633 PetscCall(MatDenseGetArrayRead(A,&aa)); 2634 PetscCheck(p == A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 2635 for (i=0; i<m; i++) { 2636 x[i] = aa[i]; if (idx) idx[i] = 0; 2637 for (j=1; j<n; j++) { 2638 if (PetscRealPart(x[i]) < PetscRealPart(aa[i+a->lda*j])) {x[i] = aa[i + a->lda*j]; if (idx) idx[i] = j;} 2639 } 2640 } 2641 PetscCall(MatDenseRestoreArrayRead(A,&aa)); 2642 PetscCall(VecRestoreArray(v,&x)); 2643 PetscFunctionReturn(0); 2644 } 2645 2646 static PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[]) 2647 { 2648 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 2649 PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 2650 PetscScalar *x; 2651 PetscReal atmp; 2652 const PetscScalar *aa; 2653 2654 PetscFunctionBegin; 2655 PetscCheck(!A->factortype,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 2656 PetscCall(VecGetArray(v,&x)); 2657 PetscCall(VecGetLocalSize(v,&p)); 2658 PetscCall(MatDenseGetArrayRead(A,&aa)); 2659 PetscCheck(p == A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 2660 for (i=0; i<m; i++) { 2661 x[i] = PetscAbsScalar(aa[i]); 2662 for (j=1; j<n; j++) { 2663 atmp = PetscAbsScalar(aa[i+a->lda*j]); 2664 if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;} 2665 } 2666 } 2667 PetscCall(MatDenseRestoreArrayRead(A,&aa)); 2668 PetscCall(VecRestoreArray(v,&x)); 2669 PetscFunctionReturn(0); 2670 } 2671 2672 static PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[]) 2673 { 2674 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 2675 PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 2676 PetscScalar *x; 2677 const PetscScalar *aa; 2678 2679 PetscFunctionBegin; 2680 PetscCheck(!A->factortype,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 2681 PetscCall(MatDenseGetArrayRead(A,&aa)); 2682 PetscCall(VecGetArray(v,&x)); 2683 PetscCall(VecGetLocalSize(v,&p)); 2684 PetscCheck(p == A->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 2685 for (i=0; i<m; i++) { 2686 x[i] = aa[i]; if (idx) idx[i] = 0; 2687 for (j=1; j<n; j++) { 2688 if (PetscRealPart(x[i]) > PetscRealPart(aa[i+a->lda*j])) {x[i] = aa[i + a->lda*j]; if (idx) idx[i] = j;} 2689 } 2690 } 2691 PetscCall(VecRestoreArray(v,&x)); 2692 PetscCall(MatDenseRestoreArrayRead(A,&aa)); 2693 PetscFunctionReturn(0); 2694 } 2695 2696 PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col) 2697 { 2698 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 2699 PetscScalar *x; 2700 const PetscScalar *aa; 2701 2702 PetscFunctionBegin; 2703 PetscCheck(!A->factortype,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 2704 PetscCall(MatDenseGetArrayRead(A,&aa)); 2705 PetscCall(VecGetArray(v,&x)); 2706 PetscCall(PetscArraycpy(x,aa+col*a->lda,A->rmap->n)); 2707 PetscCall(VecRestoreArray(v,&x)); 2708 PetscCall(MatDenseRestoreArrayRead(A,&aa)); 2709 PetscFunctionReturn(0); 2710 } 2711 2712 PETSC_INTERN PetscErrorCode MatGetColumnReductions_SeqDense(Mat A,PetscInt type,PetscReal *reductions) 2713 { 2714 PetscInt i,j,m,n; 2715 const PetscScalar *a; 2716 2717 PetscFunctionBegin; 2718 PetscCall(MatGetSize(A,&m,&n)); 2719 PetscCall(PetscArrayzero(reductions,n)); 2720 PetscCall(MatDenseGetArrayRead(A,&a)); 2721 if (type == NORM_2) { 2722 for (i=0; i<n; i++) { 2723 for (j=0; j<m; j++) { 2724 reductions[i] += PetscAbsScalar(a[j]*a[j]); 2725 } 2726 a += m; 2727 } 2728 } else if (type == NORM_1) { 2729 for (i=0; i<n; i++) { 2730 for (j=0; j<m; j++) { 2731 reductions[i] += PetscAbsScalar(a[j]); 2732 } 2733 a += m; 2734 } 2735 } else if (type == NORM_INFINITY) { 2736 for (i=0; i<n; i++) { 2737 for (j=0; j<m; j++) { 2738 reductions[i] = PetscMax(PetscAbsScalar(a[j]),reductions[i]); 2739 } 2740 a += m; 2741 } 2742 } else if (type == REDUCTION_SUM_REALPART || type == REDUCTION_MEAN_REALPART) { 2743 for (i=0; i<n; i++) { 2744 for (j=0; j<m; j++) { 2745 reductions[i] += PetscRealPart(a[j]); 2746 } 2747 a += m; 2748 } 2749 } else if (type == REDUCTION_SUM_IMAGINARYPART || type == REDUCTION_MEAN_IMAGINARYPART) { 2750 for (i=0; i<n; i++) { 2751 for (j=0; j<m; j++) { 2752 reductions[i] += PetscImaginaryPart(a[j]); 2753 } 2754 a += m; 2755 } 2756 } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown reduction type"); 2757 PetscCall(MatDenseRestoreArrayRead(A,&a)); 2758 if (type == NORM_2) { 2759 for (i=0; i<n; i++) reductions[i] = PetscSqrtReal(reductions[i]); 2760 } else if (type == REDUCTION_MEAN_REALPART || type == REDUCTION_MEAN_IMAGINARYPART) { 2761 for (i=0; i<n; i++) reductions[i] /= m; 2762 } 2763 PetscFunctionReturn(0); 2764 } 2765 2766 static PetscErrorCode MatSetRandom_SeqDense(Mat x,PetscRandom rctx) 2767 { 2768 PetscScalar *a; 2769 PetscInt lda,m,n,i,j; 2770 2771 PetscFunctionBegin; 2772 PetscCall(MatGetSize(x,&m,&n)); 2773 PetscCall(MatDenseGetLDA(x,&lda)); 2774 PetscCall(MatDenseGetArray(x,&a)); 2775 for (j=0; j<n; j++) { 2776 for (i=0; i<m; i++) { 2777 PetscCall(PetscRandomGetValue(rctx,a+j*lda+i)); 2778 } 2779 } 2780 PetscCall(MatDenseRestoreArray(x,&a)); 2781 PetscFunctionReturn(0); 2782 } 2783 2784 static PetscErrorCode MatMissingDiagonal_SeqDense(Mat A,PetscBool *missing,PetscInt *d) 2785 { 2786 PetscFunctionBegin; 2787 *missing = PETSC_FALSE; 2788 PetscFunctionReturn(0); 2789 } 2790 2791 /* vals is not const */ 2792 static PetscErrorCode MatDenseGetColumn_SeqDense(Mat A,PetscInt col,PetscScalar **vals) 2793 { 2794 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 2795 PetscScalar *v; 2796 2797 PetscFunctionBegin; 2798 PetscCheck(!A->factortype,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 2799 PetscCall(MatDenseGetArray(A,&v)); 2800 *vals = v+col*a->lda; 2801 PetscCall(MatDenseRestoreArray(A,&v)); 2802 PetscFunctionReturn(0); 2803 } 2804 2805 static PetscErrorCode MatDenseRestoreColumn_SeqDense(Mat A,PetscScalar **vals) 2806 { 2807 PetscFunctionBegin; 2808 *vals = NULL; /* user cannot accidentally use the array later */ 2809 PetscFunctionReturn(0); 2810 } 2811 2812 /* -------------------------------------------------------------------*/ 2813 static struct _MatOps MatOps_Values = { MatSetValues_SeqDense, 2814 MatGetRow_SeqDense, 2815 MatRestoreRow_SeqDense, 2816 MatMult_SeqDense, 2817 /* 4*/ MatMultAdd_SeqDense, 2818 MatMultTranspose_SeqDense, 2819 MatMultTransposeAdd_SeqDense, 2820 NULL, 2821 NULL, 2822 NULL, 2823 /* 10*/ NULL, 2824 MatLUFactor_SeqDense, 2825 MatCholeskyFactor_SeqDense, 2826 MatSOR_SeqDense, 2827 MatTranspose_SeqDense, 2828 /* 15*/ MatGetInfo_SeqDense, 2829 MatEqual_SeqDense, 2830 MatGetDiagonal_SeqDense, 2831 MatDiagonalScale_SeqDense, 2832 MatNorm_SeqDense, 2833 /* 20*/ MatAssemblyBegin_SeqDense, 2834 MatAssemblyEnd_SeqDense, 2835 MatSetOption_SeqDense, 2836 MatZeroEntries_SeqDense, 2837 /* 24*/ MatZeroRows_SeqDense, 2838 NULL, 2839 NULL, 2840 NULL, 2841 NULL, 2842 /* 29*/ MatSetUp_SeqDense, 2843 NULL, 2844 NULL, 2845 NULL, 2846 NULL, 2847 /* 34*/ MatDuplicate_SeqDense, 2848 NULL, 2849 NULL, 2850 NULL, 2851 NULL, 2852 /* 39*/ MatAXPY_SeqDense, 2853 MatCreateSubMatrices_SeqDense, 2854 NULL, 2855 MatGetValues_SeqDense, 2856 MatCopy_SeqDense, 2857 /* 44*/ MatGetRowMax_SeqDense, 2858 MatScale_SeqDense, 2859 MatShift_SeqDense, 2860 NULL, 2861 MatZeroRowsColumns_SeqDense, 2862 /* 49*/ MatSetRandom_SeqDense, 2863 NULL, 2864 NULL, 2865 NULL, 2866 NULL, 2867 /* 54*/ NULL, 2868 NULL, 2869 NULL, 2870 NULL, 2871 NULL, 2872 /* 59*/ MatCreateSubMatrix_SeqDense, 2873 MatDestroy_SeqDense, 2874 MatView_SeqDense, 2875 NULL, 2876 NULL, 2877 /* 64*/ NULL, 2878 NULL, 2879 NULL, 2880 NULL, 2881 NULL, 2882 /* 69*/ MatGetRowMaxAbs_SeqDense, 2883 NULL, 2884 NULL, 2885 NULL, 2886 NULL, 2887 /* 74*/ NULL, 2888 NULL, 2889 NULL, 2890 NULL, 2891 NULL, 2892 /* 79*/ NULL, 2893 NULL, 2894 NULL, 2895 NULL, 2896 /* 83*/ MatLoad_SeqDense, 2897 MatIsSymmetric_SeqDense, 2898 MatIsHermitian_SeqDense, 2899 NULL, 2900 NULL, 2901 NULL, 2902 /* 89*/ NULL, 2903 NULL, 2904 MatMatMultNumeric_SeqDense_SeqDense, 2905 NULL, 2906 NULL, 2907 /* 94*/ NULL, 2908 NULL, 2909 NULL, 2910 MatMatTransposeMultNumeric_SeqDense_SeqDense, 2911 NULL, 2912 /* 99*/ MatProductSetFromOptions_SeqDense, 2913 NULL, 2914 NULL, 2915 MatConjugate_SeqDense, 2916 NULL, 2917 /*104*/ NULL, 2918 MatRealPart_SeqDense, 2919 MatImaginaryPart_SeqDense, 2920 NULL, 2921 NULL, 2922 /*109*/ NULL, 2923 NULL, 2924 MatGetRowMin_SeqDense, 2925 MatGetColumnVector_SeqDense, 2926 MatMissingDiagonal_SeqDense, 2927 /*114*/ NULL, 2928 NULL, 2929 NULL, 2930 NULL, 2931 NULL, 2932 /*119*/ NULL, 2933 NULL, 2934 NULL, 2935 NULL, 2936 NULL, 2937 /*124*/ NULL, 2938 MatGetColumnReductions_SeqDense, 2939 NULL, 2940 NULL, 2941 NULL, 2942 /*129*/ NULL, 2943 NULL, 2944 NULL, 2945 MatTransposeMatMultNumeric_SeqDense_SeqDense, 2946 NULL, 2947 /*134*/ NULL, 2948 NULL, 2949 NULL, 2950 NULL, 2951 NULL, 2952 /*139*/ NULL, 2953 NULL, 2954 NULL, 2955 NULL, 2956 NULL, 2957 MatCreateMPIMatConcatenateSeqMat_SeqDense, 2958 /*145*/ NULL, 2959 NULL, 2960 NULL 2961 }; 2962 2963 /*@C 2964 MatCreateSeqDense - Creates a sequential dense matrix that 2965 is stored in column major order (the usual Fortran 77 manner). Many 2966 of the matrix operations use the BLAS and LAPACK routines. 2967 2968 Collective 2969 2970 Input Parameters: 2971 + comm - MPI communicator, set to PETSC_COMM_SELF 2972 . m - number of rows 2973 . n - number of columns 2974 - data - optional location of matrix data in column major order. Set data=NULL for PETSc 2975 to control all matrix memory allocation. 2976 2977 Output Parameter: 2978 . A - the matrix 2979 2980 Notes: 2981 The data input variable is intended primarily for Fortran programmers 2982 who wish to allocate their own matrix memory space. Most users should 2983 set data=NULL. 2984 2985 Level: intermediate 2986 2987 .seealso: `MatCreate()`, `MatCreateDense()`, `MatSetValues()` 2988 @*/ 2989 PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 2990 { 2991 PetscFunctionBegin; 2992 PetscCall(MatCreate(comm,A)); 2993 PetscCall(MatSetSizes(*A,m,n,m,n)); 2994 PetscCall(MatSetType(*A,MATSEQDENSE)); 2995 PetscCall(MatSeqDenseSetPreallocation(*A,data)); 2996 PetscFunctionReturn(0); 2997 } 2998 2999 /*@C 3000 MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 3001 3002 Collective 3003 3004 Input Parameters: 3005 + B - the matrix 3006 - data - the array (or NULL) 3007 3008 Notes: 3009 The data input variable is intended primarily for Fortran programmers 3010 who wish to allocate their own matrix memory space. Most users should 3011 need not call this routine. 3012 3013 Level: intermediate 3014 3015 .seealso: `MatCreate()`, `MatCreateDense()`, `MatSetValues()`, `MatDenseSetLDA()` 3016 3017 @*/ 3018 PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 3019 { 3020 PetscFunctionBegin; 3021 PetscValidHeaderSpecific(B,MAT_CLASSID,1); 3022 PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data)); 3023 PetscFunctionReturn(0); 3024 } 3025 3026 PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 3027 { 3028 Mat_SeqDense *b = (Mat_SeqDense*)B->data; 3029 3030 PetscFunctionBegin; 3031 PetscCheck(!b->matinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first"); 3032 B->preallocated = PETSC_TRUE; 3033 3034 PetscCall(PetscLayoutSetUp(B->rmap)); 3035 PetscCall(PetscLayoutSetUp(B->cmap)); 3036 3037 if (b->lda <= 0) b->lda = B->rmap->n; 3038 3039 if (!data) { /* petsc-allocated storage */ 3040 if (!b->user_alloc) PetscCall(PetscFree(b->v)); 3041 PetscCall(PetscCalloc1((size_t)b->lda*B->cmap->n,&b->v)); 3042 PetscCall(PetscLogObjectMemory((PetscObject)B,b->lda*B->cmap->n*sizeof(PetscScalar))); 3043 3044 b->user_alloc = PETSC_FALSE; 3045 } else { /* user-allocated storage */ 3046 if (!b->user_alloc) PetscCall(PetscFree(b->v)); 3047 b->v = data; 3048 b->user_alloc = PETSC_TRUE; 3049 } 3050 B->assembled = PETSC_TRUE; 3051 PetscFunctionReturn(0); 3052 } 3053 3054 #if defined(PETSC_HAVE_ELEMENTAL) 3055 PETSC_INTERN PetscErrorCode MatConvert_SeqDense_Elemental(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 3056 { 3057 Mat mat_elemental; 3058 const PetscScalar *array; 3059 PetscScalar *v_colwise; 3060 PetscInt M=A->rmap->N,N=A->cmap->N,i,j,k,*rows,*cols; 3061 3062 PetscFunctionBegin; 3063 PetscCall(PetscMalloc3(M*N,&v_colwise,M,&rows,N,&cols)); 3064 PetscCall(MatDenseGetArrayRead(A,&array)); 3065 /* convert column-wise array into row-wise v_colwise, see MatSetValues_Elemental() */ 3066 k = 0; 3067 for (j=0; j<N; j++) { 3068 cols[j] = j; 3069 for (i=0; i<M; i++) { 3070 v_colwise[j*M+i] = array[k++]; 3071 } 3072 } 3073 for (i=0; i<M; i++) { 3074 rows[i] = i; 3075 } 3076 PetscCall(MatDenseRestoreArrayRead(A,&array)); 3077 3078 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &mat_elemental)); 3079 PetscCall(MatSetSizes(mat_elemental,PETSC_DECIDE,PETSC_DECIDE,M,N)); 3080 PetscCall(MatSetType(mat_elemental,MATELEMENTAL)); 3081 PetscCall(MatSetUp(mat_elemental)); 3082 3083 /* PETSc-Elemental interaface uses axpy for setting off-processor entries, only ADD_VALUES is allowed */ 3084 PetscCall(MatSetValues(mat_elemental,M,rows,N,cols,v_colwise,ADD_VALUES)); 3085 PetscCall(MatAssemblyBegin(mat_elemental, MAT_FINAL_ASSEMBLY)); 3086 PetscCall(MatAssemblyEnd(mat_elemental, MAT_FINAL_ASSEMBLY)); 3087 PetscCall(PetscFree3(v_colwise,rows,cols)); 3088 3089 if (reuse == MAT_INPLACE_MATRIX) { 3090 PetscCall(MatHeaderReplace(A,&mat_elemental)); 3091 } else { 3092 *newmat = mat_elemental; 3093 } 3094 PetscFunctionReturn(0); 3095 } 3096 #endif 3097 3098 PetscErrorCode MatDenseSetLDA_SeqDense(Mat B,PetscInt lda) 3099 { 3100 Mat_SeqDense *b = (Mat_SeqDense*)B->data; 3101 PetscBool data; 3102 3103 PetscFunctionBegin; 3104 data = (PetscBool)((B->rmap->n > 0 && B->cmap->n > 0) ? (b->v ? PETSC_TRUE : PETSC_FALSE) : PETSC_FALSE); 3105 PetscCheck(b->user_alloc || !data || b->lda == lda,PETSC_COMM_SELF,PETSC_ERR_ORDER,"LDA cannot be changed after allocation of internal storage"); 3106 PetscCheck(lda >= B->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"LDA %" PetscInt_FMT " must be at least matrix dimension %" PetscInt_FMT,lda,B->rmap->n); 3107 b->lda = lda; 3108 PetscFunctionReturn(0); 3109 } 3110 3111 PetscErrorCode MatCreateMPIMatConcatenateSeqMat_SeqDense(MPI_Comm comm,Mat inmat,PetscInt n,MatReuse scall,Mat *outmat) 3112 { 3113 PetscMPIInt size; 3114 3115 PetscFunctionBegin; 3116 PetscCallMPI(MPI_Comm_size(comm,&size)); 3117 if (size == 1) { 3118 if (scall == MAT_INITIAL_MATRIX) { 3119 PetscCall(MatDuplicate(inmat,MAT_COPY_VALUES,outmat)); 3120 } else { 3121 PetscCall(MatCopy(inmat,*outmat,SAME_NONZERO_PATTERN)); 3122 } 3123 } else { 3124 PetscCall(MatCreateMPIMatConcatenateSeqMat_MPIDense(comm,inmat,n,scall,outmat)); 3125 } 3126 PetscFunctionReturn(0); 3127 } 3128 3129 PetscErrorCode MatDenseGetColumnVec_SeqDense(Mat A,PetscInt col,Vec *v) 3130 { 3131 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 3132 3133 PetscFunctionBegin; 3134 PetscCheck(!a->vecinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first"); 3135 PetscCheck(!a->matinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first"); 3136 if (!a->cvec) { 3137 PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)A),A->rmap->bs,A->rmap->n,NULL,&a->cvec)); 3138 PetscCall(PetscLogObjectParent((PetscObject)A,(PetscObject)a->cvec)); 3139 } 3140 a->vecinuse = col + 1; 3141 PetscCall(MatDenseGetArray(A,(PetscScalar**)&a->ptrinuse)); 3142 PetscCall(VecPlaceArray(a->cvec,a->ptrinuse + (size_t)col * (size_t)a->lda)); 3143 *v = a->cvec; 3144 PetscFunctionReturn(0); 3145 } 3146 3147 PetscErrorCode MatDenseRestoreColumnVec_SeqDense(Mat A,PetscInt col,Vec *v) 3148 { 3149 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 3150 3151 PetscFunctionBegin; 3152 PetscCheck(a->vecinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseGetColumnVec() first"); 3153 PetscCheck(a->cvec,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing internal column vector"); 3154 a->vecinuse = 0; 3155 PetscCall(MatDenseRestoreArray(A,(PetscScalar**)&a->ptrinuse)); 3156 PetscCall(VecResetArray(a->cvec)); 3157 if (v) *v = NULL; 3158 PetscFunctionReturn(0); 3159 } 3160 3161 PetscErrorCode MatDenseGetColumnVecRead_SeqDense(Mat A,PetscInt col,Vec *v) 3162 { 3163 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 3164 3165 PetscFunctionBegin; 3166 PetscCheck(!a->vecinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first"); 3167 PetscCheck(!a->matinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first"); 3168 if (!a->cvec) { 3169 PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)A),A->rmap->bs,A->rmap->n,NULL,&a->cvec)); 3170 PetscCall(PetscLogObjectParent((PetscObject)A,(PetscObject)a->cvec)); 3171 } 3172 a->vecinuse = col + 1; 3173 PetscCall(MatDenseGetArrayRead(A,&a->ptrinuse)); 3174 PetscCall(VecPlaceArray(a->cvec,a->ptrinuse + (size_t)col * (size_t)a->lda)); 3175 PetscCall(VecLockReadPush(a->cvec)); 3176 *v = a->cvec; 3177 PetscFunctionReturn(0); 3178 } 3179 3180 PetscErrorCode MatDenseRestoreColumnVecRead_SeqDense(Mat A,PetscInt col,Vec *v) 3181 { 3182 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 3183 3184 PetscFunctionBegin; 3185 PetscCheck(a->vecinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseGetColumnVec() first"); 3186 PetscCheck(a->cvec,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing internal column vector"); 3187 a->vecinuse = 0; 3188 PetscCall(MatDenseRestoreArrayRead(A,&a->ptrinuse)); 3189 PetscCall(VecLockReadPop(a->cvec)); 3190 PetscCall(VecResetArray(a->cvec)); 3191 if (v) *v = NULL; 3192 PetscFunctionReturn(0); 3193 } 3194 3195 PetscErrorCode MatDenseGetColumnVecWrite_SeqDense(Mat A,PetscInt col,Vec *v) 3196 { 3197 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 3198 3199 PetscFunctionBegin; 3200 PetscCheck(!a->vecinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first"); 3201 PetscCheck(!a->matinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first"); 3202 if (!a->cvec) { 3203 PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)A),A->rmap->bs,A->rmap->n,NULL,&a->cvec)); 3204 PetscCall(PetscLogObjectParent((PetscObject)A,(PetscObject)a->cvec)); 3205 } 3206 a->vecinuse = col + 1; 3207 PetscCall(MatDenseGetArrayWrite(A,(PetscScalar**)&a->ptrinuse)); 3208 PetscCall(VecPlaceArray(a->cvec,a->ptrinuse + (size_t)col * (size_t)a->lda)); 3209 *v = a->cvec; 3210 PetscFunctionReturn(0); 3211 } 3212 3213 PetscErrorCode MatDenseRestoreColumnVecWrite_SeqDense(Mat A,PetscInt col,Vec *v) 3214 { 3215 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 3216 3217 PetscFunctionBegin; 3218 PetscCheck(a->vecinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseGetColumnVec() first"); 3219 PetscCheck(a->cvec,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing internal column vector"); 3220 a->vecinuse = 0; 3221 PetscCall(MatDenseRestoreArrayWrite(A,(PetscScalar**)&a->ptrinuse)); 3222 PetscCall(VecResetArray(a->cvec)); 3223 if (v) *v = NULL; 3224 PetscFunctionReturn(0); 3225 } 3226 3227 PetscErrorCode MatDenseGetSubMatrix_SeqDense(Mat A,PetscInt cbegin,PetscInt cend,Mat *v) 3228 { 3229 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 3230 3231 PetscFunctionBegin; 3232 PetscCheck(!a->vecinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreColumnVec() first"); 3233 PetscCheck(!a->matinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseRestoreSubMatrix() first"); 3234 if (a->cmat && cend-cbegin != a->cmat->cmap->N) { 3235 PetscCall(MatDestroy(&a->cmat)); 3236 } 3237 if (!a->cmat) { 3238 PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A),A->rmap->n,PETSC_DECIDE,A->rmap->N,cend-cbegin,a->v+(size_t)cbegin*a->lda,&a->cmat)); 3239 PetscCall(PetscLogObjectParent((PetscObject)A,(PetscObject)a->cmat)); 3240 } else { 3241 PetscCall(MatDensePlaceArray(a->cmat,a->v+(size_t)cbegin*a->lda)); 3242 } 3243 PetscCall(MatDenseSetLDA(a->cmat,a->lda)); 3244 a->matinuse = cbegin + 1; 3245 *v = a->cmat; 3246 #if defined(PETSC_HAVE_CUDA) 3247 A->offloadmask = PETSC_OFFLOAD_CPU; 3248 #endif 3249 PetscFunctionReturn(0); 3250 } 3251 3252 PetscErrorCode MatDenseRestoreSubMatrix_SeqDense(Mat A,Mat *v) 3253 { 3254 Mat_SeqDense *a = (Mat_SeqDense*)A->data; 3255 3256 PetscFunctionBegin; 3257 PetscCheck(a->matinuse,PETSC_COMM_SELF,PETSC_ERR_ORDER,"Need to call MatDenseGetSubMatrix() first"); 3258 PetscCheck(a->cmat,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing internal column matrix"); 3259 PetscCheck(*v == a->cmat,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not the matrix obtained from MatDenseGetSubMatrix()"); 3260 a->matinuse = 0; 3261 PetscCall(MatDenseResetArray(a->cmat)); 3262 *v = NULL; 3263 PetscFunctionReturn(0); 3264 } 3265 3266 /*MC 3267 MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 3268 3269 Options Database Keys: 3270 . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 3271 3272 Level: beginner 3273 3274 .seealso: `MatCreateSeqDense()` 3275 3276 M*/ 3277 PetscErrorCode MatCreate_SeqDense(Mat B) 3278 { 3279 Mat_SeqDense *b; 3280 PetscMPIInt size; 3281 3282 PetscFunctionBegin; 3283 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)B),&size)); 3284 PetscCheck(size <= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 3285 3286 PetscCall(PetscNewLog(B,&b)); 3287 PetscCall(PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps))); 3288 B->data = (void*)b; 3289 3290 b->roworiented = PETSC_TRUE; 3291 3292 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatQRFactor_C",MatQRFactor_SeqDense)); 3293 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseGetLDA_C",MatDenseGetLDA_SeqDense)); 3294 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseSetLDA_C",MatDenseSetLDA_SeqDense)); 3295 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArray_C",MatDenseGetArray_SeqDense)); 3296 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArray_C",MatDenseRestoreArray_SeqDense)); 3297 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDensePlaceArray_C",MatDensePlaceArray_SeqDense)); 3298 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseResetArray_C",MatDenseResetArray_SeqDense)); 3299 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseReplaceArray_C",MatDenseReplaceArray_SeqDense)); 3300 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArrayRead_C",MatDenseGetArray_SeqDense)); 3301 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArrayRead_C",MatDenseRestoreArray_SeqDense)); 3302 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArrayWrite_C",MatDenseGetArray_SeqDense)); 3303 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArrayWrite_C",MatDenseRestoreArray_SeqDense)); 3304 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqaij_C",MatConvert_SeqDense_SeqAIJ)); 3305 #if defined(PETSC_HAVE_ELEMENTAL) 3306 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_elemental_C",MatConvert_SeqDense_Elemental)); 3307 #endif 3308 #if defined(PETSC_HAVE_SCALAPACK) 3309 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_scalapack_C",MatConvert_Dense_ScaLAPACK)); 3310 #endif 3311 #if defined(PETSC_HAVE_CUDA) 3312 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqdensecuda_C",MatConvert_SeqDense_SeqDenseCUDA)); 3313 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqdensecuda_seqdensecuda_C",MatProductSetFromOptions_SeqDense)); 3314 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqdensecuda_seqdense_C",MatProductSetFromOptions_SeqDense)); 3315 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqdense_seqdensecuda_C",MatProductSetFromOptions_SeqDense)); 3316 #endif 3317 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",MatSeqDenseSetPreallocation_SeqDense)); 3318 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqaij_seqdense_C",MatProductSetFromOptions_SeqAIJ_SeqDense)); 3319 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqdense_seqdense_C",MatProductSetFromOptions_SeqDense)); 3320 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqbaij_seqdense_C",MatProductSetFromOptions_SeqXBAIJ_SeqDense)); 3321 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatProductSetFromOptions_seqsbaij_seqdense_C",MatProductSetFromOptions_SeqXBAIJ_SeqDense)); 3322 3323 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseGetColumn_C",MatDenseGetColumn_SeqDense)); 3324 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreColumn_C",MatDenseRestoreColumn_SeqDense)); 3325 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseGetColumnVec_C",MatDenseGetColumnVec_SeqDense)); 3326 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreColumnVec_C",MatDenseRestoreColumnVec_SeqDense)); 3327 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseGetColumnVecRead_C",MatDenseGetColumnVecRead_SeqDense)); 3328 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreColumnVecRead_C",MatDenseRestoreColumnVecRead_SeqDense)); 3329 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseGetColumnVecWrite_C",MatDenseGetColumnVecWrite_SeqDense)); 3330 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreColumnVecWrite_C",MatDenseRestoreColumnVecWrite_SeqDense)); 3331 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseGetSubMatrix_C",MatDenseGetSubMatrix_SeqDense)); 3332 PetscCall(PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreSubMatrix_C",MatDenseRestoreSubMatrix_SeqDense)); 3333 PetscCall(PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE)); 3334 PetscFunctionReturn(0); 3335 } 3336 3337 /*@C 3338 MatDenseGetColumn - gives access to a column of a dense matrix. This is only the local part of the column. You MUST call MatDenseRestoreColumn() to avoid memory bleeding. 3339 3340 Not Collective 3341 3342 Input Parameters: 3343 + mat - a MATSEQDENSE or MATMPIDENSE matrix 3344 - col - column index 3345 3346 Output Parameter: 3347 . vals - pointer to the data 3348 3349 Level: intermediate 3350 3351 .seealso: `MatDenseRestoreColumn()` 3352 @*/ 3353 PetscErrorCode MatDenseGetColumn(Mat A,PetscInt col,PetscScalar **vals) 3354 { 3355 PetscFunctionBegin; 3356 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3357 PetscValidLogicalCollectiveInt(A,col,2); 3358 PetscValidPointer(vals,3); 3359 PetscUseMethod(A,"MatDenseGetColumn_C",(Mat,PetscInt,PetscScalar**),(A,col,vals)); 3360 PetscFunctionReturn(0); 3361 } 3362 3363 /*@C 3364 MatDenseRestoreColumn - returns access to a column of a dense matrix which is returned by MatDenseGetColumn(). 3365 3366 Not Collective 3367 3368 Input Parameter: 3369 . mat - a MATSEQDENSE or MATMPIDENSE matrix 3370 3371 Output Parameter: 3372 . vals - pointer to the data 3373 3374 Level: intermediate 3375 3376 .seealso: `MatDenseGetColumn()` 3377 @*/ 3378 PetscErrorCode MatDenseRestoreColumn(Mat A,PetscScalar **vals) 3379 { 3380 PetscFunctionBegin; 3381 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3382 PetscValidPointer(vals,2); 3383 PetscUseMethod(A,"MatDenseRestoreColumn_C",(Mat,PetscScalar**),(A,vals)); 3384 PetscFunctionReturn(0); 3385 } 3386 3387 /*@ 3388 MatDenseGetColumnVec - Gives read-write access to a column of a dense matrix, represented as a Vec. 3389 3390 Collective 3391 3392 Input Parameters: 3393 + mat - the Mat object 3394 - col - the column index 3395 3396 Output Parameter: 3397 . v - the vector 3398 3399 Notes: 3400 The vector is owned by PETSc. Users need to call MatDenseRestoreColumnVec() when the vector is no longer needed. 3401 Use MatDenseGetColumnVecRead() to obtain read-only access or MatDenseGetColumnVecWrite() for write-only access. 3402 3403 Level: intermediate 3404 3405 .seealso: `MATDENSE`, `MATDENSECUDA`, `MatDenseGetColumnVecRead()`, `MatDenseGetColumnVecWrite()`, `MatDenseRestoreColumnVec()`, `MatDenseRestoreColumnVecRead()`, `MatDenseRestoreColumnVecWrite()` 3406 @*/ 3407 PetscErrorCode MatDenseGetColumnVec(Mat A,PetscInt col,Vec *v) 3408 { 3409 PetscFunctionBegin; 3410 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3411 PetscValidType(A,1); 3412 PetscValidLogicalCollectiveInt(A,col,2); 3413 PetscValidPointer(v,3); 3414 PetscCheck(A->preallocated,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated"); 3415 PetscCheck(col >= 0 && col < A->cmap->N,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid col %" PetscInt_FMT ", should be in [0,%" PetscInt_FMT ")",col,A->cmap->N); 3416 PetscUseMethod(A,"MatDenseGetColumnVec_C",(Mat,PetscInt,Vec*),(A,col,v)); 3417 PetscFunctionReturn(0); 3418 } 3419 3420 /*@ 3421 MatDenseRestoreColumnVec - Returns access to a column of a dense matrix obtained from MatDenseGetColumnVec(). 3422 3423 Collective 3424 3425 Input Parameters: 3426 + mat - the Mat object 3427 . col - the column index 3428 - v - the Vec object 3429 3430 Level: intermediate 3431 3432 .seealso: `MATDENSE`, `MATDENSECUDA`, `MatDenseGetColumnVec()`, `MatDenseGetColumnVecRead()`, `MatDenseGetColumnVecWrite()`, `MatDenseRestoreColumnVecRead()`, `MatDenseRestoreColumnVecWrite()` 3433 @*/ 3434 PetscErrorCode MatDenseRestoreColumnVec(Mat A,PetscInt col,Vec *v) 3435 { 3436 PetscFunctionBegin; 3437 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3438 PetscValidType(A,1); 3439 PetscValidLogicalCollectiveInt(A,col,2); 3440 PetscCheck(A->preallocated,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated"); 3441 PetscCheck(col >= 0 && col < A->cmap->N,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid col %" PetscInt_FMT ", should be in [0,%" PetscInt_FMT ")",col,A->cmap->N); 3442 PetscUseMethod(A,"MatDenseRestoreColumnVec_C",(Mat,PetscInt,Vec*),(A,col,v)); 3443 PetscFunctionReturn(0); 3444 } 3445 3446 /*@ 3447 MatDenseGetColumnVecRead - Gives read-only access to a column of a dense matrix, represented as a Vec. 3448 3449 Collective 3450 3451 Input Parameters: 3452 + mat - the Mat object 3453 - col - the column index 3454 3455 Output Parameter: 3456 . v - the vector 3457 3458 Notes: 3459 The vector is owned by PETSc and users cannot modify it. 3460 Users need to call MatDenseRestoreColumnVecRead() when the vector is no longer needed. 3461 Use MatDenseGetColumnVec() to obtain read-write access or MatDenseGetColumnVecWrite() for write-only access. 3462 3463 Level: intermediate 3464 3465 .seealso: `MATDENSE`, `MATDENSECUDA`, `MatDenseGetColumnVec()`, `MatDenseGetColumnVecWrite()`, `MatDenseRestoreColumnVec()`, `MatDenseRestoreColumnVecRead()`, `MatDenseRestoreColumnVecWrite()` 3466 @*/ 3467 PetscErrorCode MatDenseGetColumnVecRead(Mat A,PetscInt col,Vec *v) 3468 { 3469 PetscFunctionBegin; 3470 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3471 PetscValidType(A,1); 3472 PetscValidLogicalCollectiveInt(A,col,2); 3473 PetscValidPointer(v,3); 3474 PetscCheck(A->preallocated,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated"); 3475 PetscCheck(col >= 0 && col < A->cmap->N,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid col %" PetscInt_FMT ", should be in [0,%" PetscInt_FMT ")",col,A->cmap->N); 3476 PetscUseMethod(A,"MatDenseGetColumnVecRead_C",(Mat,PetscInt,Vec*),(A,col,v)); 3477 PetscFunctionReturn(0); 3478 } 3479 3480 /*@ 3481 MatDenseRestoreColumnVecRead - Returns access to a column of a dense matrix obtained from MatDenseGetColumnVecRead(). 3482 3483 Collective 3484 3485 Input Parameters: 3486 + mat - the Mat object 3487 . col - the column index 3488 - v - the Vec object 3489 3490 Level: intermediate 3491 3492 .seealso: `MATDENSE`, `MATDENSECUDA`, `MatDenseGetColumnVec()`, `MatDenseGetColumnVecRead()`, `MatDenseGetColumnVecWrite()`, `MatDenseRestoreColumnVec()`, `MatDenseRestoreColumnVecWrite()` 3493 @*/ 3494 PetscErrorCode MatDenseRestoreColumnVecRead(Mat A,PetscInt col,Vec *v) 3495 { 3496 PetscFunctionBegin; 3497 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3498 PetscValidType(A,1); 3499 PetscValidLogicalCollectiveInt(A,col,2); 3500 PetscCheck(A->preallocated,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated"); 3501 PetscCheck(col >= 0 && col < A->cmap->N,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid col %" PetscInt_FMT ", should be in [0,%" PetscInt_FMT ")",col,A->cmap->N); 3502 PetscUseMethod(A,"MatDenseRestoreColumnVecRead_C",(Mat,PetscInt,Vec*),(A,col,v)); 3503 PetscFunctionReturn(0); 3504 } 3505 3506 /*@ 3507 MatDenseGetColumnVecWrite - Gives write-only access to a column of a dense matrix, represented as a Vec. 3508 3509 Collective 3510 3511 Input Parameters: 3512 + mat - the Mat object 3513 - col - the column index 3514 3515 Output Parameter: 3516 . v - the vector 3517 3518 Notes: 3519 The vector is owned by PETSc. Users need to call MatDenseRestoreColumnVecWrite() when the vector is no longer needed. 3520 Use MatDenseGetColumnVec() to obtain read-write access or MatDenseGetColumnVecRead() for read-only access. 3521 3522 Level: intermediate 3523 3524 .seealso: `MATDENSE`, `MATDENSECUDA`, `MatDenseGetColumnVec()`, `MatDenseGetColumnVecRead()`, `MatDenseRestoreColumnVec()`, `MatDenseRestoreColumnVecRead()`, `MatDenseRestoreColumnVecWrite()` 3525 @*/ 3526 PetscErrorCode MatDenseGetColumnVecWrite(Mat A,PetscInt col,Vec *v) 3527 { 3528 PetscFunctionBegin; 3529 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3530 PetscValidType(A,1); 3531 PetscValidLogicalCollectiveInt(A,col,2); 3532 PetscValidPointer(v,3); 3533 PetscCheck(A->preallocated,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated"); 3534 PetscCheck(col >= 0 && col < A->cmap->N,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid col %" PetscInt_FMT ", should be in [0,%" PetscInt_FMT ")",col,A->cmap->N); 3535 PetscUseMethod(A,"MatDenseGetColumnVecWrite_C",(Mat,PetscInt,Vec*),(A,col,v)); 3536 PetscFunctionReturn(0); 3537 } 3538 3539 /*@ 3540 MatDenseRestoreColumnVecWrite - Returns access to a column of a dense matrix obtained from MatDenseGetColumnVecWrite(). 3541 3542 Collective 3543 3544 Input Parameters: 3545 + mat - the Mat object 3546 . col - the column index 3547 - v - the Vec object 3548 3549 Level: intermediate 3550 3551 .seealso: `MATDENSE`, `MATDENSECUDA`, `MatDenseGetColumnVec()`, `MatDenseGetColumnVecRead()`, `MatDenseGetColumnVecWrite()`, `MatDenseRestoreColumnVec()`, `MatDenseRestoreColumnVecRead()` 3552 @*/ 3553 PetscErrorCode MatDenseRestoreColumnVecWrite(Mat A,PetscInt col,Vec *v) 3554 { 3555 PetscFunctionBegin; 3556 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3557 PetscValidType(A,1); 3558 PetscValidLogicalCollectiveInt(A,col,2); 3559 PetscCheck(A->preallocated,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated"); 3560 PetscCheck(col >= 0 && col < A->cmap->N,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid col %" PetscInt_FMT ", should be in [0,%" PetscInt_FMT ")",col,A->cmap->N); 3561 PetscUseMethod(A,"MatDenseRestoreColumnVecWrite_C",(Mat,PetscInt,Vec*),(A,col,v)); 3562 PetscFunctionReturn(0); 3563 } 3564 3565 /*@ 3566 MatDenseGetSubMatrix - Gives access to a block of columns of a dense matrix, represented as a Mat. 3567 3568 Collective 3569 3570 Input Parameters: 3571 + mat - the Mat object 3572 . cbegin - the first index in the block 3573 - cend - the index past the last one in the block 3574 3575 Output Parameter: 3576 . v - the matrix 3577 3578 Notes: 3579 The matrix is owned by PETSc. Users need to call MatDenseRestoreSubMatrix() when the matrix is no longer needed. 3580 3581 Level: intermediate 3582 3583 .seealso: `MATDENSE`, `MATDENSECUDA`, `MatDenseGetColumnVec()`, `MatDenseRestoreColumnVec()`, `MatDenseRestoreSubMatrix()` 3584 @*/ 3585 PetscErrorCode MatDenseGetSubMatrix(Mat A,PetscInt cbegin,PetscInt cend,Mat *v) 3586 { 3587 PetscFunctionBegin; 3588 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3589 PetscValidType(A,1); 3590 PetscValidLogicalCollectiveInt(A,cbegin,2); 3591 PetscValidLogicalCollectiveInt(A,cend,3); 3592 PetscValidPointer(v,4); 3593 PetscCheck(A->preallocated,PetscObjectComm((PetscObject)A),PETSC_ERR_ORDER,"Matrix not preallocated"); 3594 PetscCheck(cbegin >= 0 && cbegin < A->cmap->N,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid cbegin %" PetscInt_FMT ", should be in [0,%" PetscInt_FMT ")",cbegin,A->cmap->N); 3595 PetscCheck(cend > cbegin && cend <= A->cmap->N,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Invalid cend %" PetscInt_FMT ", should be in (%" PetscInt_FMT ",%" PetscInt_FMT "]",cend,cbegin,A->cmap->N); 3596 PetscUseMethod(A,"MatDenseGetSubMatrix_C",(Mat,PetscInt,PetscInt,Mat*),(A,cbegin,cend,v)); 3597 PetscFunctionReturn(0); 3598 } 3599 3600 /*@ 3601 MatDenseRestoreSubMatrix - Returns access to a block of columns of a dense matrix obtained from MatDenseGetSubMatrix(). 3602 3603 Collective 3604 3605 Input Parameters: 3606 + mat - the Mat object 3607 - v - the Mat object 3608 3609 Level: intermediate 3610 3611 .seealso: `MATDENSE`, `MATDENSECUDA`, `MatDenseGetColumnVec()`, `MatDenseRestoreColumnVec()`, `MatDenseGetSubMatrix()` 3612 @*/ 3613 PetscErrorCode MatDenseRestoreSubMatrix(Mat A,Mat *v) 3614 { 3615 PetscFunctionBegin; 3616 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 3617 PetscValidType(A,1); 3618 PetscValidPointer(v,2); 3619 PetscUseMethod(A,"MatDenseRestoreSubMatrix_C",(Mat,Mat*),(A,v)); 3620 PetscFunctionReturn(0); 3621 } 3622 3623 #include <petscblaslapack.h> 3624 #include <petsc/private/kernels/blockinvert.h> 3625 3626 PetscErrorCode MatSeqDenseInvert(Mat A) 3627 { 3628 Mat_SeqDense *a = (Mat_SeqDense*) A->data; 3629 PetscInt bs = A->rmap->n; 3630 MatScalar *values = a->v; 3631 const PetscReal shift = 0.0; 3632 PetscBool allowzeropivot = PetscNot(A->erroriffailure),zeropivotdetected=PETSC_FALSE; 3633 3634 PetscFunctionBegin; 3635 /* factor and invert each block */ 3636 switch (bs) { 3637 case 1: 3638 values[0] = (PetscScalar)1.0 / (values[0] + shift); 3639 break; 3640 case 2: 3641 PetscCall(PetscKernel_A_gets_inverse_A_2(values,shift,allowzeropivot,&zeropivotdetected)); 3642 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3643 break; 3644 case 3: 3645 PetscCall(PetscKernel_A_gets_inverse_A_3(values,shift,allowzeropivot,&zeropivotdetected)); 3646 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3647 break; 3648 case 4: 3649 PetscCall(PetscKernel_A_gets_inverse_A_4(values,shift,allowzeropivot,&zeropivotdetected)); 3650 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3651 break; 3652 case 5: 3653 { 3654 PetscScalar work[25]; 3655 PetscInt ipvt[5]; 3656 3657 PetscCall(PetscKernel_A_gets_inverse_A_5(values,ipvt,work,shift,allowzeropivot,&zeropivotdetected)); 3658 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3659 } 3660 break; 3661 case 6: 3662 PetscCall(PetscKernel_A_gets_inverse_A_6(values,shift,allowzeropivot,&zeropivotdetected)); 3663 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3664 break; 3665 case 7: 3666 PetscCall(PetscKernel_A_gets_inverse_A_7(values,shift,allowzeropivot,&zeropivotdetected)); 3667 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3668 break; 3669 default: 3670 { 3671 PetscInt *v_pivots,*IJ,j; 3672 PetscScalar *v_work; 3673 3674 PetscCall(PetscMalloc3(bs,&v_work,bs,&v_pivots,bs,&IJ)); 3675 for (j=0; j<bs; j++) { 3676 IJ[j] = j; 3677 } 3678 PetscCall(PetscKernel_A_gets_inverse_A(bs,values,v_pivots,v_work,allowzeropivot,&zeropivotdetected)); 3679 if (zeropivotdetected) A->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3680 PetscCall(PetscFree3(v_work,v_pivots,IJ)); 3681 } 3682 } 3683 PetscFunctionReturn(0); 3684 } 3685