Lines Matching +full:- +full:b

5   MatMAIJGetAIJ - Get the `MATAIJ` matrix describing the blockwise action of the `MATMAIJ` matrix
10 . A - the `MATMAIJ` matrix
13 . B - the `MATAIJ` matrix
22 PetscErrorCode MatMAIJGetAIJ(Mat A, Mat *B) in MatMAIJGetAIJ() argument
30 Mat_MPIMAIJ *b = (Mat_MPIMAIJ *)A->data; in MatMAIJGetAIJ() local
32 *B = b->A; in MatMAIJGetAIJ()
34 Mat_SeqMAIJ *b = (Mat_SeqMAIJ *)A->data; in MatMAIJGetAIJ() local
36 *B = b->AIJ; in MatMAIJGetAIJ()
38 *B = A; in MatMAIJGetAIJ()
44 …MatMAIJRedimension - Get a new `MATMAIJ` matrix with the same action, but for a different block si…
49 + A - the `MATMAIJ` matrix
50 - dof - the block size for the new matrix
53 . B - the new `MATMAIJ` matrix
59 PetscErrorCode MatMAIJRedimension(Mat A, PetscInt dof, Mat *B) in MatMAIJRedimension() argument
66 PetscCall(MatCreateMAIJ(Aij, dof, B)); in MatMAIJRedimension()
72 Mat_SeqMAIJ *b = (Mat_SeqMAIJ *)A->data; in MatDestroy_SeqMAIJ() local
75 PetscCall(MatDestroy(&b->AIJ)); in MatDestroy_SeqMAIJ()
76 PetscCall(PetscFree(A->data)); in MatDestroy_SeqMAIJ()
91 Mat B; in MatView_SeqMAIJ() local
94 PetscCall(MatConvert(A, MATSEQAIJ, MAT_INITIAL_MATRIX, &B)); in MatView_SeqMAIJ()
95 PetscCall(MatView(B, viewer)); in MatView_SeqMAIJ()
96 PetscCall(MatDestroy(&B)); in MatView_SeqMAIJ()
102 Mat B; in MatView_MPIMAIJ() local
105 PetscCall(MatConvert(A, MATMPIAIJ, MAT_INITIAL_MATRIX, &B)); in MatView_MPIMAIJ()
106 PetscCall(MatView(B, viewer)); in MatView_MPIMAIJ()
107 PetscCall(MatDestroy(&B)); in MatView_MPIMAIJ()
113 Mat_MPIMAIJ *b = (Mat_MPIMAIJ *)A->data; in MatDestroy_MPIMAIJ() local
116 PetscCall(MatDestroy(&b->AIJ)); in MatDestroy_MPIMAIJ()
117 PetscCall(MatDestroy(&b->OAIJ)); in MatDestroy_MPIMAIJ()
118 PetscCall(MatDestroy(&b->A)); in MatDestroy_MPIMAIJ()
119 PetscCall(VecScatterDestroy(&b->ctx)); in MatDestroy_MPIMAIJ()
120 PetscCall(VecDestroy(&b->w)); in MatDestroy_MPIMAIJ()
121 PetscCall(PetscFree(A->data)); in MatDestroy_MPIMAIJ()
130 …MATMAIJ - MATMAIJ = "maij" - A matrix type to be used for restriction and interpolation operations…
149 Mat_MPIMAIJ *b; in MatCreate_MAIJ() local
153 PetscCall(PetscNew(&b)); in MatCreate_MAIJ()
154 A->data = (void *)b; in MatCreate_MAIJ()
156 PetscCall(PetscMemzero(A->ops, sizeof(struct _MatOps))); in MatCreate_MAIJ()
158 A->ops->setup = MatSetUp_MAIJ; in MatCreate_MAIJ()
160 b->AIJ = NULL; in MatCreate_MAIJ()
161 b->dof = 0; in MatCreate_MAIJ()
162 b->OAIJ = NULL; in MatCreate_MAIJ()
163 b->ctx = NULL; in MatCreate_MAIJ()
164 b->w = NULL; in MatCreate_MAIJ()
171 A->preallocated = PETSC_TRUE; in MatCreate_MAIJ()
172 A->assembled = PETSC_TRUE; in MatCreate_MAIJ()
197 const Mat_SeqMAIJ *b = (Mat_SeqMAIJ *)A->data; in MatMult_MatMultAdd_SeqMAIJ_Template() local
198 const Mat baij = b->AIJ; in MatMult_MatMultAdd_SeqMAIJ_Template()
199 const Mat_SeqAIJ *a = (Mat_SeqAIJ *)baij->data; in MatMult_MatMultAdd_SeqMAIJ_Template()
200 const PetscInt m = baij->rmap->n; in MatMult_MatMultAdd_SeqMAIJ_Template()
201 const PetscInt nz = a->nz; in MatMult_MatMultAdd_SeqMAIJ_Template()
202 const PetscInt *idx = a->j; in MatMult_MatMultAdd_SeqMAIJ_Template()
203 const PetscInt *ii = a->i; in MatMult_MatMultAdd_SeqMAIJ_Template()
204 const PetscScalar *v = a->a; in MatMult_MatMultAdd_SeqMAIJ_Template()
221 const PetscInt n = ii[i + 1] - jrow; in MatMult_MatMultAdd_SeqMAIJ_Template()
222 // leave a line so clang-format does not align these decls in MatMult_MatMultAdd_SeqMAIJ_Template()
245 PetscCall(PetscLogFlops(2 * N * nz - (mult_add ? 0 : (N * nonzerorow)))); in MatMult_MatMultAdd_SeqMAIJ_Template()
258 const Mat_SeqMAIJ *b = (Mat_SeqMAIJ *)A->data; in MatMultTranspose_MatMultTransposeAdd_SeqMAIJ_Template() local
259 const Mat baij = b->AIJ; in MatMultTranspose_MatMultTransposeAdd_SeqMAIJ_Template()
260 const Mat_SeqAIJ *a = (Mat_SeqAIJ *)baij->data; in MatMultTranspose_MatMultTransposeAdd_SeqMAIJ_Template()
261 const PetscInt m = baij->rmap->n; in MatMultTranspose_MatMultTransposeAdd_SeqMAIJ_Template()
262 const PetscInt nz = a->nz; in MatMultTranspose_MatMultTransposeAdd_SeqMAIJ_Template()
263 const PetscInt *a_j = a->j; in MatMultTranspose_MatMultTransposeAdd_SeqMAIJ_Template()
264 const PetscInt *a_i = a->i; in MatMultTranspose_MatMultTransposeAdd_SeqMAIJ_Template()
265 const PetscScalar *a_a = a->a; in MatMultTranspose_MatMultTransposeAdd_SeqMAIJ_Template()
283 const PetscInt n = a_i[i + 1] - a_ii; in MatMultTranspose_MatMultTransposeAdd_SeqMAIJ_Template()
346 Mat_SeqMAIJ *b = (Mat_SeqMAIJ *)A->data; in MatMult_SeqMAIJ_N() local
347 Mat_SeqAIJ *a = (Mat_SeqAIJ *)b->AIJ->data; in MatMult_SeqMAIJ_N()
350 const PetscInt m = b->AIJ->rmap->n, *idx, *ii; in MatMult_SeqMAIJ_N()
351 PetscInt n, i, jrow, j, dof = b->dof, k; in MatMult_SeqMAIJ_N()
357 idx = a->j; in MatMult_SeqMAIJ_N()
358 v = a->a; in MatMult_SeqMAIJ_N()
359 ii = a->i; in MatMult_SeqMAIJ_N()
363 n = ii[i + 1] - jrow; in MatMult_SeqMAIJ_N()
371 PetscCall(PetscLogFlops(2.0 * dof * a->nz)); in MatMult_SeqMAIJ_N()
379 Mat_SeqMAIJ *b = (Mat_SeqMAIJ *)A->data; in MatMultAdd_SeqMAIJ_N() local
380 Mat_SeqAIJ *a = (Mat_SeqAIJ *)b->AIJ->data; in MatMultAdd_SeqMAIJ_N()
383 const PetscInt m = b->AIJ->rmap->n, *idx, *ii; in MatMultAdd_SeqMAIJ_N()
384 PetscInt n, i, jrow, j, dof = b->dof, k; in MatMultAdd_SeqMAIJ_N()
390 idx = a->j; in MatMultAdd_SeqMAIJ_N()
391 v = a->a; in MatMultAdd_SeqMAIJ_N()
392 ii = a->i; in MatMultAdd_SeqMAIJ_N()
396 n = ii[i + 1] - jrow; in MatMultAdd_SeqMAIJ_N()
404 PetscCall(PetscLogFlops(2.0 * dof * a->nz)); in MatMultAdd_SeqMAIJ_N()
412 Mat_SeqMAIJ *b = (Mat_SeqMAIJ *)A->data; in MatMultTranspose_SeqMAIJ_N() local
413 Mat_SeqAIJ *a = (Mat_SeqAIJ *)b->AIJ->data; in MatMultTranspose_SeqMAIJ_N()
416 const PetscInt m = b->AIJ->rmap->n, *idx, dof = b->dof; in MatMultTranspose_SeqMAIJ_N()
424 idx = PetscSafePointerPlusOffset(a->j, a->i[i]); in MatMultTranspose_SeqMAIJ_N()
425 v = PetscSafePointerPlusOffset(a->a, a->i[i]); in MatMultTranspose_SeqMAIJ_N()
426 n = a->i[i + 1] - a->i[i]; in MatMultTranspose_SeqMAIJ_N()
428 while (n-- > 0) { in MatMultTranspose_SeqMAIJ_N()
434 PetscCall(PetscLogFlops(2.0 * dof * a->nz)); in MatMultTranspose_SeqMAIJ_N()
442 Mat_SeqMAIJ *b = (Mat_SeqMAIJ *)A->data; in MatMultTransposeAdd_SeqMAIJ_N() local
443 Mat_SeqAIJ *a = (Mat_SeqAIJ *)b->AIJ->data; in MatMultTransposeAdd_SeqMAIJ_N()
446 const PetscInt m = b->AIJ->rmap->n, *idx, dof = b->dof; in MatMultTransposeAdd_SeqMAIJ_N()
454 idx = a->j + a->i[i]; in MatMultTransposeAdd_SeqMAIJ_N()
455 v = a->a + a->i[i]; in MatMultTransposeAdd_SeqMAIJ_N()
456 n = a->i[i + 1] - a->i[i]; in MatMultTransposeAdd_SeqMAIJ_N()
458 while (n-- > 0) { in MatMultTransposeAdd_SeqMAIJ_N()
464 PetscCall(PetscLogFlops(2.0 * dof * a->nz)); in MatMultTransposeAdd_SeqMAIJ_N()
472 Mat_MPIMAIJ *b = (Mat_MPIMAIJ *)A->data; in MatMult_MPIMAIJ_dof() local
476 PetscCall(VecScatterBegin(b->ctx, xx, b->w, INSERT_VALUES, SCATTER_FORWARD)); in MatMult_MPIMAIJ_dof()
477 PetscCall((*b->AIJ->ops->mult)(b->AIJ, xx, yy)); in MatMult_MPIMAIJ_dof()
478 PetscCall(VecScatterEnd(b->ctx, xx, b->w, INSERT_VALUES, SCATTER_FORWARD)); in MatMult_MPIMAIJ_dof()
479 PetscCall((*b->OAIJ->ops->multadd)(b->OAIJ, b->w, yy, yy)); in MatMult_MPIMAIJ_dof()
485 Mat_MPIMAIJ *b = (Mat_MPIMAIJ *)A->data; in MatMultTranspose_MPIMAIJ_dof() local
488 PetscCall((*b->OAIJ->ops->multtranspose)(b->OAIJ, xx, b->w)); in MatMultTranspose_MPIMAIJ_dof()
489 PetscCall((*b->AIJ->ops->multtranspose)(b->AIJ, xx, yy)); in MatMultTranspose_MPIMAIJ_dof()
490 PetscCall(VecScatterBegin(b->ctx, b->w, yy, ADD_VALUES, SCATTER_REVERSE)); in MatMultTranspose_MPIMAIJ_dof()
491 PetscCall(VecScatterEnd(b->ctx, b->w, yy, ADD_VALUES, SCATTER_REVERSE)); in MatMultTranspose_MPIMAIJ_dof()
497 Mat_MPIMAIJ *b = (Mat_MPIMAIJ *)A->data; in MatMultAdd_MPIMAIJ_dof() local
501 PetscCall(VecScatterBegin(b->ctx, xx, b->w, INSERT_VALUES, SCATTER_FORWARD)); in MatMultAdd_MPIMAIJ_dof()
502 PetscCall((*b->AIJ->ops->multadd)(b->AIJ, xx, yy, zz)); in MatMultAdd_MPIMAIJ_dof()
503 PetscCall(VecScatterEnd(b->ctx, xx, b->w, INSERT_VALUES, SCATTER_FORWARD)); in MatMultAdd_MPIMAIJ_dof()
504 PetscCall((*b->OAIJ->ops->multadd)(b->OAIJ, b->w, zz, zz)); in MatMultAdd_MPIMAIJ_dof()
510 Mat_MPIMAIJ *b = (Mat_MPIMAIJ *)A->data; in MatMultTransposeAdd_MPIMAIJ_dof() local
513 PetscCall((*b->OAIJ->ops->multtranspose)(b->OAIJ, xx, b->w)); in MatMultTransposeAdd_MPIMAIJ_dof()
514 PetscCall((*b->AIJ->ops->multtransposeadd)(b->AIJ, xx, yy, zz)); in MatMultTransposeAdd_MPIMAIJ_dof()
515 PetscCall(VecScatterBegin(b->ctx, b->w, zz, ADD_VALUES, SCATTER_REVERSE)); in MatMultTransposeAdd_MPIMAIJ_dof()
516 PetscCall(VecScatterEnd(b->ctx, b->w, zz, ADD_VALUES, SCATTER_REVERSE)); in MatMultTransposeAdd_MPIMAIJ_dof()
522 Mat_Product *product = C->product; in MatProductSetFromOptions_SeqAIJ_SeqMAIJ()
525->type == MATPRODUCT_PtAP, PETSC_COMM_SELF, PETSC_ERR_SUP, "Mat Product type %s is not supported f… in MatProductSetFromOptions_SeqAIJ_SeqMAIJ()
526 C->ops->productsymbolic = MatProductSymbolic_PtAP_SeqAIJ_SeqMAIJ; in MatProductSetFromOptions_SeqAIJ_SeqMAIJ()
532 Mat_Product *product = C->product; in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
534 Mat A = product->A, P = product->B; in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
545->type == MATPRODUCT_PtAP, PETSC_COMM_SELF, PETSC_ERR_SUP, "Mat Product type %s is not supported f… in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
549 …PetscCheck(A->rmap->rstart == P->rmap->rstart && A->rmap->rend == P->rmap->rend, PETSC_COMM_SELF, … in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
550 A->rmap->rstart, A->rmap->rend, P->rmap->rstart, P->rmap->rend); in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
551 …PetscCheck(A->cmap->rstart == P->rmap->rstart && A->cmap->rend == P->rmap->rend, PETSC_COMM_SELF, … in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
552 A->cmap->rstart, A->cmap->rend, P->rmap->rstart, P->rmap->rend); in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
555 PetscCall(PetscStrcmp(C->product->alg, "default", &flg)); in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
559 …PetscOptionsBegin(PetscObjectComm((PetscObject)C), ((PetscObject)C)->prefix, "MatProduct_PtAP", "M… in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
560 …PetscCall(PetscOptionsEList("-mat_product_algorithm", "Algorithmic approach", "MatPtAP", algTypes,… in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
564 PetscCall(PetscStrcmp(C->product->alg, "allatonce", &flg)); in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
566 C->ops->productsymbolic = MatProductSymbolic_PtAP_MPIAIJ_MPIMAIJ; in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
570 PetscCall(PetscStrcmp(C->product->alg, "allatonce_merged", &flg)); in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
572 C->ops->productsymbolic = MatProductSymbolic_PtAP_MPIAIJ_MPIMAIJ; in MatProductSetFromOptions_MPIAIJ_MPIMAIJ()
585 /* This routine requires testing -- first draft only */ in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
586 Mat_SeqMAIJ *pp = (Mat_SeqMAIJ *)PP->data; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
587 Mat P = pp->AIJ; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
588 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
589 Mat_SeqAIJ *p = (Mat_SeqAIJ *)P->data; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
590 Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
591 const PetscInt *ai = a->i, *aj = a->j, *pi = p->i, *pj = p->j, *pJ, *pjj; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
592 const PetscInt *ci = c->i, *cj = c->j, *cjj; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
593 const PetscInt am = A->rmap->N, cn = C->cmap->N, cm = C->rmap->N, ppdof = pp->dof; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
595 const MatScalar *aa = a->a, *pa = p->a, *pA, *paj; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
596 MatScalar *ca = c->a, *caj, *apa; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
607 anzi = ai[i + 1] - ai[i]; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
614 pnzj = pi[prow + 1] - pi[prow]; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
620 apjdense[poffset] = -1; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
637 pnzi = pi[prow + 1] - poffset; in MatPtAPNumeric_SeqAIJ_SeqMAIJ()
671 Mat_SeqMAIJ *pp = (Mat_SeqMAIJ *)PP->data; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
672 Mat P = pp->AIJ; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
673 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *p = (Mat_SeqAIJ *)P->data, *c; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
676 …const PetscInt an = A->cmap->N, am = A->rmap->N, pn = P->cmap->N, pm = P->rmap->N, ppdof = pp- in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
679 const PetscInt *pi = p->i, *pj = p->j, *pjj, *ai = a->i, *aj = a->j, *ajj; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
704 ptnzi = pti[i + 1] - pti[i]; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
713 anzj = ai[arow + 1] - ai[arow]; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
717 ptadenserow[ajj[k]] = -1; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
731 pnzj = pi[prow + 1] - pi[prow]; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
737 denserow[pjj[k] * ppdof + pshift] = -1; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
748 …if (current_space->local_remaining < cnzi) PetscCall(PetscFreeSpaceGet(PetscIntSumTruncate(cnzi, c… in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
751 PetscCall(PetscArraycpy(current_space->array, sparserow, cnzi)); in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
753 current_space->array += cnzi; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
754 current_space->local_used += cnzi; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
755 current_space->local_remaining -= cnzi; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
777 PetscCall(MatSetBlockSize(C, pp->dof)); in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
781 c = (Mat_SeqAIJ *)C->data; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
782 c->free_a = PETSC_TRUE; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
783 c->free_ij = PETSC_TRUE; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
784 c->nonew = 0; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
786 C->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqMAIJ; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
787 C->ops->productnumeric = MatProductNumeric_PtAP; in MatPtAPSymbolic_SeqAIJ_SeqMAIJ()
796 Mat_Product *product = C->product; in MatProductSymbolic_PtAP_SeqAIJ_SeqMAIJ()
797 Mat A = product->A, P = product->B; in MatProductSymbolic_PtAP_SeqAIJ_SeqMAIJ()
800 PetscCall(MatPtAPSymbolic_SeqAIJ_SeqMAIJ(A, P, product->fill, C)); in MatProductSymbolic_PtAP_SeqAIJ_SeqMAIJ()
808 Mat_MPIMAIJ *maij = (Mat_MPIMAIJ *)P->data; in MatPtAPNumeric_MPIAIJ_MPIMAIJ_allatonce()
811 PetscCall(MatPtAPNumeric_MPIAIJ_MPIXAIJ_allatonce(A, maij->A, maij->dof, C)); in MatPtAPNumeric_MPIAIJ_MPIMAIJ_allatonce()
819 Mat_MPIMAIJ *maij = (Mat_MPIMAIJ *)P->data; in MatPtAPSymbolic_MPIAIJ_MPIMAIJ_allatonce()
822 PetscCall(MatPtAPSymbolic_MPIAIJ_MPIXAIJ_allatonce(A, maij->A, maij->dof, fill, C)); in MatPtAPSymbolic_MPIAIJ_MPIMAIJ_allatonce()
823 C->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIMAIJ_allatonce; in MatPtAPSymbolic_MPIAIJ_MPIMAIJ_allatonce()
831 Mat_MPIMAIJ *maij = (Mat_MPIMAIJ *)P->data; in MatPtAPNumeric_MPIAIJ_MPIMAIJ_allatonce_merged()
834 PetscCall(MatPtAPNumeric_MPIAIJ_MPIXAIJ_allatonce_merged(A, maij->A, maij->dof, C)); in MatPtAPNumeric_MPIAIJ_MPIMAIJ_allatonce_merged()
842 Mat_MPIMAIJ *maij = (Mat_MPIMAIJ *)P->data; in MatPtAPSymbolic_MPIAIJ_MPIMAIJ_allatonce_merged()
845 PetscCall(MatPtAPSymbolic_MPIAIJ_MPIXAIJ_allatonce_merged(A, maij->A, maij->dof, fill, C)); in MatPtAPSymbolic_MPIAIJ_MPIMAIJ_allatonce_merged()
846 C->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIMAIJ_allatonce_merged; in MatPtAPSymbolic_MPIAIJ_MPIMAIJ_allatonce_merged()
852 Mat_Product *product = C->product; in MatProductSymbolic_PtAP_MPIAIJ_MPIMAIJ()
853 Mat A = product->A, P = product->B; in MatProductSymbolic_PtAP_MPIAIJ_MPIMAIJ()
857 PetscCall(PetscStrcmp(product->alg, "allatonce", &flg)); in MatProductSymbolic_PtAP_MPIAIJ_MPIMAIJ()
859 PetscCall(MatPtAPSymbolic_MPIAIJ_MPIMAIJ_allatonce(A, P, product->fill, C)); in MatProductSymbolic_PtAP_MPIAIJ_MPIMAIJ()
860 C->ops->productnumeric = MatProductNumeric_PtAP; in MatProductSymbolic_PtAP_MPIAIJ_MPIMAIJ()
864 PetscCall(PetscStrcmp(product->alg, "allatonce_merged", &flg)); in MatProductSymbolic_PtAP_MPIAIJ_MPIMAIJ()
866 PetscCall(MatPtAPSymbolic_MPIAIJ_MPIMAIJ_allatonce_merged(A, P, product->fill, C)); in MatProductSymbolic_PtAP_MPIAIJ_MPIMAIJ()
867 C->ops->productnumeric = MatProductNumeric_PtAP; in MatProductSymbolic_PtAP_MPIAIJ_MPIMAIJ()
873 Mat_SeqMAIJ *b = (Mat_SeqMAIJ *)A->data; in MatConvert_SeqMAIJ_SeqAIJ() local
874 Mat a = b->AIJ, B; in MatConvert_SeqMAIJ_SeqAIJ() local
875 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)a->data; in MatConvert_SeqMAIJ_SeqAIJ()
876 PetscInt m, n, i, ncols, *ilen, nmax = 0, *icols, j, k, ii, dof = b->dof; in MatConvert_SeqMAIJ_SeqAIJ()
884 nmax = PetscMax(nmax, aij->ilen[i]); in MatConvert_SeqMAIJ_SeqAIJ()
885 for (j = 0; j < dof; j++) ilen[dof * i + j] = aij->ilen[i]; in MatConvert_SeqMAIJ_SeqAIJ()
887 PetscCall(MatCreate(PETSC_COMM_SELF, &B)); in MatConvert_SeqMAIJ_SeqAIJ()
888 PetscCall(MatSetSizes(B, dof * m, dof * n, dof * m, dof * n)); in MatConvert_SeqMAIJ_SeqAIJ()
889 PetscCall(MatSetType(B, newtype)); in MatConvert_SeqMAIJ_SeqAIJ()
890 PetscCall(MatSeqAIJSetPreallocation(B, 0, ilen)); in MatConvert_SeqMAIJ_SeqAIJ()
898 PetscCall(MatSetValues_SeqAIJ(B, 1, &ii, ncols, icols, vals, INSERT_VALUES)); in MatConvert_SeqMAIJ_SeqAIJ()
904 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); in MatConvert_SeqMAIJ_SeqAIJ()
905 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); in MatConvert_SeqMAIJ_SeqAIJ()
908 PetscCall(MatHeaderReplace(A, &B)); in MatConvert_SeqMAIJ_SeqAIJ()
910 *newmat = B; in MatConvert_SeqMAIJ_SeqAIJ()
919 Mat_MPIMAIJ *maij = (Mat_MPIMAIJ *)A->data; in MatConvert_MPIMAIJ_MPIAIJ()
920 Mat MatAIJ = ((Mat_SeqMAIJ *)maij->AIJ->data)->AIJ, B; in MatConvert_MPIMAIJ_MPIAIJ() local
921 Mat MatOAIJ = ((Mat_SeqMAIJ *)maij->OAIJ->data)->AIJ; in MatConvert_MPIMAIJ_MPIAIJ()
922 Mat_SeqAIJ *AIJ = (Mat_SeqAIJ *)MatAIJ->data; in MatConvert_MPIMAIJ_MPIAIJ()
923 Mat_SeqAIJ *OAIJ = (Mat_SeqAIJ *)MatOAIJ->data; in MatConvert_MPIMAIJ_MPIAIJ()
924 Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)maij->A->data; in MatConvert_MPIMAIJ_MPIAIJ()
925 PetscInt dof = maij->dof, i, j, *dnz = NULL, *onz = NULL, nmax = 0, onmax = 0; in MatConvert_MPIMAIJ_MPIAIJ()
931 PetscCall(PetscMalloc2(A->rmap->n, &dnz, A->rmap->n, &onz)); in MatConvert_MPIMAIJ_MPIAIJ()
932 for (i = 0; i < A->rmap->n / dof; i++) { in MatConvert_MPIMAIJ_MPIAIJ()
933 nmax = PetscMax(nmax, AIJ->ilen[i]); in MatConvert_MPIMAIJ_MPIAIJ()
934 onmax = PetscMax(onmax, OAIJ->ilen[i]); in MatConvert_MPIMAIJ_MPIAIJ()
936 dnz[dof * i + j] = AIJ->ilen[i]; in MatConvert_MPIMAIJ_MPIAIJ()
937 onz[dof * i + j] = OAIJ->ilen[i]; in MatConvert_MPIMAIJ_MPIAIJ()
940 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); in MatConvert_MPIMAIJ_MPIAIJ()
941 PetscCall(MatSetSizes(B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); in MatConvert_MPIMAIJ_MPIAIJ()
942 PetscCall(MatSetType(B, newtype)); in MatConvert_MPIMAIJ_MPIAIJ()
943 PetscCall(MatMPIAIJSetPreallocation(B, 0, dnz, 0, onz)); in MatConvert_MPIMAIJ_MPIAIJ()
944 PetscCall(MatSetBlockSize(B, dof)); in MatConvert_MPIMAIJ_MPIAIJ()
948 rstart = dof * maij->A->rmap->rstart; in MatConvert_MPIMAIJ_MPIAIJ()
949 cstart = dof * maij->A->cmap->rstart; in MatConvert_MPIMAIJ_MPIAIJ()
950 garray = mpiaij->garray; in MatConvert_MPIMAIJ_MPIAIJ()
953 for (i = 0; i < A->rmap->n / dof; i++) { in MatConvert_MPIMAIJ_MPIAIJ()
959 PetscCall(MatSetValues_MPIAIJ(B, 1, &ii, ncols, icols, vals, INSERT_VALUES)); in MatConvert_MPIMAIJ_MPIAIJ()
960 PetscCall(MatSetValues_MPIAIJ(B, 1, &ii, oncols, oicols, ovals, INSERT_VALUES)); in MatConvert_MPIMAIJ_MPIAIJ()
968 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); in MatConvert_MPIMAIJ_MPIAIJ()
969 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); in MatConvert_MPIMAIJ_MPIAIJ()
972 PetscInt refct = ((PetscObject)A)->refct; /* save ((PetscObject)A)->refct */ in MatConvert_MPIMAIJ_MPIAIJ()
973 ((PetscObject)A)->refct = 1; in MatConvert_MPIMAIJ_MPIAIJ()
975 PetscCall(MatHeaderReplace(A, &B)); in MatConvert_MPIMAIJ_MPIAIJ()
977 ((PetscObject)A)->refct = refct; /* restore ((PetscObject)A)->refct */ in MatConvert_MPIMAIJ_MPIAIJ()
979 *newmat = B; in MatConvert_MPIMAIJ_MPIAIJ()
1007 MatCreateMAIJ - Creates a matrix type providing restriction and interpolation
1015 + A - the `MATAIJ` matrix describing the action on blocks
1016 - dof - the block size (number of components per node)
1019 . maij - the new `MATMAIJ` matrix
1028 Mat B; in MatCreateMAIJ() local
1041 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); in MatCreateMAIJ()
1043 PetscCall(MatSetVecType(B, A->defaultvectype)); in MatCreateMAIJ()
1044 … PetscCall(MatSetSizes(B, dof * A->rmap->n, dof * A->cmap->n, dof * A->rmap->N, dof * A->cmap->N)); in MatCreateMAIJ()
1045 PetscCall(PetscLayoutSetBlockSize(B->rmap, dof)); in MatCreateMAIJ()
1046 PetscCall(PetscLayoutSetBlockSize(B->cmap, dof)); in MatCreateMAIJ()
1047 PetscCall(PetscLayoutSetUp(B->rmap)); in MatCreateMAIJ()
1048 PetscCall(PetscLayoutSetUp(B->cmap)); in MatCreateMAIJ()
1050 B->assembled = PETSC_TRUE; in MatCreateMAIJ()
1054 Mat_SeqMAIJ *b; in MatCreateMAIJ() local
1056 PetscCall(MatSetType(B, MATSEQMAIJ)); in MatCreateMAIJ()
1058 B->ops->setup = NULL; in MatCreateMAIJ()
1059 B->ops->destroy = MatDestroy_SeqMAIJ; in MatCreateMAIJ()
1060 B->ops->view = MatView_SeqMAIJ; in MatCreateMAIJ()
1062 b = (Mat_SeqMAIJ *)B->data; in MatCreateMAIJ()
1063 b->dof = dof; in MatCreateMAIJ()
1064 b->AIJ = A; in MatCreateMAIJ()
1067 B->ops->mult = MatMult_SeqMAIJ_2; in MatCreateMAIJ()
1068 B->ops->multadd = MatMultAdd_SeqMAIJ_2; in MatCreateMAIJ()
1069 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_2; in MatCreateMAIJ()
1070 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_2; in MatCreateMAIJ()
1072 B->ops->mult = MatMult_SeqMAIJ_3; in MatCreateMAIJ()
1073 B->ops->multadd = MatMultAdd_SeqMAIJ_3; in MatCreateMAIJ()
1074 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_3; in MatCreateMAIJ()
1075 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_3; in MatCreateMAIJ()
1077 B->ops->mult = MatMult_SeqMAIJ_4; in MatCreateMAIJ()
1078 B->ops->multadd = MatMultAdd_SeqMAIJ_4; in MatCreateMAIJ()
1079 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_4; in MatCreateMAIJ()
1080 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_4; in MatCreateMAIJ()
1082 B->ops->mult = MatMult_SeqMAIJ_5; in MatCreateMAIJ()
1083 B->ops->multadd = MatMultAdd_SeqMAIJ_5; in MatCreateMAIJ()
1084 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_5; in MatCreateMAIJ()
1085 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_5; in MatCreateMAIJ()
1087 B->ops->mult = MatMult_SeqMAIJ_6; in MatCreateMAIJ()
1088 B->ops->multadd = MatMultAdd_SeqMAIJ_6; in MatCreateMAIJ()
1089 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_6; in MatCreateMAIJ()
1090 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_6; in MatCreateMAIJ()
1092 B->ops->mult = MatMult_SeqMAIJ_7; in MatCreateMAIJ()
1093 B->ops->multadd = MatMultAdd_SeqMAIJ_7; in MatCreateMAIJ()
1094 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_7; in MatCreateMAIJ()
1095 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_7; in MatCreateMAIJ()
1097 B->ops->mult = MatMult_SeqMAIJ_8; in MatCreateMAIJ()
1098 B->ops->multadd = MatMultAdd_SeqMAIJ_8; in MatCreateMAIJ()
1099 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_8; in MatCreateMAIJ()
1100 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_8; in MatCreateMAIJ()
1102 B->ops->mult = MatMult_SeqMAIJ_9; in MatCreateMAIJ()
1103 B->ops->multadd = MatMultAdd_SeqMAIJ_9; in MatCreateMAIJ()
1104 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_9; in MatCreateMAIJ()
1105 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_9; in MatCreateMAIJ()
1107 B->ops->mult = MatMult_SeqMAIJ_10; in MatCreateMAIJ()
1108 B->ops->multadd = MatMultAdd_SeqMAIJ_10; in MatCreateMAIJ()
1109 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_10; in MatCreateMAIJ()
1110 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_10; in MatCreateMAIJ()
1112 B->ops->mult = MatMult_SeqMAIJ_11; in MatCreateMAIJ()
1113 B->ops->multadd = MatMultAdd_SeqMAIJ_11; in MatCreateMAIJ()
1114 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_11; in MatCreateMAIJ()
1115 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_11; in MatCreateMAIJ()
1117 B->ops->mult = MatMult_SeqMAIJ_16; in MatCreateMAIJ()
1118 B->ops->multadd = MatMultAdd_SeqMAIJ_16; in MatCreateMAIJ()
1119 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_16; in MatCreateMAIJ()
1120 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_16; in MatCreateMAIJ()
1122 B->ops->mult = MatMult_SeqMAIJ_18; in MatCreateMAIJ()
1123 B->ops->multadd = MatMultAdd_SeqMAIJ_18; in MatCreateMAIJ()
1124 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_18; in MatCreateMAIJ()
1125 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_18; in MatCreateMAIJ()
1127 B->ops->mult = MatMult_SeqMAIJ_N; in MatCreateMAIJ()
1128 B->ops->multadd = MatMultAdd_SeqMAIJ_N; in MatCreateMAIJ()
1129 B->ops->multtranspose = MatMultTranspose_SeqMAIJ_N; in MatCreateMAIJ()
1130 B->ops->multtransposeadd = MatMultTransposeAdd_SeqMAIJ_N; in MatCreateMAIJ()
1133 …PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqmaij_seqaijcusparse_C", MatCon… in MatCreateMAIJ()
1135 …PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_seqmaij_seqaij_C", MatConvert_Seq… in MatCreateMAIJ()
1136 …PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_seqaij_seqmaij_C", … in MatCreateMAIJ()
1138 Mat_MPIAIJ *mpiaij = (Mat_MPIAIJ *)A->data; in MatCreateMAIJ()
1139 Mat_MPIMAIJ *b; in MatCreateMAIJ() local
1143 PetscCall(MatSetType(B, MATMPIMAIJ)); in MatCreateMAIJ()
1145 B->ops->setup = NULL; in MatCreateMAIJ()
1146 B->ops->destroy = MatDestroy_MPIMAIJ; in MatCreateMAIJ()
1147 B->ops->view = MatView_MPIMAIJ; in MatCreateMAIJ()
1149 b = (Mat_MPIMAIJ *)B->data; in MatCreateMAIJ()
1150 b->dof = dof; in MatCreateMAIJ()
1151 b->A = A; in MatCreateMAIJ()
1153 PetscCall(MatCreateMAIJ(mpiaij->A, -dof, &b->AIJ)); in MatCreateMAIJ()
1154 PetscCall(MatCreateMAIJ(mpiaij->B, -dof, &b->OAIJ)); in MatCreateMAIJ()
1156 PetscCall(VecGetSize(mpiaij->lvec, &n)); in MatCreateMAIJ()
1157 PetscCall(VecCreate(PETSC_COMM_SELF, &b->w)); in MatCreateMAIJ()
1158 PetscCall(VecSetSizes(b->w, n * dof, n * dof)); in MatCreateMAIJ()
1159 PetscCall(VecSetBlockSize(b->w, dof)); in MatCreateMAIJ()
1160 PetscCall(VecSetType(b->w, VECSEQ)); in MatCreateMAIJ()
1163 …PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)A), dof, n, mpiaij->garray, PETSC_COPY_VALUES… in MatCreateMAIJ()
1167 …CreateMPIWithArray(PetscObjectComm((PetscObject)A), dof, dof * A->cmap->n, dof * A->cmap->N, NULL,… in MatCreateMAIJ()
1170 PetscCall(VecScatterCreate(gvec, from, b->w, to, &b->ctx)); in MatCreateMAIJ()
1176 B->ops->mult = MatMult_MPIMAIJ_dof; in MatCreateMAIJ()
1177 B->ops->multtranspose = MatMultTranspose_MPIMAIJ_dof; in MatCreateMAIJ()
1178 B->ops->multadd = MatMultAdd_MPIMAIJ_dof; in MatCreateMAIJ()
1179 B->ops->multtransposeadd = MatMultTransposeAdd_MPIMAIJ_dof; in MatCreateMAIJ()
1182 …PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpimaij_mpiaijcusparse_C", MatCon… in MatCreateMAIJ()
1184 …PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatConvert_mpimaij_mpiaij_C", MatConvert_MPI… in MatCreateMAIJ()
1185 …PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatProductSetFromOptions_mpiaij_mpimaij_C", … in MatCreateMAIJ()
1187 B->ops->createsubmatrix = MatCreateSubMatrix_MAIJ; in MatCreateMAIJ()
1188 B->ops->createsubmatrices = MatCreateSubMatrices_MAIJ; in MatCreateMAIJ()
1189 PetscCall(MatSetUp(B)); in MatCreateMAIJ()
1196 if (flg) PetscCall(MatConvert(B, ((PetscObject)A)->type_name, MAT_INPLACE_MATRIX, &B)); in MatCreateMAIJ()
1200 *maij = B; in MatCreateMAIJ()