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