xref: /petsc/src/mat/impls/nest/matnest.c (revision e9d3347a9a7c60c1d9de7be8d2006267b8246ad9)
1aaa7dc30SBarry Smith #include <../src/mat/impls/nest/matnestimpl.h> /*I   "petscmat.h"   I*/
2b68353e5Sstefano_zampini #include <../src/mat/impls/aij/seq/aij.h>
30c312b8eSJed Brown #include <petscsf.h>
4d8588912SDave May 
5c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat, PetscInt, const IS[], PetscInt, const IS[]);
606a1af2fSStefano Zampini static PetscErrorCode MatCreateVecs_Nest(Mat, Vec *, Vec *);
706a1af2fSStefano Zampini static PetscErrorCode MatReset_Nest(Mat);
806a1af2fSStefano Zampini 
95e3038f0Sstefano_zampini PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat, MatType, MatReuse, Mat *);
10c8883902SJed Brown 
11d8588912SDave May /* private functions */
12d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNestGetSizes_Private(Mat A, PetscInt *m, PetscInt *n, PetscInt *M, PetscInt *N)
13d71ae5a4SJacob Faibussowitsch {
14d8588912SDave May   Mat_Nest *bA = (Mat_Nest *)A->data;
158188e55aSJed Brown   PetscInt  i, j;
16d8588912SDave May 
17d8588912SDave May   PetscFunctionBegin;
188188e55aSJed Brown   *m = *n = *M = *N = 0;
198188e55aSJed Brown   for (i = 0; i < bA->nr; i++) { /* rows */
208188e55aSJed Brown     PetscInt sm, sM;
219566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(bA->isglobal.row[i], &sm));
229566063dSJacob Faibussowitsch     PetscCall(ISGetSize(bA->isglobal.row[i], &sM));
238188e55aSJed Brown     *m += sm;
248188e55aSJed Brown     *M += sM;
25d8588912SDave May   }
268188e55aSJed Brown   for (j = 0; j < bA->nc; j++) { /* cols */
278188e55aSJed Brown     PetscInt sn, sN;
289566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(bA->isglobal.col[j], &sn));
299566063dSJacob Faibussowitsch     PetscCall(ISGetSize(bA->isglobal.col[j], &sN));
308188e55aSJed Brown     *n += sn;
318188e55aSJed Brown     *N += sN;
32d8588912SDave May   }
333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34d8588912SDave May }
35d8588912SDave May 
36d8588912SDave May /* operations */
37d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_Nest(Mat A, Vec x, Vec y)
38d71ae5a4SJacob Faibussowitsch {
39d8588912SDave May   Mat_Nest *bA = (Mat_Nest *)A->data;
40207556f9SJed Brown   Vec      *bx = bA->right, *by = bA->left;
41207556f9SJed Brown   PetscInt  i, j, nr = bA->nr, nc = bA->nc;
42d8588912SDave May 
43d8588912SDave May   PetscFunctionBegin;
449566063dSJacob Faibussowitsch   for (i = 0; i < nr; i++) PetscCall(VecGetSubVector(y, bA->isglobal.row[i], &by[i]));
459566063dSJacob Faibussowitsch   for (i = 0; i < nc; i++) PetscCall(VecGetSubVector(x, bA->isglobal.col[i], &bx[i]));
46207556f9SJed Brown   for (i = 0; i < nr; i++) {
479566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(by[i]));
48207556f9SJed Brown     for (j = 0; j < nc; j++) {
49207556f9SJed Brown       if (!bA->m[i][j]) continue;
50d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
519566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(bA->m[i][j], bx[j], by[i], by[i]));
52d8588912SDave May     }
53d8588912SDave May   }
549566063dSJacob Faibussowitsch   for (i = 0; i < nr; i++) PetscCall(VecRestoreSubVector(y, bA->isglobal.row[i], &by[i]));
559566063dSJacob Faibussowitsch   for (i = 0; i < nc; i++) PetscCall(VecRestoreSubVector(x, bA->isglobal.col[i], &bx[i]));
563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57d8588912SDave May }
58d8588912SDave May 
59d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultAdd_Nest(Mat A, Vec x, Vec y, Vec z)
60d71ae5a4SJacob Faibussowitsch {
619194d70fSJed Brown   Mat_Nest *bA = (Mat_Nest *)A->data;
629194d70fSJed Brown   Vec      *bx = bA->right, *bz = bA->left;
639194d70fSJed Brown   PetscInt  i, j, nr = bA->nr, nc = bA->nc;
649194d70fSJed Brown 
659194d70fSJed Brown   PetscFunctionBegin;
669566063dSJacob Faibussowitsch   for (i = 0; i < nr; i++) PetscCall(VecGetSubVector(z, bA->isglobal.row[i], &bz[i]));
679566063dSJacob Faibussowitsch   for (i = 0; i < nc; i++) PetscCall(VecGetSubVector(x, bA->isglobal.col[i], &bx[i]));
689194d70fSJed Brown   for (i = 0; i < nr; i++) {
699194d70fSJed Brown     if (y != z) {
709194d70fSJed Brown       Vec by;
719566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(y, bA->isglobal.row[i], &by));
729566063dSJacob Faibussowitsch       PetscCall(VecCopy(by, bz[i]));
739566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(y, bA->isglobal.row[i], &by));
749194d70fSJed Brown     }
759194d70fSJed Brown     for (j = 0; j < nc; j++) {
769194d70fSJed Brown       if (!bA->m[i][j]) continue;
779194d70fSJed Brown       /* y[i] <- y[i] + A[i][j] * x[j] */
789566063dSJacob Faibussowitsch       PetscCall(MatMultAdd(bA->m[i][j], bx[j], bz[i], bz[i]));
799194d70fSJed Brown     }
809194d70fSJed Brown   }
819566063dSJacob Faibussowitsch   for (i = 0; i < nr; i++) PetscCall(VecRestoreSubVector(z, bA->isglobal.row[i], &bz[i]));
829566063dSJacob Faibussowitsch   for (i = 0; i < nc; i++) PetscCall(VecRestoreSubVector(x, bA->isglobal.col[i], &bx[i]));
833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
849194d70fSJed Brown }
859194d70fSJed Brown 
8652c5f739Sprj- typedef struct {
8752c5f739Sprj-   Mat         *workC;      /* array of Mat with specific containers depending on the underlying MatMatMult implementation */
8852c5f739Sprj-   PetscScalar *tarray;     /* buffer for storing all temporary products A[i][j] B[j] */
8952c5f739Sprj-   PetscInt    *dm, *dn, k; /* displacements and number of submatrices */
9052c5f739Sprj- } Nest_Dense;
9152c5f739Sprj- 
92a678f235SPierre Jolivet static PetscErrorCode MatProductNumeric_Nest_Dense(Mat C)
93d71ae5a4SJacob Faibussowitsch {
946718818eSStefano Zampini   Mat_Nest          *bA;
9552c5f739Sprj-   Nest_Dense        *contents;
966718818eSStefano Zampini   Mat                viewB, viewC, productB, workC;
9752c5f739Sprj-   const PetscScalar *barray;
9852c5f739Sprj-   PetscScalar       *carray;
996718818eSStefano Zampini   PetscInt           i, j, M, N, nr, nc, ldb, ldc;
1006718818eSStefano Zampini   Mat                A, B;
10152c5f739Sprj- 
10252c5f739Sprj-   PetscFunctionBegin;
1030d6f747bSJacob Faibussowitsch   MatCheckProduct(C, 1);
1046718818eSStefano Zampini   A = C->product->A;
1056718818eSStefano Zampini   B = C->product->B;
1069566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, NULL, &N));
1076718818eSStefano Zampini   if (!N) {
1089566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
1099566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
1103ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1116718818eSStefano Zampini   }
1126718818eSStefano Zampini   contents = (Nest_Dense *)C->product->data;
11328b400f6SJacob Faibussowitsch   PetscCheck(contents, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data empty");
1146718818eSStefano Zampini   bA = (Mat_Nest *)A->data;
1156718818eSStefano Zampini   nr = bA->nr;
1166718818eSStefano Zampini   nc = bA->nc;
1179566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(B, &ldb));
1189566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(C, &ldc));
1199566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(C));
1209566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(B, &barray));
1219566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(C, &carray));
12252c5f739Sprj-   for (i = 0; i < nr; i++) {
1239566063dSJacob Faibussowitsch     PetscCall(ISGetSize(bA->isglobal.row[i], &M));
1249566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), contents->dm[i + 1] - contents->dm[i], PETSC_DECIDE, M, N, carray + contents->dm[i], &viewC));
1259566063dSJacob Faibussowitsch     PetscCall(MatDenseSetLDA(viewC, ldc));
12652c5f739Sprj-     for (j = 0; j < nc; j++) {
12752c5f739Sprj-       if (!bA->m[i][j]) continue;
1289566063dSJacob Faibussowitsch       PetscCall(ISGetSize(bA->isglobal.col[j], &M));
1299566063dSJacob Faibussowitsch       PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), contents->dn[j + 1] - contents->dn[j], PETSC_DECIDE, M, N, (PetscScalar *)(barray + contents->dn[j]), &viewB));
1309566063dSJacob Faibussowitsch       PetscCall(MatDenseSetLDA(viewB, ldb));
1314222ddf1SHong Zhang 
1324222ddf1SHong Zhang       /* MatMatMultNumeric(bA->m[i][j],viewB,contents->workC[i*nc + j]); */
1334222ddf1SHong Zhang       workC             = contents->workC[i * nc + j];
1344222ddf1SHong Zhang       productB          = workC->product->B;
1354222ddf1SHong Zhang       workC->product->B = viewB; /* use newly created dense matrix viewB */
1369566063dSJacob Faibussowitsch       PetscCall(MatProductNumeric(workC));
1379566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&viewB));
1384222ddf1SHong Zhang       workC->product->B = productB; /* resume original B */
1394222ddf1SHong Zhang 
14052c5f739Sprj-       /* C[i] <- workC + C[i] */
1419566063dSJacob Faibussowitsch       PetscCall(MatAXPY(viewC, 1.0, contents->workC[i * nc + j], SAME_NONZERO_PATTERN));
14252c5f739Sprj-     }
1439566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&viewC));
14452c5f739Sprj-   }
1459566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(C, &carray));
1469566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(B, &barray));
1474222ddf1SHong Zhang 
1489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
1499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
1503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15152c5f739Sprj- }
15252c5f739Sprj- 
15366976f2fSJacob Faibussowitsch static PetscErrorCode MatNest_DenseDestroy(void *ctx)
154d71ae5a4SJacob Faibussowitsch {
15552c5f739Sprj-   Nest_Dense *contents = (Nest_Dense *)ctx;
15652c5f739Sprj-   PetscInt    i;
15752c5f739Sprj- 
15852c5f739Sprj-   PetscFunctionBegin;
1599566063dSJacob Faibussowitsch   PetscCall(PetscFree(contents->tarray));
16048a46eb9SPierre Jolivet   for (i = 0; i < contents->k; i++) PetscCall(MatDestroy(contents->workC + i));
1619566063dSJacob Faibussowitsch   PetscCall(PetscFree3(contents->dm, contents->dn, contents->workC));
1629566063dSJacob Faibussowitsch   PetscCall(PetscFree(contents));
1633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16452c5f739Sprj- }
16552c5f739Sprj- 
166a678f235SPierre Jolivet static PetscErrorCode MatProductSymbolic_Nest_Dense(Mat C)
167d71ae5a4SJacob Faibussowitsch {
1686718818eSStefano Zampini   Mat_Nest          *bA;
1696718818eSStefano Zampini   Mat                viewB, workC;
17052c5f739Sprj-   const PetscScalar *barray;
1716718818eSStefano Zampini   PetscInt           i, j, M, N, m, n, nr, nc, maxm = 0, ldb;
1724222ddf1SHong Zhang   Nest_Dense        *contents = NULL;
1736718818eSStefano Zampini   PetscBool          cisdense;
1746718818eSStefano Zampini   Mat                A, B;
1756718818eSStefano Zampini   PetscReal          fill;
17652c5f739Sprj- 
17752c5f739Sprj-   PetscFunctionBegin;
1780d6f747bSJacob Faibussowitsch   MatCheckProduct(C, 1);
17928b400f6SJacob Faibussowitsch   PetscCheck(!C->product->data, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Product data not empty");
1806718818eSStefano Zampini   A    = C->product->A;
1816718818eSStefano Zampini   B    = C->product->B;
1826718818eSStefano Zampini   fill = C->product->fill;
1836718818eSStefano Zampini   bA   = (Mat_Nest *)A->data;
1846718818eSStefano Zampini   nr   = bA->nr;
1856718818eSStefano Zampini   nc   = bA->nc;
1869566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(C, &m, &n));
1879566063dSJacob Faibussowitsch   PetscCall(MatGetSize(C, &M, &N));
1880572eedcSPierre Jolivet   if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) {
1899566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(B, NULL, &n));
1909566063dSJacob Faibussowitsch     PetscCall(MatGetSize(B, NULL, &N));
1919566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(A, &m, NULL));
1929566063dSJacob Faibussowitsch     PetscCall(MatGetSize(A, &M, NULL));
1939566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, m, n, M, N));
1940572eedcSPierre Jolivet   }
1959566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATMPIDENSE, MATSEQDENSECUDA, MATMPIDENSECUDA, ""));
19648a46eb9SPierre Jolivet   if (!cisdense) PetscCall(MatSetType(C, ((PetscObject)B)->type_name));
1979566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
1986718818eSStefano Zampini   if (!N) {
1996718818eSStefano Zampini     C->ops->productnumeric = MatProductNumeric_Nest_Dense;
2003ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
20152c5f739Sprj-   }
20252c5f739Sprj- 
2039566063dSJacob Faibussowitsch   PetscCall(PetscNew(&contents));
2046718818eSStefano Zampini   C->product->data    = contents;
2056718818eSStefano Zampini   C->product->destroy = MatNest_DenseDestroy;
2069566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(nr + 1, &contents->dm, nc + 1, &contents->dn, nr * nc, &contents->workC));
20752c5f739Sprj-   contents->k = nr * nc;
20852c5f739Sprj-   for (i = 0; i < nr; i++) {
2099566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(bA->isglobal.row[i], contents->dm + i + 1));
21052c5f739Sprj-     maxm = PetscMax(maxm, contents->dm[i + 1]);
21152c5f739Sprj-     contents->dm[i + 1] += contents->dm[i];
21252c5f739Sprj-   }
21352c5f739Sprj-   for (i = 0; i < nc; i++) {
2149566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(bA->isglobal.col[i], contents->dn + i + 1));
21552c5f739Sprj-     contents->dn[i + 1] += contents->dn[i];
21652c5f739Sprj-   }
2179566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(maxm * N, &contents->tarray));
2189566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(B, &ldb));
2199566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, NULL, &N));
2209566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(B, &barray));
22152c5f739Sprj-   /* loops are permuted compared to MatMatMultNumeric so that viewB is created only once per column of A */
22252c5f739Sprj-   for (j = 0; j < nc; j++) {
2239566063dSJacob Faibussowitsch     PetscCall(ISGetSize(bA->isglobal.col[j], &M));
2249566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), contents->dn[j + 1] - contents->dn[j], PETSC_DECIDE, M, N, (PetscScalar *)(barray + contents->dn[j]), &viewB));
2259566063dSJacob Faibussowitsch     PetscCall(MatDenseSetLDA(viewB, ldb));
22652c5f739Sprj-     for (i = 0; i < nr; i++) {
22752c5f739Sprj-       if (!bA->m[i][j]) continue;
22852c5f739Sprj-       /* MatMatMultSymbolic may attach a specific container (depending on MatType of bA->m[i][j]) to workC[i][j] */
2294222ddf1SHong Zhang 
2309566063dSJacob Faibussowitsch       PetscCall(MatProductCreate(bA->m[i][j], viewB, NULL, &contents->workC[i * nc + j]));
2314222ddf1SHong Zhang       workC = contents->workC[i * nc + j];
2329566063dSJacob Faibussowitsch       PetscCall(MatProductSetType(workC, MATPRODUCT_AB));
2339566063dSJacob Faibussowitsch       PetscCall(MatProductSetAlgorithm(workC, "default"));
2349566063dSJacob Faibussowitsch       PetscCall(MatProductSetFill(workC, fill));
2359566063dSJacob Faibussowitsch       PetscCall(MatProductSetFromOptions(workC));
2369566063dSJacob Faibussowitsch       PetscCall(MatProductSymbolic(workC));
2374222ddf1SHong Zhang 
2386718818eSStefano Zampini       /* since tarray will be shared by all Mat */
2399566063dSJacob Faibussowitsch       PetscCall(MatSeqDenseSetPreallocation(workC, contents->tarray));
2409566063dSJacob Faibussowitsch       PetscCall(MatMPIDenseSetPreallocation(workC, contents->tarray));
24152c5f739Sprj-     }
2429566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&viewB));
24352c5f739Sprj-   }
2449566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(B, &barray));
24552c5f739Sprj- 
2466718818eSStefano Zampini   C->ops->productnumeric = MatProductNumeric_Nest_Dense;
2473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24852c5f739Sprj- }
24952c5f739Sprj- 
250d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_Nest_Dense_AB(Mat C)
251d71ae5a4SJacob Faibussowitsch {
2524222ddf1SHong Zhang   PetscFunctionBegin;
2536718818eSStefano Zampini   C->ops->productsymbolic = MatProductSymbolic_Nest_Dense;
2543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2554222ddf1SHong Zhang }
2564222ddf1SHong Zhang 
257a678f235SPierre Jolivet static PetscErrorCode MatProductSetFromOptions_Nest_Dense(Mat C)
258d71ae5a4SJacob Faibussowitsch {
2594222ddf1SHong Zhang   Mat_Product *product = C->product;
26052c5f739Sprj- 
26152c5f739Sprj-   PetscFunctionBegin;
26248a46eb9SPierre Jolivet   if (product->type == MATPRODUCT_AB) PetscCall(MatProductSetFromOptions_Nest_Dense_AB(C));
2633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26452c5f739Sprj- }
26552c5f739Sprj- 
266d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_Nest(Mat A, Vec x, Vec y)
267d71ae5a4SJacob Faibussowitsch {
268d8588912SDave May   Mat_Nest *bA = (Mat_Nest *)A->data;
269207556f9SJed Brown   Vec      *bx = bA->left, *by = bA->right;
270207556f9SJed Brown   PetscInt  i, j, nr = bA->nr, nc = bA->nc;
271d8588912SDave May 
272d8588912SDave May   PetscFunctionBegin;
2739566063dSJacob Faibussowitsch   for (i = 0; i < nr; i++) PetscCall(VecGetSubVector(x, bA->isglobal.row[i], &bx[i]));
2749566063dSJacob Faibussowitsch   for (i = 0; i < nc; i++) PetscCall(VecGetSubVector(y, bA->isglobal.col[i], &by[i]));
275207556f9SJed Brown   for (j = 0; j < nc; j++) {
2769566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(by[j]));
277609e31cbSJed Brown     for (i = 0; i < nr; i++) {
2786c75ac25SJed Brown       if (!bA->m[i][j]) continue;
279609e31cbSJed Brown       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
2809566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(bA->m[i][j], bx[i], by[j], by[j]));
281d8588912SDave May     }
282d8588912SDave May   }
2839566063dSJacob Faibussowitsch   for (i = 0; i < nr; i++) PetscCall(VecRestoreSubVector(x, bA->isglobal.row[i], &bx[i]));
2849566063dSJacob Faibussowitsch   for (i = 0; i < nc; i++) PetscCall(VecRestoreSubVector(y, bA->isglobal.col[i], &by[i]));
2853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
286d8588912SDave May }
287d8588912SDave May 
288d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTransposeAdd_Nest(Mat A, Vec x, Vec y, Vec z)
289d71ae5a4SJacob Faibussowitsch {
2909194d70fSJed Brown   Mat_Nest *bA = (Mat_Nest *)A->data;
2919194d70fSJed Brown   Vec      *bx = bA->left, *bz = bA->right;
2929194d70fSJed Brown   PetscInt  i, j, nr = bA->nr, nc = bA->nc;
2939194d70fSJed Brown 
2949194d70fSJed Brown   PetscFunctionBegin;
2959566063dSJacob Faibussowitsch   for (i = 0; i < nr; i++) PetscCall(VecGetSubVector(x, bA->isglobal.row[i], &bx[i]));
2969566063dSJacob Faibussowitsch   for (i = 0; i < nc; i++) PetscCall(VecGetSubVector(z, bA->isglobal.col[i], &bz[i]));
2979194d70fSJed Brown   for (j = 0; j < nc; j++) {
2989194d70fSJed Brown     if (y != z) {
2999194d70fSJed Brown       Vec by;
3009566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(y, bA->isglobal.col[j], &by));
3019566063dSJacob Faibussowitsch       PetscCall(VecCopy(by, bz[j]));
3029566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(y, bA->isglobal.col[j], &by));
3039194d70fSJed Brown     }
3049194d70fSJed Brown     for (i = 0; i < nr; i++) {
3056c75ac25SJed Brown       if (!bA->m[i][j]) continue;
3069194d70fSJed Brown       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
3079566063dSJacob Faibussowitsch       PetscCall(MatMultTransposeAdd(bA->m[i][j], bx[i], bz[j], bz[j]));
3089194d70fSJed Brown     }
3099194d70fSJed Brown   }
3109566063dSJacob Faibussowitsch   for (i = 0; i < nr; i++) PetscCall(VecRestoreSubVector(x, bA->isglobal.row[i], &bx[i]));
3119566063dSJacob Faibussowitsch   for (i = 0; i < nc; i++) PetscCall(VecRestoreSubVector(z, bA->isglobal.col[i], &bz[i]));
3123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3139194d70fSJed Brown }
3149194d70fSJed Brown 
315d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatTranspose_Nest(Mat A, MatReuse reuse, Mat *B)
316d71ae5a4SJacob Faibussowitsch {
317f8170845SAlex Fikl   Mat_Nest *bA = (Mat_Nest *)A->data, *bC;
318f8170845SAlex Fikl   Mat       C;
319f8170845SAlex Fikl   PetscInt  i, j, nr = bA->nr, nc = bA->nc;
320f8170845SAlex Fikl 
321f8170845SAlex Fikl   PetscFunctionBegin;
3227fb60732SBarry Smith   if (reuse == MAT_REUSE_MATRIX) PetscCall(MatTransposeCheckNonzeroState_Private(A, *B));
323aed4548fSBarry Smith   PetscCheck(reuse != MAT_INPLACE_MATRIX || nr == nc, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Square nested matrix only for in-place");
324f8170845SAlex Fikl 
325cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
326f8170845SAlex Fikl     Mat *subs;
327f8170845SAlex Fikl     IS  *is_row, *is_col;
328f8170845SAlex Fikl 
3299566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(nr * nc, &subs));
3309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(nr, &is_row, nc, &is_col));
3319566063dSJacob Faibussowitsch     PetscCall(MatNestGetISs(A, is_row, is_col));
332cf37664fSBarry Smith     if (reuse == MAT_INPLACE_MATRIX) {
333ddeb9bd8SAlex Fikl       for (i = 0; i < nr; i++) {
334ad540459SPierre Jolivet         for (j = 0; j < nc; j++) subs[i + nr * j] = bA->m[i][j];
335ddeb9bd8SAlex Fikl       }
336ddeb9bd8SAlex Fikl     }
337ddeb9bd8SAlex Fikl 
3389566063dSJacob Faibussowitsch     PetscCall(MatCreateNest(PetscObjectComm((PetscObject)A), nc, is_col, nr, is_row, subs, &C));
3399566063dSJacob Faibussowitsch     PetscCall(PetscFree(subs));
3409566063dSJacob Faibussowitsch     PetscCall(PetscFree2(is_row, is_col));
341f8170845SAlex Fikl   } else {
342f8170845SAlex Fikl     C = *B;
343f8170845SAlex Fikl   }
344f8170845SAlex Fikl 
345f8170845SAlex Fikl   bC = (Mat_Nest *)C->data;
346f8170845SAlex Fikl   for (i = 0; i < nr; i++) {
347f8170845SAlex Fikl     for (j = 0; j < nc; j++) {
348f8170845SAlex Fikl       if (bA->m[i][j]) {
3499566063dSJacob Faibussowitsch         PetscCall(MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i])));
350f8170845SAlex Fikl       } else {
351f8170845SAlex Fikl         bC->m[j][i] = NULL;
352f8170845SAlex Fikl       }
353f8170845SAlex Fikl     }
354f8170845SAlex Fikl   }
355f8170845SAlex Fikl 
356cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
357f8170845SAlex Fikl     *B = C;
358f8170845SAlex Fikl   } else {
3599566063dSJacob Faibussowitsch     PetscCall(MatHeaderMerge(A, &C));
360f8170845SAlex Fikl   }
3613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
362f8170845SAlex Fikl }
363f8170845SAlex Fikl 
364d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNestDestroyISList(PetscInt n, IS **list)
365d71ae5a4SJacob Faibussowitsch {
366e2d7f03fSJed Brown   IS      *lst = *list;
367e2d7f03fSJed Brown   PetscInt i;
368e2d7f03fSJed Brown 
369e2d7f03fSJed Brown   PetscFunctionBegin;
3703ba16761SJacob Faibussowitsch   if (!lst) PetscFunctionReturn(PETSC_SUCCESS);
3719371c9d4SSatish Balay   for (i = 0; i < n; i++)
3729371c9d4SSatish Balay     if (lst[i]) PetscCall(ISDestroy(&lst[i]));
3739566063dSJacob Faibussowitsch   PetscCall(PetscFree(lst));
3740298fd71SBarry Smith   *list = NULL;
3753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
376e2d7f03fSJed Brown }
377e2d7f03fSJed Brown 
378d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatReset_Nest(Mat A)
379d71ae5a4SJacob Faibussowitsch {
380d8588912SDave May   Mat_Nest *vs = (Mat_Nest *)A->data;
381d8588912SDave May   PetscInt  i, j;
382d8588912SDave May 
383d8588912SDave May   PetscFunctionBegin;
384d8588912SDave May   /* release the matrices and the place holders */
3859566063dSJacob Faibussowitsch   PetscCall(MatNestDestroyISList(vs->nr, &vs->isglobal.row));
3869566063dSJacob Faibussowitsch   PetscCall(MatNestDestroyISList(vs->nc, &vs->isglobal.col));
3879566063dSJacob Faibussowitsch   PetscCall(MatNestDestroyISList(vs->nr, &vs->islocal.row));
3889566063dSJacob Faibussowitsch   PetscCall(MatNestDestroyISList(vs->nc, &vs->islocal.col));
389d8588912SDave May 
3909566063dSJacob Faibussowitsch   PetscCall(PetscFree(vs->row_len));
3919566063dSJacob Faibussowitsch   PetscCall(PetscFree(vs->col_len));
3929566063dSJacob Faibussowitsch   PetscCall(PetscFree(vs->nnzstate));
393d8588912SDave May 
3949566063dSJacob Faibussowitsch   PetscCall(PetscFree2(vs->left, vs->right));
395207556f9SJed Brown 
396d8588912SDave May   /* release the matrices and the place holders */
397d8588912SDave May   if (vs->m) {
398d8588912SDave May     for (i = 0; i < vs->nr; i++) {
39948a46eb9SPierre Jolivet       for (j = 0; j < vs->nc; j++) PetscCall(MatDestroy(&vs->m[i][j]));
4009566063dSJacob Faibussowitsch       PetscCall(PetscFree(vs->m[i]));
401d8588912SDave May     }
4029566063dSJacob Faibussowitsch     PetscCall(PetscFree(vs->m));
403d8588912SDave May   }
40406a1af2fSStefano Zampini 
40506a1af2fSStefano Zampini   /* restore defaults */
40606a1af2fSStefano Zampini   vs->nr            = 0;
40706a1af2fSStefano Zampini   vs->nc            = 0;
40806a1af2fSStefano Zampini   vs->splitassembly = PETSC_FALSE;
4093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41006a1af2fSStefano Zampini }
41106a1af2fSStefano Zampini 
412d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_Nest(Mat A)
413d71ae5a4SJacob Faibussowitsch {
414362febeeSStefano Zampini   PetscFunctionBegin;
4159566063dSJacob Faibussowitsch   PetscCall(MatReset_Nest(A));
4169566063dSJacob Faibussowitsch   PetscCall(PetscFree(A->data));
4179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetSubMat_C", NULL));
4189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestSetSubMat_C", NULL));
4199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetSubMats_C", NULL));
4209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetSize_C", NULL));
4219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetISs_C", NULL));
4229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetLocalISs_C", NULL));
4239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestSetVecType_C", NULL));
4249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestSetSubMats_C", NULL));
4259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_mpiaij_C", NULL));
4269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_seqaij_C", NULL));
4279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_aij_C", NULL));
4289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_is_C", NULL));
4299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_mpidense_C", NULL));
4309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_seqdense_C", NULL));
4319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_nest_seqdense_C", NULL));
4329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_nest_mpidense_C", NULL));
4339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_nest_dense_C", NULL));
4343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
435d8588912SDave May }
436d8588912SDave May 
437d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMissingDiagonal_Nest(Mat mat, PetscBool *missing, PetscInt *dd)
438d71ae5a4SJacob Faibussowitsch {
439381b8e50SStefano Zampini   Mat_Nest *vs = (Mat_Nest *)mat->data;
440381b8e50SStefano Zampini   PetscInt  i;
441381b8e50SStefano Zampini 
442381b8e50SStefano Zampini   PetscFunctionBegin;
443381b8e50SStefano Zampini   if (dd) *dd = 0;
444381b8e50SStefano Zampini   if (!vs->nr) {
445381b8e50SStefano Zampini     *missing = PETSC_TRUE;
4463ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
447381b8e50SStefano Zampini   }
448381b8e50SStefano Zampini   *missing = PETSC_FALSE;
449381b8e50SStefano Zampini   for (i = 0; i < vs->nr && !(*missing); i++) {
450381b8e50SStefano Zampini     *missing = PETSC_TRUE;
451381b8e50SStefano Zampini     if (vs->m[i][i]) {
4529566063dSJacob Faibussowitsch       PetscCall(MatMissingDiagonal(vs->m[i][i], missing, NULL));
45308401ef6SPierre Jolivet       PetscCheck(!*missing || !dd, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "First missing entry not yet implemented");
454381b8e50SStefano Zampini     }
455381b8e50SStefano Zampini   }
4563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
457381b8e50SStefano Zampini }
458381b8e50SStefano Zampini 
459d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatAssemblyBegin_Nest(Mat A, MatAssemblyType type)
460d71ae5a4SJacob Faibussowitsch {
461d8588912SDave May   Mat_Nest *vs = (Mat_Nest *)A->data;
462d8588912SDave May   PetscInt  i, j;
46306a1af2fSStefano Zampini   PetscBool nnzstate = PETSC_FALSE;
464d8588912SDave May 
465d8588912SDave May   PetscFunctionBegin;
466d8588912SDave May   for (i = 0; i < vs->nr; i++) {
467d8588912SDave May     for (j = 0; j < vs->nc; j++) {
46806a1af2fSStefano Zampini       PetscObjectState subnnzstate = 0;
469e7c19651SJed Brown       if (vs->m[i][j]) {
4709566063dSJacob Faibussowitsch         PetscCall(MatAssemblyBegin(vs->m[i][j], type));
471e7c19651SJed Brown         if (!vs->splitassembly) {
472e7c19651SJed Brown           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
473e7c19651SJed Brown            * sub-block). This could be fixed by adding a flag to Mat so that there was a way to check if a Mat was
474e7c19651SJed Brown            * already performing an assembly, but the result would by more complicated and appears to offer less
475e7c19651SJed Brown            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
476e7c19651SJed Brown            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
477e7c19651SJed Brown            */
4789566063dSJacob Faibussowitsch           PetscCall(MatAssemblyEnd(vs->m[i][j], type));
4799566063dSJacob Faibussowitsch           PetscCall(MatGetNonzeroState(vs->m[i][j], &subnnzstate));
480e7c19651SJed Brown         }
481e7c19651SJed Brown       }
48206a1af2fSStefano Zampini       nnzstate                     = (PetscBool)(nnzstate || vs->nnzstate[i * vs->nc + j] != subnnzstate);
48306a1af2fSStefano Zampini       vs->nnzstate[i * vs->nc + j] = subnnzstate;
484d8588912SDave May     }
485d8588912SDave May   }
48606a1af2fSStefano Zampini   if (nnzstate) A->nonzerostate++;
4873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
488d8588912SDave May }
489d8588912SDave May 
490d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
491d71ae5a4SJacob Faibussowitsch {
492d8588912SDave May   Mat_Nest *vs = (Mat_Nest *)A->data;
493d8588912SDave May   PetscInt  i, j;
494d8588912SDave May 
495d8588912SDave May   PetscFunctionBegin;
496d8588912SDave May   for (i = 0; i < vs->nr; i++) {
497d8588912SDave May     for (j = 0; j < vs->nc; j++) {
498e7c19651SJed Brown       if (vs->m[i][j]) {
49948a46eb9SPierre Jolivet         if (vs->splitassembly) PetscCall(MatAssemblyEnd(vs->m[i][j], type));
500e7c19651SJed Brown       }
501d8588912SDave May     }
502d8588912SDave May   }
5033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
504d8588912SDave May }
505d8588912SDave May 
506d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A, PetscInt row, Mat *B)
507d71ae5a4SJacob Faibussowitsch {
508f349c1fdSJed Brown   Mat_Nest *vs = (Mat_Nest *)A->data;
509f349c1fdSJed Brown   PetscInt  j;
510f349c1fdSJed Brown   Mat       sub;
511d8588912SDave May 
512d8588912SDave May   PetscFunctionBegin;
5130298fd71SBarry Smith   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
514f349c1fdSJed Brown   for (j = 0; !sub && j < vs->nc; j++) sub = vs->m[row][j];
5159566063dSJacob Faibussowitsch   if (sub) PetscCall(MatSetUp(sub)); /* Ensure that the sizes are available */
516f349c1fdSJed Brown   *B = sub;
5173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
518d8588912SDave May }
519d8588912SDave May 
520d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A, PetscInt col, Mat *B)
521d71ae5a4SJacob Faibussowitsch {
522f349c1fdSJed Brown   Mat_Nest *vs = (Mat_Nest *)A->data;
523f349c1fdSJed Brown   PetscInt  i;
524f349c1fdSJed Brown   Mat       sub;
525f349c1fdSJed Brown 
526f349c1fdSJed Brown   PetscFunctionBegin;
5270298fd71SBarry Smith   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
528f349c1fdSJed Brown   for (i = 0; !sub && i < vs->nr; i++) sub = vs->m[i][col];
5299566063dSJacob Faibussowitsch   if (sub) PetscCall(MatSetUp(sub)); /* Ensure that the sizes are available */
530f349c1fdSJed Brown   *B = sub;
5313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
532d8588912SDave May }
533d8588912SDave May 
534d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNestFindISRange(Mat A, PetscInt n, const IS list[], IS is, PetscInt *begin, PetscInt *end)
535d71ae5a4SJacob Faibussowitsch {
53618d228c0SPierre Jolivet   PetscInt  i, j, size, m;
537f349c1fdSJed Brown   PetscBool flg;
53818d228c0SPierre Jolivet   IS        out, concatenate[2];
539f349c1fdSJed Brown 
540f349c1fdSJed Brown   PetscFunctionBegin;
5414f572ea9SToby Isaac   PetscAssertPointer(list, 3);
542f349c1fdSJed Brown   PetscValidHeaderSpecific(is, IS_CLASSID, 4);
54318d228c0SPierre Jolivet   if (begin) {
5444f572ea9SToby Isaac     PetscAssertPointer(begin, 5);
54518d228c0SPierre Jolivet     *begin = -1;
54618d228c0SPierre Jolivet   }
54718d228c0SPierre Jolivet   if (end) {
5484f572ea9SToby Isaac     PetscAssertPointer(end, 6);
54918d228c0SPierre Jolivet     *end = -1;
55018d228c0SPierre Jolivet   }
551f349c1fdSJed Brown   for (i = 0; i < n; i++) {
552207556f9SJed Brown     if (!list[i]) continue;
5539566063dSJacob Faibussowitsch     PetscCall(ISEqualUnsorted(list[i], is, &flg));
554f349c1fdSJed Brown     if (flg) {
55518d228c0SPierre Jolivet       if (begin) *begin = i;
55618d228c0SPierre Jolivet       if (end) *end = i + 1;
5573ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
558f349c1fdSJed Brown     }
559f349c1fdSJed Brown   }
5609566063dSJacob Faibussowitsch   PetscCall(ISGetSize(is, &size));
56118d228c0SPierre Jolivet   for (i = 0; i < n - 1; i++) {
56218d228c0SPierre Jolivet     if (!list[i]) continue;
56318d228c0SPierre Jolivet     m = 0;
5649566063dSJacob Faibussowitsch     PetscCall(ISConcatenate(PetscObjectComm((PetscObject)A), 2, list + i, &out));
5659566063dSJacob Faibussowitsch     PetscCall(ISGetSize(out, &m));
56618d228c0SPierre Jolivet     for (j = i + 2; j < n && m < size; j++) {
56718d228c0SPierre Jolivet       if (list[j]) {
56818d228c0SPierre Jolivet         concatenate[0] = out;
56918d228c0SPierre Jolivet         concatenate[1] = list[j];
5709566063dSJacob Faibussowitsch         PetscCall(ISConcatenate(PetscObjectComm((PetscObject)A), 2, concatenate, &out));
5719566063dSJacob Faibussowitsch         PetscCall(ISDestroy(concatenate));
5729566063dSJacob Faibussowitsch         PetscCall(ISGetSize(out, &m));
57318d228c0SPierre Jolivet       }
57418d228c0SPierre Jolivet     }
57518d228c0SPierre Jolivet     if (m == size) {
5769566063dSJacob Faibussowitsch       PetscCall(ISEqualUnsorted(out, is, &flg));
57718d228c0SPierre Jolivet       if (flg) {
57818d228c0SPierre Jolivet         if (begin) *begin = i;
57918d228c0SPierre Jolivet         if (end) *end = j;
5809566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&out));
5813ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
58218d228c0SPierre Jolivet       }
58318d228c0SPierre Jolivet     }
5849566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&out));
58518d228c0SPierre Jolivet   }
5863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
587f349c1fdSJed Brown }
588f349c1fdSJed Brown 
589d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNestFillEmptyMat_Private(Mat A, PetscInt i, PetscInt j, Mat *B)
590d71ae5a4SJacob Faibussowitsch {
5918188e55aSJed Brown   Mat_Nest *vs = (Mat_Nest *)A->data;
59218d228c0SPierre Jolivet   PetscInt  lr, lc;
59318d228c0SPierre Jolivet 
59418d228c0SPierre Jolivet   PetscFunctionBegin;
5959566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
5969566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(vs->isglobal.row[i], &lr));
5979566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(vs->isglobal.col[j], &lc));
5989566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*B, lr, lc, PETSC_DECIDE, PETSC_DECIDE));
5999566063dSJacob Faibussowitsch   PetscCall(MatSetType(*B, MATAIJ));
6009566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*B, 0, NULL));
6019566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*B, 0, NULL, 0, NULL));
6029566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*B));
6039566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*B, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE));
6049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*B, MAT_FINAL_ASSEMBLY));
6059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*B, MAT_FINAL_ASSEMBLY));
6063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
60718d228c0SPierre Jolivet }
60818d228c0SPierre Jolivet 
609d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNestGetBlock_Private(Mat A, PetscInt rbegin, PetscInt rend, PetscInt cbegin, PetscInt cend, Mat *B)
610d71ae5a4SJacob Faibussowitsch {
61118d228c0SPierre Jolivet   Mat_Nest  *vs = (Mat_Nest *)A->data;
61218d228c0SPierre Jolivet   Mat       *a;
61318d228c0SPierre Jolivet   PetscInt   i, j, k, l, nr = rend - rbegin, nc = cend - cbegin;
6148188e55aSJed Brown   char       keyname[256];
61518d228c0SPierre Jolivet   PetscBool *b;
61618d228c0SPierre Jolivet   PetscBool  flg;
6178188e55aSJed Brown 
6188188e55aSJed Brown   PetscFunctionBegin;
6190298fd71SBarry Smith   *B = NULL;
6209566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(keyname, sizeof(keyname), "NestBlock_%" PetscInt_FMT "-%" PetscInt_FMT "x%" PetscInt_FMT "-%" PetscInt_FMT, rbegin, rend, cbegin, cend));
6219566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)A, keyname, (PetscObject *)B));
6223ba16761SJacob Faibussowitsch   if (*B) PetscFunctionReturn(PETSC_SUCCESS);
6238188e55aSJed Brown 
6249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nr * nc, &a, nr * nc, &b));
62518d228c0SPierre Jolivet   for (i = 0; i < nr; i++) {
62618d228c0SPierre Jolivet     for (j = 0; j < nc; j++) {
62718d228c0SPierre Jolivet       a[i * nc + j] = vs->m[rbegin + i][cbegin + j];
62818d228c0SPierre Jolivet       b[i * nc + j] = PETSC_FALSE;
62918d228c0SPierre Jolivet     }
63018d228c0SPierre Jolivet   }
63118d228c0SPierre Jolivet   if (nc != vs->nc && nr != vs->nr) {
63218d228c0SPierre Jolivet     for (i = 0; i < nr; i++) {
63318d228c0SPierre Jolivet       for (j = 0; j < nc; j++) {
63418d228c0SPierre Jolivet         flg = PETSC_FALSE;
63518d228c0SPierre Jolivet         for (k = 0; (k < nr && !flg); k++) {
63618d228c0SPierre Jolivet           if (a[j + k * nc]) flg = PETSC_TRUE;
63718d228c0SPierre Jolivet         }
63818d228c0SPierre Jolivet         if (flg) {
63918d228c0SPierre Jolivet           flg = PETSC_FALSE;
64018d228c0SPierre Jolivet           for (l = 0; (l < nc && !flg); l++) {
64118d228c0SPierre Jolivet             if (a[i * nc + l]) flg = PETSC_TRUE;
64218d228c0SPierre Jolivet           }
64318d228c0SPierre Jolivet         }
64418d228c0SPierre Jolivet         if (!flg) {
64518d228c0SPierre Jolivet           b[i * nc + j] = PETSC_TRUE;
6469566063dSJacob Faibussowitsch           PetscCall(MatNestFillEmptyMat_Private(A, rbegin + i, cbegin + j, a + i * nc + j));
64718d228c0SPierre Jolivet         }
64818d228c0SPierre Jolivet       }
64918d228c0SPierre Jolivet     }
65018d228c0SPierre Jolivet   }
6519566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PetscObjectComm((PetscObject)A), nr, nr != vs->nr ? NULL : vs->isglobal.row, nc, nc != vs->nc ? NULL : vs->isglobal.col, a, B));
65218d228c0SPierre Jolivet   for (i = 0; i < nr; i++) {
65318d228c0SPierre Jolivet     for (j = 0; j < nc; j++) {
65448a46eb9SPierre Jolivet       if (b[i * nc + j]) PetscCall(MatDestroy(a + i * nc + j));
65518d228c0SPierre Jolivet     }
65618d228c0SPierre Jolivet   }
6579566063dSJacob Faibussowitsch   PetscCall(PetscFree2(a, b));
6588188e55aSJed Brown   (*B)->assembled = A->assembled;
6599566063dSJacob Faibussowitsch   PetscCall(PetscObjectCompose((PetscObject)A, keyname, (PetscObject)*B));
6609566063dSJacob Faibussowitsch   PetscCall(PetscObjectDereference((PetscObject)*B)); /* Leave the only remaining reference in the composition */
6613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6628188e55aSJed Brown }
6638188e55aSJed Brown 
664d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNestFindSubMat(Mat A, struct MatNestISPair *is, IS isrow, IS iscol, Mat *B)
665d71ae5a4SJacob Faibussowitsch {
666f349c1fdSJed Brown   Mat_Nest *vs = (Mat_Nest *)A->data;
66718d228c0SPierre Jolivet   PetscInt  rbegin, rend, cbegin, cend;
668f349c1fdSJed Brown 
669f349c1fdSJed Brown   PetscFunctionBegin;
6709566063dSJacob Faibussowitsch   PetscCall(MatNestFindISRange(A, vs->nr, is->row, isrow, &rbegin, &rend));
6719566063dSJacob Faibussowitsch   PetscCall(MatNestFindISRange(A, vs->nc, is->col, iscol, &cbegin, &cend));
67218d228c0SPierre Jolivet   if (rend == rbegin + 1 && cend == cbegin + 1) {
67348a46eb9SPierre Jolivet     if (!vs->m[rbegin][cbegin]) PetscCall(MatNestFillEmptyMat_Private(A, rbegin, cbegin, vs->m[rbegin] + cbegin));
67418d228c0SPierre Jolivet     *B = vs->m[rbegin][cbegin];
67518d228c0SPierre Jolivet   } else if (rbegin != -1 && cbegin != -1) {
6769566063dSJacob Faibussowitsch     PetscCall(MatNestGetBlock_Private(A, rbegin, rend, cbegin, cend, B));
67718d228c0SPierre Jolivet   } else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Could not find index set");
6783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
679f349c1fdSJed Brown }
680f349c1fdSJed Brown 
68106a1af2fSStefano Zampini /*
68206a1af2fSStefano Zampini    TODO: This does not actually returns a submatrix we can modify
68306a1af2fSStefano Zampini */
684d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateSubMatrix_Nest(Mat A, IS isrow, IS iscol, MatReuse reuse, Mat *B)
685d71ae5a4SJacob Faibussowitsch {
686f349c1fdSJed Brown   Mat_Nest *vs = (Mat_Nest *)A->data;
687f349c1fdSJed Brown   Mat       sub;
688f349c1fdSJed Brown 
689f349c1fdSJed Brown   PetscFunctionBegin;
6909566063dSJacob Faibussowitsch   PetscCall(MatNestFindSubMat(A, &vs->isglobal, isrow, iscol, &sub));
691f349c1fdSJed Brown   switch (reuse) {
692f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
6939566063dSJacob Faibussowitsch     if (sub) PetscCall(PetscObjectReference((PetscObject)sub));
694f349c1fdSJed Brown     *B = sub;
695f349c1fdSJed Brown     break;
696d71ae5a4SJacob Faibussowitsch   case MAT_REUSE_MATRIX:
697d71ae5a4SJacob Faibussowitsch     PetscCheck(sub == *B, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Submatrix was not used before in this call");
698d71ae5a4SJacob Faibussowitsch     break;
699d71ae5a4SJacob Faibussowitsch   case MAT_IGNORE_MATRIX: /* Nothing to do */
700d71ae5a4SJacob Faibussowitsch     break;
701d71ae5a4SJacob Faibussowitsch   case MAT_INPLACE_MATRIX: /* Nothing to do */
702d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MAT_INPLACE_MATRIX is not supported yet");
703f349c1fdSJed Brown   }
7043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
705f349c1fdSJed Brown }
706f349c1fdSJed Brown 
70766976f2fSJacob Faibussowitsch static PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A, IS isrow, IS iscol, Mat *B)
708d71ae5a4SJacob Faibussowitsch {
709f349c1fdSJed Brown   Mat_Nest *vs = (Mat_Nest *)A->data;
710f349c1fdSJed Brown   Mat       sub;
711f349c1fdSJed Brown 
712f349c1fdSJed Brown   PetscFunctionBegin;
7139566063dSJacob Faibussowitsch   PetscCall(MatNestFindSubMat(A, &vs->islocal, isrow, iscol, &sub));
714f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
7159566063dSJacob Faibussowitsch   if (sub) PetscCall(PetscObjectReference((PetscObject)sub));
716f349c1fdSJed Brown   *B = sub;
7173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
718d8588912SDave May }
719d8588912SDave May 
720d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A, IS isrow, IS iscol, Mat *B)
721d71ae5a4SJacob Faibussowitsch {
722f349c1fdSJed Brown   Mat_Nest *vs = (Mat_Nest *)A->data;
723f349c1fdSJed Brown   Mat       sub;
724d8588912SDave May 
725d8588912SDave May   PetscFunctionBegin;
7269566063dSJacob Faibussowitsch   PetscCall(MatNestFindSubMat(A, &vs->islocal, isrow, iscol, &sub));
72708401ef6SPierre Jolivet   PetscCheck(*B == sub, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Local submatrix has not been gotten");
728f349c1fdSJed Brown   if (sub) {
729aed4548fSBarry Smith     PetscCheck(((PetscObject)sub)->refct > 1, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Local submatrix has had reference count decremented too many times");
7309566063dSJacob Faibussowitsch     PetscCall(MatDestroy(B));
731d8588912SDave May   }
7323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
733d8588912SDave May }
734d8588912SDave May 
735d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetDiagonal_Nest(Mat A, Vec v)
736d71ae5a4SJacob Faibussowitsch {
7377874fa86SDave May   Mat_Nest *bA = (Mat_Nest *)A->data;
7387874fa86SDave May   PetscInt  i;
7397874fa86SDave May 
7407874fa86SDave May   PetscFunctionBegin;
7417874fa86SDave May   for (i = 0; i < bA->nr; i++) {
742429bac76SJed Brown     Vec bv;
7439566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(v, bA->isglobal.row[i], &bv));
7447874fa86SDave May     if (bA->m[i][i]) {
7459566063dSJacob Faibussowitsch       PetscCall(MatGetDiagonal(bA->m[i][i], bv));
7467874fa86SDave May     } else {
7479566063dSJacob Faibussowitsch       PetscCall(VecSet(bv, 0.0));
7487874fa86SDave May     }
7499566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(v, bA->isglobal.row[i], &bv));
7507874fa86SDave May   }
7513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7527874fa86SDave May }
7537874fa86SDave May 
754d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDiagonalScale_Nest(Mat A, Vec l, Vec r)
755d71ae5a4SJacob Faibussowitsch {
7567874fa86SDave May   Mat_Nest *bA = (Mat_Nest *)A->data;
757429bac76SJed Brown   Vec       bl, *br;
7587874fa86SDave May   PetscInt  i, j;
7597874fa86SDave May 
7607874fa86SDave May   PetscFunctionBegin;
7619566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(bA->nc, &br));
7622e6472ebSElliott Sales de Andrade   if (r) {
7639566063dSJacob Faibussowitsch     for (j = 0; j < bA->nc; j++) PetscCall(VecGetSubVector(r, bA->isglobal.col[j], &br[j]));
7642e6472ebSElliott Sales de Andrade   }
7652e6472ebSElliott Sales de Andrade   bl = NULL;
7667874fa86SDave May   for (i = 0; i < bA->nr; i++) {
76748a46eb9SPierre Jolivet     if (l) PetscCall(VecGetSubVector(l, bA->isglobal.row[i], &bl));
7687874fa86SDave May     for (j = 0; j < bA->nc; j++) {
76948a46eb9SPierre Jolivet       if (bA->m[i][j]) PetscCall(MatDiagonalScale(bA->m[i][j], bl, br[j]));
7707874fa86SDave May     }
77148a46eb9SPierre Jolivet     if (l) PetscCall(VecRestoreSubVector(l, bA->isglobal.row[i], &bl));
7722e6472ebSElliott Sales de Andrade   }
7732e6472ebSElliott Sales de Andrade   if (r) {
7749566063dSJacob Faibussowitsch     for (j = 0; j < bA->nc; j++) PetscCall(VecRestoreSubVector(r, bA->isglobal.col[j], &br[j]));
7752e6472ebSElliott Sales de Andrade   }
7769566063dSJacob Faibussowitsch   PetscCall(PetscFree(br));
7773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7787874fa86SDave May }
7797874fa86SDave May 
780d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_Nest(Mat A, PetscScalar a)
781d71ae5a4SJacob Faibussowitsch {
782a061e289SJed Brown   Mat_Nest *bA = (Mat_Nest *)A->data;
783a061e289SJed Brown   PetscInt  i, j;
784a061e289SJed Brown 
785a061e289SJed Brown   PetscFunctionBegin;
786a061e289SJed Brown   for (i = 0; i < bA->nr; i++) {
787a061e289SJed Brown     for (j = 0; j < bA->nc; j++) {
78848a46eb9SPierre Jolivet       if (bA->m[i][j]) PetscCall(MatScale(bA->m[i][j], a));
789a061e289SJed Brown     }
790a061e289SJed Brown   }
7913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
792a061e289SJed Brown }
793a061e289SJed Brown 
794d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatShift_Nest(Mat A, PetscScalar a)
795d71ae5a4SJacob Faibussowitsch {
796a061e289SJed Brown   Mat_Nest *bA = (Mat_Nest *)A->data;
797a061e289SJed Brown   PetscInt  i;
79806a1af2fSStefano Zampini   PetscBool nnzstate = PETSC_FALSE;
799a061e289SJed Brown 
800a061e289SJed Brown   PetscFunctionBegin;
801a061e289SJed Brown   for (i = 0; i < bA->nr; i++) {
80206a1af2fSStefano Zampini     PetscObjectState subnnzstate = 0;
80308401ef6SPierre Jolivet     PetscCheck(bA->m[i][i], PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "No support for shifting an empty diagonal block, insert a matrix in block (%" PetscInt_FMT ",%" PetscInt_FMT ")", i, i);
8049566063dSJacob Faibussowitsch     PetscCall(MatShift(bA->m[i][i], a));
8059566063dSJacob Faibussowitsch     PetscCall(MatGetNonzeroState(bA->m[i][i], &subnnzstate));
80606a1af2fSStefano Zampini     nnzstate                     = (PetscBool)(nnzstate || bA->nnzstate[i * bA->nc + i] != subnnzstate);
80706a1af2fSStefano Zampini     bA->nnzstate[i * bA->nc + i] = subnnzstate;
808a061e289SJed Brown   }
80906a1af2fSStefano Zampini   if (nnzstate) A->nonzerostate++;
8103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
811a061e289SJed Brown }
812a061e289SJed Brown 
813d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDiagonalSet_Nest(Mat A, Vec D, InsertMode is)
814d71ae5a4SJacob Faibussowitsch {
81513135bc6SAlex Fikl   Mat_Nest *bA = (Mat_Nest *)A->data;
81613135bc6SAlex Fikl   PetscInt  i;
81706a1af2fSStefano Zampini   PetscBool nnzstate = PETSC_FALSE;
81813135bc6SAlex Fikl 
81913135bc6SAlex Fikl   PetscFunctionBegin;
82013135bc6SAlex Fikl   for (i = 0; i < bA->nr; i++) {
82106a1af2fSStefano Zampini     PetscObjectState subnnzstate = 0;
82213135bc6SAlex Fikl     Vec              bv;
8239566063dSJacob Faibussowitsch     PetscCall(VecGetSubVector(D, bA->isglobal.row[i], &bv));
82413135bc6SAlex Fikl     if (bA->m[i][i]) {
8259566063dSJacob Faibussowitsch       PetscCall(MatDiagonalSet(bA->m[i][i], bv, is));
8269566063dSJacob Faibussowitsch       PetscCall(MatGetNonzeroState(bA->m[i][i], &subnnzstate));
82713135bc6SAlex Fikl     }
8289566063dSJacob Faibussowitsch     PetscCall(VecRestoreSubVector(D, bA->isglobal.row[i], &bv));
82906a1af2fSStefano Zampini     nnzstate                     = (PetscBool)(nnzstate || bA->nnzstate[i * bA->nc + i] != subnnzstate);
83006a1af2fSStefano Zampini     bA->nnzstate[i * bA->nc + i] = subnnzstate;
83113135bc6SAlex Fikl   }
83206a1af2fSStefano Zampini   if (nnzstate) A->nonzerostate++;
8333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83413135bc6SAlex Fikl }
83513135bc6SAlex Fikl 
836d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetRandom_Nest(Mat A, PetscRandom rctx)
837d71ae5a4SJacob Faibussowitsch {
838f8170845SAlex Fikl   Mat_Nest *bA = (Mat_Nest *)A->data;
839f8170845SAlex Fikl   PetscInt  i, j;
840f8170845SAlex Fikl 
841f8170845SAlex Fikl   PetscFunctionBegin;
842f8170845SAlex Fikl   for (i = 0; i < bA->nr; i++) {
843f8170845SAlex Fikl     for (j = 0; j < bA->nc; j++) {
84448a46eb9SPierre Jolivet       if (bA->m[i][j]) PetscCall(MatSetRandom(bA->m[i][j], rctx));
845f8170845SAlex Fikl     }
846f8170845SAlex Fikl   }
8473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
848f8170845SAlex Fikl }
849f8170845SAlex Fikl 
850d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCreateVecs_Nest(Mat A, Vec *right, Vec *left)
851d71ae5a4SJacob Faibussowitsch {
852d8588912SDave May   Mat_Nest *bA = (Mat_Nest *)A->data;
853d8588912SDave May   Vec      *L, *R;
854d8588912SDave May   MPI_Comm  comm;
855d8588912SDave May   PetscInt  i, j;
856d8588912SDave May 
857d8588912SDave May   PetscFunctionBegin;
8589566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
859d8588912SDave May   if (right) {
860d8588912SDave May     /* allocate R */
8619566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bA->nc, &R));
862d8588912SDave May     /* Create the right vectors */
863d8588912SDave May     for (j = 0; j < bA->nc; j++) {
864d8588912SDave May       for (i = 0; i < bA->nr; i++) {
865d8588912SDave May         if (bA->m[i][j]) {
8669566063dSJacob Faibussowitsch           PetscCall(MatCreateVecs(bA->m[i][j], &R[j], NULL));
867d8588912SDave May           break;
868d8588912SDave May         }
869d8588912SDave May       }
87008401ef6SPierre Jolivet       PetscCheck(i != bA->nr, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
871d8588912SDave May     }
8729566063dSJacob Faibussowitsch     PetscCall(VecCreateNest(comm, bA->nc, bA->isglobal.col, R, right));
873d8588912SDave May     /* hand back control to the nest vector */
87448a46eb9SPierre Jolivet     for (j = 0; j < bA->nc; j++) PetscCall(VecDestroy(&R[j]));
8759566063dSJacob Faibussowitsch     PetscCall(PetscFree(R));
876d8588912SDave May   }
877d8588912SDave May 
878d8588912SDave May   if (left) {
879d8588912SDave May     /* allocate L */
8809566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(bA->nr, &L));
881d8588912SDave May     /* Create the left vectors */
882d8588912SDave May     for (i = 0; i < bA->nr; i++) {
883d8588912SDave May       for (j = 0; j < bA->nc; j++) {
884d8588912SDave May         if (bA->m[i][j]) {
8859566063dSJacob Faibussowitsch           PetscCall(MatCreateVecs(bA->m[i][j], NULL, &L[i]));
886d8588912SDave May           break;
887d8588912SDave May         }
888d8588912SDave May       }
88908401ef6SPierre Jolivet       PetscCheck(j != bA->nc, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
890d8588912SDave May     }
891d8588912SDave May 
8929566063dSJacob Faibussowitsch     PetscCall(VecCreateNest(comm, bA->nr, bA->isglobal.row, L, left));
89348a46eb9SPierre Jolivet     for (i = 0; i < bA->nr; i++) PetscCall(VecDestroy(&L[i]));
894d8588912SDave May 
8959566063dSJacob Faibussowitsch     PetscCall(PetscFree(L));
896d8588912SDave May   }
8973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
898d8588912SDave May }
899d8588912SDave May 
900d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Nest(Mat A, PetscViewer viewer)
901d71ae5a4SJacob Faibussowitsch {
902d8588912SDave May   Mat_Nest *bA = (Mat_Nest *)A->data;
90329e60adbSStefano Zampini   PetscBool isascii, viewSub = PETSC_FALSE;
904d8588912SDave May   PetscInt  i, j;
905d8588912SDave May 
906d8588912SDave May   PetscFunctionBegin;
9079566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
908d8588912SDave May   if (isascii) {
9099566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_view_nest_sub", &viewSub, NULL));
9109566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix object:\n"));
9119566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
9129566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "type=nest, rows=%" PetscInt_FMT ", cols=%" PetscInt_FMT "\n", bA->nr, bA->nc));
913d8588912SDave May 
9149566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "MatNest structure:\n"));
915d8588912SDave May     for (i = 0; i < bA->nr; i++) {
916d8588912SDave May       for (j = 0; j < bA->nc; j++) {
91719fd82e9SBarry Smith         MatType   type;
918270f95d7SJed Brown         char      name[256] = "", prefix[256] = "";
919d8588912SDave May         PetscInt  NR, NC;
920d8588912SDave May         PetscBool isNest = PETSC_FALSE;
921d8588912SDave May 
922d8588912SDave May         if (!bA->m[i][j]) {
9239566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ",%" PetscInt_FMT ") : NULL\n", i, j));
924d8588912SDave May           continue;
925d8588912SDave May         }
9269566063dSJacob Faibussowitsch         PetscCall(MatGetSize(bA->m[i][j], &NR, &NC));
9279566063dSJacob Faibussowitsch         PetscCall(MatGetType(bA->m[i][j], &type));
9289566063dSJacob Faibussowitsch         if (((PetscObject)bA->m[i][j])->name) PetscCall(PetscSNPrintf(name, sizeof(name), "name=\"%s\", ", ((PetscObject)bA->m[i][j])->name));
9299566063dSJacob Faibussowitsch         if (((PetscObject)bA->m[i][j])->prefix) PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "prefix=\"%s\", ", ((PetscObject)bA->m[i][j])->prefix));
9309566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)bA->m[i][j], MATNEST, &isNest));
931d8588912SDave May 
9329566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "(%" PetscInt_FMT ",%" PetscInt_FMT ") : %s%stype=%s, rows=%" PetscInt_FMT ", cols=%" PetscInt_FMT "\n", i, j, name, prefix, type, NR, NC));
933d8588912SDave May 
93429e60adbSStefano Zampini         if (isNest || viewSub) {
9359566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer)); /* push1 */
9369566063dSJacob Faibussowitsch           PetscCall(MatView(bA->m[i][j], viewer));
9379566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer)); /* pop1 */
938d8588912SDave May         }
939d8588912SDave May       }
940d8588912SDave May     }
9419566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer)); /* pop0 */
942d8588912SDave May   }
9433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
944d8588912SDave May }
945d8588912SDave May 
946d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatZeroEntries_Nest(Mat A)
947d71ae5a4SJacob Faibussowitsch {
948d8588912SDave May   Mat_Nest *bA = (Mat_Nest *)A->data;
949d8588912SDave May   PetscInt  i, j;
950d8588912SDave May 
951d8588912SDave May   PetscFunctionBegin;
952d8588912SDave May   for (i = 0; i < bA->nr; i++) {
953d8588912SDave May     for (j = 0; j < bA->nc; j++) {
954d8588912SDave May       if (!bA->m[i][j]) continue;
9559566063dSJacob Faibussowitsch       PetscCall(MatZeroEntries(bA->m[i][j]));
956d8588912SDave May     }
957d8588912SDave May   }
9583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
959d8588912SDave May }
960d8588912SDave May 
961d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCopy_Nest(Mat A, Mat B, MatStructure str)
962d71ae5a4SJacob Faibussowitsch {
963c222c20dSDavid Ham   Mat_Nest *bA = (Mat_Nest *)A->data, *bB = (Mat_Nest *)B->data;
964c222c20dSDavid Ham   PetscInt  i, j, nr = bA->nr, nc = bA->nc;
96506a1af2fSStefano Zampini   PetscBool nnzstate = PETSC_FALSE;
966c222c20dSDavid Ham 
967c222c20dSDavid Ham   PetscFunctionBegin;
968aed4548fSBarry Smith   PetscCheck(nr == bB->nr && nc == bB->nc, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Cannot copy a Mat_Nest of block size (%" PetscInt_FMT ",%" PetscInt_FMT ") to a Mat_Nest of block size (%" PetscInt_FMT ",%" PetscInt_FMT ")", bB->nr, bB->nc, nr, nc);
969c222c20dSDavid Ham   for (i = 0; i < nr; i++) {
970c222c20dSDavid Ham     for (j = 0; j < nc; j++) {
97106a1af2fSStefano Zampini       PetscObjectState subnnzstate = 0;
97246a2b97cSJed Brown       if (bA->m[i][j] && bB->m[i][j]) {
9739566063dSJacob Faibussowitsch         PetscCall(MatCopy(bA->m[i][j], bB->m[i][j], str));
97408401ef6SPierre Jolivet       } else PetscCheck(!bA->m[i][j] && !bB->m[i][j], PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_INCOMP, "Matrix block does not exist at %" PetscInt_FMT ",%" PetscInt_FMT, i, j);
9759566063dSJacob Faibussowitsch       PetscCall(MatGetNonzeroState(bB->m[i][j], &subnnzstate));
97606a1af2fSStefano Zampini       nnzstate                 = (PetscBool)(nnzstate || bB->nnzstate[i * nc + j] != subnnzstate);
97706a1af2fSStefano Zampini       bB->nnzstate[i * nc + j] = subnnzstate;
978c222c20dSDavid Ham     }
979c222c20dSDavid Ham   }
98006a1af2fSStefano Zampini   if (nnzstate) B->nonzerostate++;
9813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
982c222c20dSDavid Ham }
983c222c20dSDavid Ham 
984d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatAXPY_Nest(Mat Y, PetscScalar a, Mat X, MatStructure str)
985d71ae5a4SJacob Faibussowitsch {
9866e76ffeaSPierre Jolivet   Mat_Nest *bY = (Mat_Nest *)Y->data, *bX = (Mat_Nest *)X->data;
9876e76ffeaSPierre Jolivet   PetscInt  i, j, nr = bY->nr, nc = bY->nc;
98806a1af2fSStefano Zampini   PetscBool nnzstate = PETSC_FALSE;
9896e76ffeaSPierre Jolivet 
9906e76ffeaSPierre Jolivet   PetscFunctionBegin;
991aed4548fSBarry Smith   PetscCheck(nr == bX->nr && nc == bX->nc, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_INCOMP, "Cannot AXPY a MatNest of block size (%" PetscInt_FMT ",%" PetscInt_FMT ") with a MatNest of block size (%" PetscInt_FMT ",%" PetscInt_FMT ")", bX->nr, bX->nc, nr, nc);
9926e76ffeaSPierre Jolivet   for (i = 0; i < nr; i++) {
9936e76ffeaSPierre Jolivet     for (j = 0; j < nc; j++) {
99406a1af2fSStefano Zampini       PetscObjectState subnnzstate = 0;
9956e76ffeaSPierre Jolivet       if (bY->m[i][j] && bX->m[i][j]) {
9969566063dSJacob Faibussowitsch         PetscCall(MatAXPY(bY->m[i][j], a, bX->m[i][j], str));
997c066aebcSStefano Zampini       } else if (bX->m[i][j]) {
998c066aebcSStefano Zampini         Mat M;
999c066aebcSStefano Zampini 
1000e75569e9SPierre Jolivet         PetscCheck(str == DIFFERENT_NONZERO_PATTERN || str == UNKNOWN_NONZERO_PATTERN, PetscObjectComm((PetscObject)Y), PETSC_ERR_ARG_INCOMP, "Matrix block does not exist at %" PetscInt_FMT ",%" PetscInt_FMT ". Use DIFFERENT_NONZERO_PATTERN or UNKNOWN_NONZERO_PATTERN", i, j);
10019566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(bX->m[i][j], MAT_COPY_VALUES, &M));
10029566063dSJacob Faibussowitsch         PetscCall(MatNestSetSubMat(Y, i, j, M));
10039566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&M));
1004c066aebcSStefano Zampini       }
10059566063dSJacob Faibussowitsch       if (bY->m[i][j]) PetscCall(MatGetNonzeroState(bY->m[i][j], &subnnzstate));
100606a1af2fSStefano Zampini       nnzstate                 = (PetscBool)(nnzstate || bY->nnzstate[i * nc + j] != subnnzstate);
100706a1af2fSStefano Zampini       bY->nnzstate[i * nc + j] = subnnzstate;
10086e76ffeaSPierre Jolivet     }
10096e76ffeaSPierre Jolivet   }
101006a1af2fSStefano Zampini   if (nnzstate) Y->nonzerostate++;
10113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10126e76ffeaSPierre Jolivet }
10136e76ffeaSPierre Jolivet 
1014d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDuplicate_Nest(Mat A, MatDuplicateOption op, Mat *B)
1015d71ae5a4SJacob Faibussowitsch {
1016d8588912SDave May   Mat_Nest *bA = (Mat_Nest *)A->data;
1017841e96a3SJed Brown   Mat      *b;
1018841e96a3SJed Brown   PetscInt  i, j, nr = bA->nr, nc = bA->nc;
1019d8588912SDave May 
1020d8588912SDave May   PetscFunctionBegin;
10219566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nr * nc, &b));
1022841e96a3SJed Brown   for (i = 0; i < nr; i++) {
1023841e96a3SJed Brown     for (j = 0; j < nc; j++) {
1024841e96a3SJed Brown       if (bA->m[i][j]) {
10259566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(bA->m[i][j], op, &b[i * nc + j]));
1026841e96a3SJed Brown       } else {
10270298fd71SBarry Smith         b[i * nc + j] = NULL;
1028d8588912SDave May       }
1029d8588912SDave May     }
1030d8588912SDave May   }
10319566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PetscObjectComm((PetscObject)A), nr, bA->isglobal.row, nc, bA->isglobal.col, b, B));
1032841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
103348a46eb9SPierre Jolivet   for (i = 0; i < nr * nc; i++) PetscCall(MatDestroy(&b[i]));
10349566063dSJacob Faibussowitsch   PetscCall(PetscFree(b));
1035d8588912SDave May 
10369566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*B, MAT_FINAL_ASSEMBLY));
10379566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*B, MAT_FINAL_ASSEMBLY));
10383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1039d8588912SDave May }
1040d8588912SDave May 
1041d8588912SDave May /* nest api */
104266976f2fSJacob Faibussowitsch static PetscErrorCode MatNestGetSubMat_Nest(Mat A, PetscInt idxm, PetscInt jdxm, Mat *mat)
1043d71ae5a4SJacob Faibussowitsch {
1044d8588912SDave May   Mat_Nest *bA = (Mat_Nest *)A->data;
10455fd66863SKarl Rupp 
1046d8588912SDave May   PetscFunctionBegin;
104708401ef6SPierre Jolivet   PetscCheck(idxm < bA->nr, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, idxm, bA->nr - 1);
104808401ef6SPierre Jolivet   PetscCheck(jdxm < bA->nc, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Col too large: row %" PetscInt_FMT " max %" PetscInt_FMT, jdxm, bA->nc - 1);
1049d8588912SDave May   *mat = bA->m[idxm][jdxm];
10503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1051d8588912SDave May }
1052d8588912SDave May 
10539ba0d327SJed Brown /*@
105411a5261eSBarry Smith   MatNestGetSubMat - Returns a single, sub-matrix from a `MATNEST`
1055d8588912SDave May 
10562ef1f0ffSBarry Smith   Not Collective
1057d8588912SDave May 
1058d8588912SDave May   Input Parameters:
105911a5261eSBarry Smith + A    - `MATNEST` matrix
1060d8588912SDave May . idxm - index of the matrix within the nest matrix
1061629881c0SJed Brown - jdxm - index of the matrix within the nest matrix
1062d8588912SDave May 
1063d8588912SDave May   Output Parameter:
10642ef1f0ffSBarry Smith . sub - matrix at index `idxm`, `jdxm` within the nest matrix
1065d8588912SDave May 
1066d8588912SDave May   Level: developer
1067d8588912SDave May 
1068fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatNestGetSize()`, `MatNestGetSubMats()`, `MatCreateNest()`, `MatNestSetSubMat()`,
1069db781477SPatrick Sanan           `MatNestGetLocalISs()`, `MatNestGetISs()`
1070d8588912SDave May @*/
1071d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNestGetSubMat(Mat A, PetscInt idxm, PetscInt jdxm, Mat *sub)
1072d71ae5a4SJacob Faibussowitsch {
1073d8588912SDave May   PetscFunctionBegin;
1074cac4c232SBarry Smith   PetscUseMethod(A, "MatNestGetSubMat_C", (Mat, PetscInt, PetscInt, Mat *), (A, idxm, jdxm, sub));
10753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1076d8588912SDave May }
1077d8588912SDave May 
107866976f2fSJacob Faibussowitsch static PetscErrorCode MatNestSetSubMat_Nest(Mat A, PetscInt idxm, PetscInt jdxm, Mat mat)
1079d71ae5a4SJacob Faibussowitsch {
10800782ca92SJed Brown   Mat_Nest *bA = (Mat_Nest *)A->data;
10810782ca92SJed Brown   PetscInt  m, n, M, N, mi, ni, Mi, Ni;
10820782ca92SJed Brown 
10830782ca92SJed Brown   PetscFunctionBegin;
108408401ef6SPierre Jolivet   PetscCheck(idxm < bA->nr, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Row too large: row %" PetscInt_FMT " max %" PetscInt_FMT, idxm, bA->nr - 1);
108508401ef6SPierre Jolivet   PetscCheck(jdxm < bA->nc, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Col too large: row %" PetscInt_FMT " max %" PetscInt_FMT, jdxm, bA->nc - 1);
10869566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(mat, &m, &n));
10879566063dSJacob Faibussowitsch   PetscCall(MatGetSize(mat, &M, &N));
10889566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(bA->isglobal.row[idxm], &mi));
10899566063dSJacob Faibussowitsch   PetscCall(ISGetSize(bA->isglobal.row[idxm], &Mi));
10909566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(bA->isglobal.col[jdxm], &ni));
10919566063dSJacob Faibussowitsch   PetscCall(ISGetSize(bA->isglobal.col[jdxm], &Ni));
1092aed4548fSBarry Smith   PetscCheck(M == Mi && N == Ni, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_INCOMP, "Submatrix dimension (%" PetscInt_FMT ",%" PetscInt_FMT ") incompatible with nest block (%" PetscInt_FMT ",%" PetscInt_FMT ")", M, N, Mi, Ni);
1093aed4548fSBarry Smith   PetscCheck(m == mi && n == ni, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_INCOMP, "Submatrix local dimension (%" PetscInt_FMT ",%" PetscInt_FMT ") incompatible with nest block (%" PetscInt_FMT ",%" PetscInt_FMT ")", m, n, mi, ni);
109426fbe8dcSKarl Rupp 
109506a1af2fSStefano Zampini   /* do not increase object state */
10963ba16761SJacob Faibussowitsch   if (mat == bA->m[idxm][jdxm]) PetscFunctionReturn(PETSC_SUCCESS);
109706a1af2fSStefano Zampini 
10989566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
10999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&bA->m[idxm][jdxm]));
11000782ca92SJed Brown   bA->m[idxm][jdxm] = mat;
11019566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)A));
11029566063dSJacob Faibussowitsch   PetscCall(MatGetNonzeroState(mat, &bA->nnzstate[idxm * bA->nc + jdxm]));
110306a1af2fSStefano Zampini   A->nonzerostate++;
11043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11050782ca92SJed Brown }
11060782ca92SJed Brown 
11079ba0d327SJed Brown /*@
110811a5261eSBarry Smith   MatNestSetSubMat - Set a single submatrix in the `MATNEST`
11090782ca92SJed Brown 
11102ef1f0ffSBarry Smith   Logically Collective
11110782ca92SJed Brown 
11120782ca92SJed Brown   Input Parameters:
111311a5261eSBarry Smith + A    - `MATNEST` matrix
11140782ca92SJed Brown . idxm - index of the matrix within the nest matrix
11150782ca92SJed Brown . jdxm - index of the matrix within the nest matrix
11162ef1f0ffSBarry Smith - sub  - matrix at index `idxm`, `jdxm` within the nest matrix
11172ef1f0ffSBarry Smith 
11182ef1f0ffSBarry Smith   Level: developer
11190782ca92SJed Brown 
11200782ca92SJed Brown   Notes:
11210782ca92SJed Brown   The new submatrix must have the same size and communicator as that block of the nest.
11220782ca92SJed Brown 
11230782ca92SJed Brown   This increments the reference count of the submatrix.
11240782ca92SJed Brown 
1125fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatNestSetSubMats()`, `MatNestGetSubMats()`, `MatNestGetLocalISs()`, `MatCreateNest()`,
1126db781477SPatrick Sanan           `MatNestGetSubMat()`, `MatNestGetISs()`, `MatNestGetSize()`
11270782ca92SJed Brown @*/
1128d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNestSetSubMat(Mat A, PetscInt idxm, PetscInt jdxm, Mat sub)
1129d71ae5a4SJacob Faibussowitsch {
11300782ca92SJed Brown   PetscFunctionBegin;
1131cac4c232SBarry Smith   PetscUseMethod(A, "MatNestSetSubMat_C", (Mat, PetscInt, PetscInt, Mat), (A, idxm, jdxm, sub));
11323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11330782ca92SJed Brown }
11340782ca92SJed Brown 
113566976f2fSJacob Faibussowitsch static PetscErrorCode MatNestGetSubMats_Nest(Mat A, PetscInt *M, PetscInt *N, Mat ***mat)
1136d71ae5a4SJacob Faibussowitsch {
1137d8588912SDave May   Mat_Nest *bA = (Mat_Nest *)A->data;
11385fd66863SKarl Rupp 
1139d8588912SDave May   PetscFunctionBegin;
114026fbe8dcSKarl Rupp   if (M) *M = bA->nr;
114126fbe8dcSKarl Rupp   if (N) *N = bA->nc;
114226fbe8dcSKarl Rupp   if (mat) *mat = bA->m;
11433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1144d8588912SDave May }
1145d8588912SDave May 
1146d8588912SDave May /*@C
114711a5261eSBarry Smith   MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a `MATNEST` matrix.
1148d8588912SDave May 
11492ef1f0ffSBarry Smith   Not Collective
1150d8588912SDave May 
1151f899ff85SJose E. Roman   Input Parameter:
1152629881c0SJed Brown . A - nest matrix
1153d8588912SDave May 
1154d8d19677SJose E. Roman   Output Parameters:
1155629881c0SJed Brown + M   - number of rows in the nest matrix
1156d8588912SDave May . N   - number of cols in the nest matrix
1157*e9d3347aSJose E. Roman - mat - array of matrices
1158d8588912SDave May 
11592ef1f0ffSBarry Smith   Level: developer
11602ef1f0ffSBarry Smith 
116111a5261eSBarry Smith   Note:
11622ef1f0ffSBarry Smith   The user should not free the array `mat`.
1163d8588912SDave May 
1164fe59aa6dSJacob Faibussowitsch   Fortran Notes:
116511a5261eSBarry Smith   This routine has a calling sequence
1166351962e3SVincent Le Chenadec $   call MatNestGetSubMats(A, M, N, mat, ierr)
116720f4b53cSBarry Smith   where the space allocated for the optional argument `mat` is assumed large enough (if provided).
1168*e9d3347aSJose E. Roman   Matrices in `mat` are returned in row-major order, see `MatCreateNest()` for an example.
1169351962e3SVincent Le Chenadec 
1170fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatNestGetSize()`, `MatNestGetSubMat()`, `MatNestGetLocalISs()`, `MatCreateNest()`,
1171db781477SPatrick Sanan           `MatNestSetSubMats()`, `MatNestGetISs()`, `MatNestSetSubMat()`
1172d8588912SDave May @*/
1173d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNestGetSubMats(Mat A, PetscInt *M, PetscInt *N, Mat ***mat)
1174d71ae5a4SJacob Faibussowitsch {
1175d8588912SDave May   PetscFunctionBegin;
1176cac4c232SBarry Smith   PetscUseMethod(A, "MatNestGetSubMats_C", (Mat, PetscInt *, PetscInt *, Mat ***), (A, M, N, mat));
11773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1178d8588912SDave May }
1179d8588912SDave May 
118066976f2fSJacob Faibussowitsch static PetscErrorCode MatNestGetSize_Nest(Mat A, PetscInt *M, PetscInt *N)
1181d71ae5a4SJacob Faibussowitsch {
1182d8588912SDave May   Mat_Nest *bA = (Mat_Nest *)A->data;
1183d8588912SDave May 
1184d8588912SDave May   PetscFunctionBegin;
118526fbe8dcSKarl Rupp   if (M) *M = bA->nr;
118626fbe8dcSKarl Rupp   if (N) *N = bA->nc;
11873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1188d8588912SDave May }
1189d8588912SDave May 
11909ba0d327SJed Brown /*@
119111a5261eSBarry Smith   MatNestGetSize - Returns the size of the `MATNEST` matrix.
1192d8588912SDave May 
11932ef1f0ffSBarry Smith   Not Collective
1194d8588912SDave May 
1195f899ff85SJose E. Roman   Input Parameter:
119611a5261eSBarry Smith . A - `MATNEST` matrix
1197d8588912SDave May 
1198d8d19677SJose E. Roman   Output Parameters:
1199629881c0SJed Brown + M - number of rows in the nested mat
1200629881c0SJed Brown - N - number of cols in the nested mat
1201d8588912SDave May 
1202d8588912SDave May   Level: developer
1203d8588912SDave May 
1204fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatNestGetSubMat()`, `MatNestGetSubMats()`, `MatCreateNest()`, `MatNestGetLocalISs()`,
1205db781477SPatrick Sanan           `MatNestGetISs()`
1206d8588912SDave May @*/
1207d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNestGetSize(Mat A, PetscInt *M, PetscInt *N)
1208d71ae5a4SJacob Faibussowitsch {
1209d8588912SDave May   PetscFunctionBegin;
1210cac4c232SBarry Smith   PetscUseMethod(A, "MatNestGetSize_C", (Mat, PetscInt *, PetscInt *), (A, M, N));
12113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1212d8588912SDave May }
1213d8588912SDave May 
1214d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNestGetISs_Nest(Mat A, IS rows[], IS cols[])
1215d71ae5a4SJacob Faibussowitsch {
1216900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest *)A->data;
1217900e7ff2SJed Brown   PetscInt  i;
1218900e7ff2SJed Brown 
1219900e7ff2SJed Brown   PetscFunctionBegin;
12209371c9d4SSatish Balay   if (rows)
12219371c9d4SSatish Balay     for (i = 0; i < vs->nr; i++) rows[i] = vs->isglobal.row[i];
12229371c9d4SSatish Balay   if (cols)
12239371c9d4SSatish Balay     for (i = 0; i < vs->nc; i++) cols[i] = vs->isglobal.col[i];
12243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1225900e7ff2SJed Brown }
1226900e7ff2SJed Brown 
12273a4d7b9aSSatish Balay /*@C
122811a5261eSBarry Smith   MatNestGetISs - Returns the index sets partitioning the row and column spaces of a `MATNEST`
1229900e7ff2SJed Brown 
12302ef1f0ffSBarry Smith   Not Collective
1231900e7ff2SJed Brown 
1232f899ff85SJose E. Roman   Input Parameter:
123311a5261eSBarry Smith . A - `MATNEST` matrix
1234900e7ff2SJed Brown 
1235d8d19677SJose E. Roman   Output Parameters:
1236900e7ff2SJed Brown + rows - array of row index sets
1237900e7ff2SJed Brown - cols - array of column index sets
1238900e7ff2SJed Brown 
1239900e7ff2SJed Brown   Level: advanced
1240900e7ff2SJed Brown 
124111a5261eSBarry Smith   Note:
12422ef1f0ffSBarry Smith   The user must have allocated arrays of the correct size. The reference count is not increased on the returned `IS`s.
1243900e7ff2SJed Brown 
1244fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatNestGetSubMat()`, `MatNestGetSubMats()`, `MatNestGetSize()`, `MatNestGetLocalISs()`,
1245fe59aa6dSJacob Faibussowitsch           `MatCreateNest()`, `MatNestSetSubMats()`
1246900e7ff2SJed Brown @*/
1247d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNestGetISs(Mat A, IS rows[], IS cols[])
1248d71ae5a4SJacob Faibussowitsch {
1249900e7ff2SJed Brown   PetscFunctionBegin;
1250900e7ff2SJed Brown   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1251cac4c232SBarry Smith   PetscUseMethod(A, "MatNestGetISs_C", (Mat, IS[], IS[]), (A, rows, cols));
12523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1253900e7ff2SJed Brown }
1254900e7ff2SJed Brown 
1255d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNestGetLocalISs_Nest(Mat A, IS rows[], IS cols[])
1256d71ae5a4SJacob Faibussowitsch {
1257900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest *)A->data;
1258900e7ff2SJed Brown   PetscInt  i;
1259900e7ff2SJed Brown 
1260900e7ff2SJed Brown   PetscFunctionBegin;
12619371c9d4SSatish Balay   if (rows)
12629371c9d4SSatish Balay     for (i = 0; i < vs->nr; i++) rows[i] = vs->islocal.row[i];
12639371c9d4SSatish Balay   if (cols)
12649371c9d4SSatish Balay     for (i = 0; i < vs->nc; i++) cols[i] = vs->islocal.col[i];
12653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1266900e7ff2SJed Brown }
1267900e7ff2SJed Brown 
1268900e7ff2SJed Brown /*@C
126911a5261eSBarry Smith   MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces of a `MATNEST`
1270900e7ff2SJed Brown 
12712ef1f0ffSBarry Smith   Not Collective
1272900e7ff2SJed Brown 
1273f899ff85SJose E. Roman   Input Parameter:
127411a5261eSBarry Smith . A - `MATNEST` matrix
1275900e7ff2SJed Brown 
1276d8d19677SJose E. Roman   Output Parameters:
12772ef1f0ffSBarry Smith + rows - array of row index sets (or `NULL` to ignore)
12782ef1f0ffSBarry Smith - cols - array of column index sets (or `NULL` to ignore)
1279900e7ff2SJed Brown 
1280900e7ff2SJed Brown   Level: advanced
1281900e7ff2SJed Brown 
128211a5261eSBarry Smith   Note:
12832ef1f0ffSBarry Smith   The user must have allocated arrays of the correct size. The reference count is not increased on the returned `IS`s.
1284900e7ff2SJed Brown 
12851cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatNestGetSubMat()`, `MatNestGetSubMats()`, `MatNestGetSize()`, `MatNestGetISs()`, `MatCreateNest()`,
1286fe59aa6dSJacob Faibussowitsch           `MatNestSetSubMats()`, `MatNestSetSubMat()`
1287900e7ff2SJed Brown @*/
1288d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNestGetLocalISs(Mat A, IS rows[], IS cols[])
1289d71ae5a4SJacob Faibussowitsch {
1290900e7ff2SJed Brown   PetscFunctionBegin;
1291900e7ff2SJed Brown   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1292cac4c232SBarry Smith   PetscUseMethod(A, "MatNestGetLocalISs_C", (Mat, IS[], IS[]), (A, rows, cols));
12933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1294900e7ff2SJed Brown }
1295900e7ff2SJed Brown 
129666976f2fSJacob Faibussowitsch static PetscErrorCode MatNestSetVecType_Nest(Mat A, VecType vtype)
1297d71ae5a4SJacob Faibussowitsch {
1298207556f9SJed Brown   PetscBool flg;
1299207556f9SJed Brown 
1300207556f9SJed Brown   PetscFunctionBegin;
13019566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(vtype, VECNEST, &flg));
1302207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
13032a7a6963SBarry Smith   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
130412b53f24SSatish Balay   else A->ops->getvecs = (PetscErrorCode(*)(Mat, Vec *, Vec *))0;
13053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1306207556f9SJed Brown }
1307207556f9SJed Brown 
1308207556f9SJed Brown /*@C
130911a5261eSBarry Smith   MatNestSetVecType - Sets the type of `Vec` returned by `MatCreateVecs()`
1310207556f9SJed Brown 
13112ef1f0ffSBarry Smith   Not Collective
1312207556f9SJed Brown 
1313207556f9SJed Brown   Input Parameters:
131411a5261eSBarry Smith + A     - `MATNEST` matrix
131511a5261eSBarry Smith - vtype - `VecType` to use for creating vectors
1316207556f9SJed Brown 
1317207556f9SJed Brown   Level: developer
1318207556f9SJed Brown 
1319fe59aa6dSJacob Faibussowitsch .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatCreateVecs()`, `MatCreateNest()`, `VecType`
1320207556f9SJed Brown @*/
1321d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNestSetVecType(Mat A, VecType vtype)
1322d71ae5a4SJacob Faibussowitsch {
1323207556f9SJed Brown   PetscFunctionBegin;
1324cac4c232SBarry Smith   PetscTryMethod(A, "MatNestSetVecType_C", (Mat, VecType), (A, vtype));
13253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1326207556f9SJed Brown }
1327207556f9SJed Brown 
132866976f2fSJacob Faibussowitsch static PetscErrorCode MatNestSetSubMats_Nest(Mat A, PetscInt nr, const IS is_row[], PetscInt nc, const IS is_col[], const Mat a[])
1329d71ae5a4SJacob Faibussowitsch {
1330c8883902SJed Brown   Mat_Nest *s = (Mat_Nest *)A->data;
1331c8883902SJed Brown   PetscInt  i, j, m, n, M, N;
133288ffe2e8SJose E. Roman   PetscBool cong, isstd, sametype = PETSC_FALSE;
133388ffe2e8SJose E. Roman   VecType   vtype, type;
1334d8588912SDave May 
1335d8588912SDave May   PetscFunctionBegin;
13369566063dSJacob Faibussowitsch   PetscCall(MatReset_Nest(A));
133706a1af2fSStefano Zampini 
1338c8883902SJed Brown   s->nr = nr;
1339c8883902SJed Brown   s->nc = nc;
1340d8588912SDave May 
1341c8883902SJed Brown   /* Create space for submatrices */
13429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nr, &s->m));
134348a46eb9SPierre Jolivet   for (i = 0; i < nr; i++) PetscCall(PetscMalloc1(nc, &s->m[i]));
1344c8883902SJed Brown   for (i = 0; i < nr; i++) {
1345c8883902SJed Brown     for (j = 0; j < nc; j++) {
1346c8883902SJed Brown       s->m[i][j] = a[i * nc + j];
134748a46eb9SPierre Jolivet       if (a[i * nc + j]) PetscCall(PetscObjectReference((PetscObject)a[i * nc + j]));
1348d8588912SDave May     }
1349d8588912SDave May   }
13509566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(A, &vtype));
13519566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(vtype, VECSTANDARD, &isstd));
135288ffe2e8SJose E. Roman   if (isstd) {
135388ffe2e8SJose E. Roman     /* check if all blocks have the same vectype */
135488ffe2e8SJose E. Roman     vtype = NULL;
135588ffe2e8SJose E. Roman     for (i = 0; i < nr; i++) {
135688ffe2e8SJose E. Roman       for (j = 0; j < nc; j++) {
135788ffe2e8SJose E. Roman         if (a[i * nc + j]) {
135888ffe2e8SJose E. Roman           if (!vtype) { /* first visited block */
13599566063dSJacob Faibussowitsch             PetscCall(MatGetVecType(a[i * nc + j], &vtype));
136088ffe2e8SJose E. Roman             sametype = PETSC_TRUE;
136188ffe2e8SJose E. Roman           } else if (sametype) {
13629566063dSJacob Faibussowitsch             PetscCall(MatGetVecType(a[i * nc + j], &type));
13639566063dSJacob Faibussowitsch             PetscCall(PetscStrcmp(vtype, type, &sametype));
136488ffe2e8SJose E. Roman           }
136588ffe2e8SJose E. Roman         }
136688ffe2e8SJose E. Roman       }
136788ffe2e8SJose E. Roman     }
136888ffe2e8SJose E. Roman     if (sametype) { /* propagate vectype */
13699566063dSJacob Faibussowitsch       PetscCall(MatSetVecType(A, vtype));
137088ffe2e8SJose E. Roman     }
137188ffe2e8SJose E. Roman   }
1372d8588912SDave May 
13739566063dSJacob Faibussowitsch   PetscCall(MatSetUp_NestIS_Private(A, nr, is_row, nc, is_col));
1374d8588912SDave May 
13759566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nr, &s->row_len));
13769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nc, &s->col_len));
1377c8883902SJed Brown   for (i = 0; i < nr; i++) s->row_len[i] = -1;
1378c8883902SJed Brown   for (j = 0; j < nc; j++) s->col_len[j] = -1;
1379d8588912SDave May 
13809566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nr * nc, &s->nnzstate));
138106a1af2fSStefano Zampini   for (i = 0; i < nr; i++) {
138206a1af2fSStefano Zampini     for (j = 0; j < nc; j++) {
138348a46eb9SPierre Jolivet       if (s->m[i][j]) PetscCall(MatGetNonzeroState(s->m[i][j], &s->nnzstate[i * nc + j]));
138406a1af2fSStefano Zampini     }
138506a1af2fSStefano Zampini   }
138606a1af2fSStefano Zampini 
13879566063dSJacob Faibussowitsch   PetscCall(MatNestGetSizes_Private(A, &m, &n, &M, &N));
1388d8588912SDave May 
13899566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetSize(A->rmap, M));
13909566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(A->rmap, m));
13919566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetSize(A->cmap, N));
13929566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(A->cmap, n));
1393c8883902SJed Brown 
13949566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->rmap));
13959566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(A->cmap));
1396c8883902SJed Brown 
139706a1af2fSStefano Zampini   /* disable operations that are not supported for non-square matrices,
139806a1af2fSStefano Zampini      or matrices for which is_row != is_col  */
13999566063dSJacob Faibussowitsch   PetscCall(MatHasCongruentLayouts(A, &cong));
140006a1af2fSStefano Zampini   if (cong && nr != nc) cong = PETSC_FALSE;
140106a1af2fSStefano Zampini   if (cong) {
140248a46eb9SPierre Jolivet     for (i = 0; cong && i < nr; i++) PetscCall(ISEqualUnsorted(s->isglobal.row[i], s->isglobal.col[i], &cong));
140306a1af2fSStefano Zampini   }
140406a1af2fSStefano Zampini   if (!cong) {
1405381b8e50SStefano Zampini     A->ops->missingdiagonal = NULL;
140606a1af2fSStefano Zampini     A->ops->getdiagonal     = NULL;
140706a1af2fSStefano Zampini     A->ops->shift           = NULL;
140806a1af2fSStefano Zampini     A->ops->diagonalset     = NULL;
140906a1af2fSStefano Zampini   }
141006a1af2fSStefano Zampini 
14119566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(nr, &s->left, nc, &s->right));
14129566063dSJacob Faibussowitsch   PetscCall(PetscObjectStateIncrease((PetscObject)A));
141306a1af2fSStefano Zampini   A->nonzerostate++;
14143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1415d8588912SDave May }
1416d8588912SDave May 
1417c8883902SJed Brown /*@
141811a5261eSBarry Smith   MatNestSetSubMats - Sets the nested submatrices in a `MATNEST`
1419c8883902SJed Brown 
1420c3339decSBarry Smith   Collective
1421c8883902SJed Brown 
1422d8d19677SJose E. Roman   Input Parameters:
142311a5261eSBarry Smith + A      - `MATNEST` matrix
1424c8883902SJed Brown . nr     - number of nested row blocks
14252ef1f0ffSBarry Smith . is_row - index sets for each nested row block, or `NULL` to make contiguous
1426c8883902SJed Brown . nc     - number of nested column blocks
14272ef1f0ffSBarry Smith . is_col - index sets for each nested column block, or `NULL` to make contiguous
1428*e9d3347aSJose E. Roman - a      - array of nr*nc submatrices, empty submatrices can be passed using `NULL`
14292ef1f0ffSBarry Smith 
14302ef1f0ffSBarry Smith   Level: advanced
1431c8883902SJed Brown 
1432*e9d3347aSJose E. Roman   Notes:
143311a5261eSBarry Smith   This always resets any submatrix information previously set
143406a1af2fSStefano Zampini 
1435*e9d3347aSJose E. Roman   In both C and Fortran, `a` must be a row-major order array containing the matrices. See
1436*e9d3347aSJose E. Roman   `MatCreateNest()` for an example.
1437*e9d3347aSJose E. Roman 
14381cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatCreateNest()`, `MatNestSetSubMat()`, `MatNestGetSubMat()`, `MatNestGetSubMats()`
1439c8883902SJed Brown @*/
1440d71ae5a4SJacob Faibussowitsch PetscErrorCode MatNestSetSubMats(Mat A, PetscInt nr, const IS is_row[], PetscInt nc, const IS is_col[], const Mat a[])
1441d71ae5a4SJacob Faibussowitsch {
144206a1af2fSStefano Zampini   PetscInt i;
1443c8883902SJed Brown 
1444c8883902SJed Brown   PetscFunctionBegin;
1445c8883902SJed Brown   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
144608401ef6SPierre Jolivet   PetscCheck(nr >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Number of rows cannot be negative");
1447c8883902SJed Brown   if (nr && is_row) {
14484f572ea9SToby Isaac     PetscAssertPointer(is_row, 3);
1449c8883902SJed Brown     for (i = 0; i < nr; i++) PetscValidHeaderSpecific(is_row[i], IS_CLASSID, 3);
1450c8883902SJed Brown   }
145108401ef6SPierre Jolivet   PetscCheck(nc >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_OUTOFRANGE, "Number of columns cannot be negative");
14521664e352SJed Brown   if (nc && is_col) {
14534f572ea9SToby Isaac     PetscAssertPointer(is_col, 5);
14549b30a8f6SBarry Smith     for (i = 0; i < nc; i++) PetscValidHeaderSpecific(is_col[i], IS_CLASSID, 5);
1455c8883902SJed Brown   }
14564f572ea9SToby Isaac   if (nr * nc > 0) PetscAssertPointer(a, 6);
1457cac4c232SBarry Smith   PetscUseMethod(A, "MatNestSetSubMats_C", (Mat, PetscInt, const IS[], PetscInt, const IS[], const Mat[]), (A, nr, is_row, nc, is_col, a));
14583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1459c8883902SJed Brown }
1460d8588912SDave May 
1461d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A, PetscInt n, const IS islocal[], const IS isglobal[], PetscBool colflg, ISLocalToGlobalMapping *ltog)
1462d71ae5a4SJacob Faibussowitsch {
146377019fcaSJed Brown   PetscBool flg;
146477019fcaSJed Brown   PetscInt  i, j, m, mi, *ix;
146577019fcaSJed Brown 
146677019fcaSJed Brown   PetscFunctionBegin;
1467aea6d515SStefano Zampini   *ltog = NULL;
146877019fcaSJed Brown   for (i = 0, m = 0, flg = PETSC_FALSE; i < n; i++) {
146977019fcaSJed Brown     if (islocal[i]) {
14709566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(islocal[i], &mi));
147177019fcaSJed Brown       flg = PETSC_TRUE; /* We found a non-trivial entry */
147277019fcaSJed Brown     } else {
14739566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(isglobal[i], &mi));
147477019fcaSJed Brown     }
147577019fcaSJed Brown     m += mi;
147677019fcaSJed Brown   }
14773ba16761SJacob Faibussowitsch   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
1478aea6d515SStefano Zampini 
14799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &ix));
1480165cd838SBarry Smith   for (i = 0, m = 0; i < n; i++) {
14810298fd71SBarry Smith     ISLocalToGlobalMapping smap = NULL;
1482e108cb99SStefano Zampini     Mat                    sub  = NULL;
1483f6d38dbbSStefano Zampini     PetscSF                sf;
1484f6d38dbbSStefano Zampini     PetscLayout            map;
1485aea6d515SStefano Zampini     const PetscInt        *ix2;
148677019fcaSJed Brown 
1487165cd838SBarry Smith     if (!colflg) {
14889566063dSJacob Faibussowitsch       PetscCall(MatNestFindNonzeroSubMatRow(A, i, &sub));
148977019fcaSJed Brown     } else {
14909566063dSJacob Faibussowitsch       PetscCall(MatNestFindNonzeroSubMatCol(A, i, &sub));
149177019fcaSJed Brown     }
1492191fd14bSBarry Smith     if (sub) {
1493191fd14bSBarry Smith       if (!colflg) {
14949566063dSJacob Faibussowitsch         PetscCall(MatGetLocalToGlobalMapping(sub, &smap, NULL));
1495191fd14bSBarry Smith       } else {
14969566063dSJacob Faibussowitsch         PetscCall(MatGetLocalToGlobalMapping(sub, NULL, &smap));
1497191fd14bSBarry Smith       }
1498191fd14bSBarry Smith     }
149977019fcaSJed Brown     /*
150077019fcaSJed Brown        Now we need to extract the monolithic global indices that correspond to the given split global indices.
150177019fcaSJed Brown        In many/most cases, we only want MatGetLocalSubMatrix() to work, in which case we only need to know the size of the local spaces.
150277019fcaSJed Brown     */
15039566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(isglobal[i], &ix2));
1504aea6d515SStefano Zampini     if (islocal[i]) {
1505aea6d515SStefano Zampini       PetscInt *ilocal, *iremote;
1506aea6d515SStefano Zampini       PetscInt  mil, nleaves;
1507aea6d515SStefano Zampini 
15089566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(islocal[i], &mi));
150928b400f6SJacob Faibussowitsch       PetscCheck(smap, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing local to global map");
1510aea6d515SStefano Zampini       for (j = 0; j < mi; j++) ix[m + j] = j;
15119566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingApply(smap, mi, ix + m, ix + m));
1512aea6d515SStefano Zampini 
1513aea6d515SStefano Zampini       /* PetscSFSetGraphLayout does not like negative indices */
15149566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(mi, &ilocal, mi, &iremote));
1515aea6d515SStefano Zampini       for (j = 0, nleaves = 0; j < mi; j++) {
1516aea6d515SStefano Zampini         if (ix[m + j] < 0) continue;
1517aea6d515SStefano Zampini         ilocal[nleaves]  = j;
1518aea6d515SStefano Zampini         iremote[nleaves] = ix[m + j];
1519aea6d515SStefano Zampini         nleaves++;
1520aea6d515SStefano Zampini       }
15219566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(isglobal[i], &mil));
15229566063dSJacob Faibussowitsch       PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &sf));
15239566063dSJacob Faibussowitsch       PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)A), &map));
15249566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetLocalSize(map, mil));
15259566063dSJacob Faibussowitsch       PetscCall(PetscLayoutSetUp(map));
15269566063dSJacob Faibussowitsch       PetscCall(PetscSFSetGraphLayout(sf, map, nleaves, ilocal, PETSC_USE_POINTER, iremote));
15279566063dSJacob Faibussowitsch       PetscCall(PetscLayoutDestroy(&map));
15289566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastBegin(sf, MPIU_INT, ix2, ix + m, MPI_REPLACE));
15299566063dSJacob Faibussowitsch       PetscCall(PetscSFBcastEnd(sf, MPIU_INT, ix2, ix + m, MPI_REPLACE));
15309566063dSJacob Faibussowitsch       PetscCall(PetscSFDestroy(&sf));
15319566063dSJacob Faibussowitsch       PetscCall(PetscFree2(ilocal, iremote));
1532aea6d515SStefano Zampini     } else {
15339566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(isglobal[i], &mi));
1534aea6d515SStefano Zampini       for (j = 0; j < mi; j++) ix[m + j] = ix2[i];
1535aea6d515SStefano Zampini     }
15369566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(isglobal[i], &ix2));
153777019fcaSJed Brown     m += mi;
153877019fcaSJed Brown   }
15399566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A), 1, m, ix, PETSC_OWN_POINTER, ltog));
15403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
154177019fcaSJed Brown }
154277019fcaSJed Brown 
1543d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1544d8588912SDave May /*
1545d8588912SDave May   nprocessors = NP
1546d8588912SDave May   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1547d8588912SDave May        proc 0: => (g_0,h_0,)
1548d8588912SDave May        proc 1: => (g_1,h_1,)
1549d8588912SDave May        ...
1550d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1551d8588912SDave May 
1552d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1553d8588912SDave May     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))
1554d8588912SDave May 
1555d8588912SDave May             proc 0:
1556d8588912SDave May     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1557d8588912SDave May             proc 1:
1558d8588912SDave May     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)
1559d8588912SDave May 
1560d8588912SDave May             proc NP-1:
1561d8588912SDave May     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1562d8588912SDave May */
1563d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetUp_NestIS_Private(Mat A, PetscInt nr, const IS is_row[], PetscInt nc, const IS is_col[])
1564d71ae5a4SJacob Faibussowitsch {
1565e2d7f03fSJed Brown   Mat_Nest *vs = (Mat_Nest *)A->data;
15668188e55aSJed Brown   PetscInt  i, j, offset, n, nsum, bs;
15670298fd71SBarry Smith   Mat       sub = NULL;
1568d8588912SDave May 
1569d8588912SDave May   PetscFunctionBegin;
15709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nr, &vs->isglobal.row));
15719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nc, &vs->isglobal.col));
1572d8588912SDave May   if (is_row) { /* valid IS is passed in */
1573a5b23f4aSJose E. Roman     /* refs on is[] are incremented */
1574e2d7f03fSJed Brown     for (i = 0; i < vs->nr; i++) {
15759566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)is_row[i]));
157626fbe8dcSKarl Rupp 
1577e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1578d8588912SDave May     }
15792ae74bdbSJed Brown   } else { /* Create the ISs by inspecting sizes of a submatrix in each row */
15808188e55aSJed Brown     nsum = 0;
15818188e55aSJed Brown     for (i = 0; i < vs->nr; i++) { /* Add up the local sizes to compute the aggregate offset */
15829566063dSJacob Faibussowitsch       PetscCall(MatNestFindNonzeroSubMatRow(A, i, &sub));
158328b400f6SJacob Faibussowitsch       PetscCheck(sub, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "No nonzero submatrix in row %" PetscInt_FMT, i);
15849566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(sub, &n, NULL));
158508401ef6SPierre Jolivet       PetscCheck(n >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Sizes have not yet been set for submatrix");
15868188e55aSJed Brown       nsum += n;
15878188e55aSJed Brown     }
15889566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&nsum, &offset, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
158930bc264bSJed Brown     offset -= nsum;
1590e2d7f03fSJed Brown     for (i = 0; i < vs->nr; i++) {
15919566063dSJacob Faibussowitsch       PetscCall(MatNestFindNonzeroSubMatRow(A, i, &sub));
15929566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(sub, &n, NULL));
15939566063dSJacob Faibussowitsch       PetscCall(MatGetBlockSizes(sub, &bs, NULL));
15949566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PetscObjectComm((PetscObject)sub), n, offset, 1, &vs->isglobal.row[i]));
15959566063dSJacob Faibussowitsch       PetscCall(ISSetBlockSize(vs->isglobal.row[i], bs));
15962ae74bdbSJed Brown       offset += n;
1597d8588912SDave May     }
1598d8588912SDave May   }
1599d8588912SDave May 
1600d8588912SDave May   if (is_col) { /* valid IS is passed in */
1601a5b23f4aSJose E. Roman     /* refs on is[] are incremented */
1602e2d7f03fSJed Brown     for (j = 0; j < vs->nc; j++) {
16039566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)is_col[j]));
160426fbe8dcSKarl Rupp 
1605e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1606d8588912SDave May     }
16072ae74bdbSJed Brown   } else { /* Create the ISs by inspecting sizes of a submatrix in each column */
16082ae74bdbSJed Brown     offset = A->cmap->rstart;
16098188e55aSJed Brown     nsum   = 0;
16108188e55aSJed Brown     for (j = 0; j < vs->nc; j++) {
16119566063dSJacob Faibussowitsch       PetscCall(MatNestFindNonzeroSubMatCol(A, j, &sub));
161228b400f6SJacob Faibussowitsch       PetscCheck(sub, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "No nonzero submatrix in column %" PetscInt_FMT, i);
16139566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(sub, NULL, &n));
161408401ef6SPierre Jolivet       PetscCheck(n >= 0, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Sizes have not yet been set for submatrix");
16158188e55aSJed Brown       nsum += n;
16168188e55aSJed Brown     }
16179566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&nsum, &offset, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)A)));
161830bc264bSJed Brown     offset -= nsum;
1619e2d7f03fSJed Brown     for (j = 0; j < vs->nc; j++) {
16209566063dSJacob Faibussowitsch       PetscCall(MatNestFindNonzeroSubMatCol(A, j, &sub));
16219566063dSJacob Faibussowitsch       PetscCall(MatGetLocalSize(sub, NULL, &n));
16229566063dSJacob Faibussowitsch       PetscCall(MatGetBlockSizes(sub, NULL, &bs));
16239566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PetscObjectComm((PetscObject)sub), n, offset, 1, &vs->isglobal.col[j]));
16249566063dSJacob Faibussowitsch       PetscCall(ISSetBlockSize(vs->isglobal.col[j], bs));
16252ae74bdbSJed Brown       offset += n;
1626d8588912SDave May     }
1627d8588912SDave May   }
1628e2d7f03fSJed Brown 
1629e2d7f03fSJed Brown   /* Set up the local ISs */
16309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(vs->nr, &vs->islocal.row));
16319566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(vs->nc, &vs->islocal.col));
1632e2d7f03fSJed Brown   for (i = 0, offset = 0; i < vs->nr; i++) {
1633e2d7f03fSJed Brown     IS                     isloc;
16340298fd71SBarry Smith     ISLocalToGlobalMapping rmap = NULL;
1635e2d7f03fSJed Brown     PetscInt               nlocal, bs;
16369566063dSJacob Faibussowitsch     PetscCall(MatNestFindNonzeroSubMatRow(A, i, &sub));
16379566063dSJacob Faibussowitsch     if (sub) PetscCall(MatGetLocalToGlobalMapping(sub, &rmap, NULL));
1638207556f9SJed Brown     if (rmap) {
16399566063dSJacob Faibussowitsch       PetscCall(MatGetBlockSizes(sub, &bs, NULL));
16409566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingGetSize(rmap, &nlocal));
16419566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, nlocal, offset, 1, &isloc));
16429566063dSJacob Faibussowitsch       PetscCall(ISSetBlockSize(isloc, bs));
1643207556f9SJed Brown     } else {
1644207556f9SJed Brown       nlocal = 0;
16450298fd71SBarry Smith       isloc  = NULL;
1646207556f9SJed Brown     }
1647e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1648e2d7f03fSJed Brown     offset += nlocal;
1649e2d7f03fSJed Brown   }
16508188e55aSJed Brown   for (i = 0, offset = 0; i < vs->nc; i++) {
1651e2d7f03fSJed Brown     IS                     isloc;
16520298fd71SBarry Smith     ISLocalToGlobalMapping cmap = NULL;
1653e2d7f03fSJed Brown     PetscInt               nlocal, bs;
16549566063dSJacob Faibussowitsch     PetscCall(MatNestFindNonzeroSubMatCol(A, i, &sub));
16559566063dSJacob Faibussowitsch     if (sub) PetscCall(MatGetLocalToGlobalMapping(sub, NULL, &cmap));
1656207556f9SJed Brown     if (cmap) {
16579566063dSJacob Faibussowitsch       PetscCall(MatGetBlockSizes(sub, NULL, &bs));
16589566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingGetSize(cmap, &nlocal));
16599566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, nlocal, offset, 1, &isloc));
16609566063dSJacob Faibussowitsch       PetscCall(ISSetBlockSize(isloc, bs));
1661207556f9SJed Brown     } else {
1662207556f9SJed Brown       nlocal = 0;
16630298fd71SBarry Smith       isloc  = NULL;
1664207556f9SJed Brown     }
1665e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1666e2d7f03fSJed Brown     offset += nlocal;
1667e2d7f03fSJed Brown   }
16680189643fSJed Brown 
166977019fcaSJed Brown   /* Set up the aggregate ISLocalToGlobalMapping */
167077019fcaSJed Brown   {
167145b6f7e9SBarry Smith     ISLocalToGlobalMapping rmap, cmap;
16729566063dSJacob Faibussowitsch     PetscCall(MatNestCreateAggregateL2G_Private(A, vs->nr, vs->islocal.row, vs->isglobal.row, PETSC_FALSE, &rmap));
16739566063dSJacob Faibussowitsch     PetscCall(MatNestCreateAggregateL2G_Private(A, vs->nc, vs->islocal.col, vs->isglobal.col, PETSC_TRUE, &cmap));
16749566063dSJacob Faibussowitsch     if (rmap && cmap) PetscCall(MatSetLocalToGlobalMapping(A, rmap, cmap));
16759566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&rmap));
16769566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&cmap));
167777019fcaSJed Brown   }
167877019fcaSJed Brown 
167976bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
16800189643fSJed Brown     for (i = 0; i < vs->nr; i++) {
16810189643fSJed Brown       for (j = 0; j < vs->nc; j++) {
16820189643fSJed Brown         PetscInt m, n, M, N, mi, ni, Mi, Ni;
16830189643fSJed Brown         Mat      B = vs->m[i][j];
16840189643fSJed Brown         if (!B) continue;
16859566063dSJacob Faibussowitsch         PetscCall(MatGetSize(B, &M, &N));
16869566063dSJacob Faibussowitsch         PetscCall(MatGetLocalSize(B, &m, &n));
16879566063dSJacob Faibussowitsch         PetscCall(ISGetSize(vs->isglobal.row[i], &Mi));
16889566063dSJacob Faibussowitsch         PetscCall(ISGetSize(vs->isglobal.col[j], &Ni));
16899566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(vs->isglobal.row[i], &mi));
16909566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(vs->isglobal.col[j], &ni));
1691aed4548fSBarry Smith         PetscCheck(M == Mi && N == Ni, PetscObjectComm((PetscObject)sub), PETSC_ERR_ARG_INCOMP, "Global sizes (%" PetscInt_FMT ",%" PetscInt_FMT ") of nested submatrix (%" PetscInt_FMT ",%" PetscInt_FMT ") do not agree with space defined by index sets (%" PetscInt_FMT ",%" PetscInt_FMT ")", M, N, i, j, Mi, Ni);
1692aed4548fSBarry Smith         PetscCheck(m == mi && n == ni, PetscObjectComm((PetscObject)sub), PETSC_ERR_ARG_INCOMP, "Local sizes (%" PetscInt_FMT ",%" PetscInt_FMT ") of nested submatrix (%" PetscInt_FMT ",%" PetscInt_FMT ") do not agree with space defined by index sets (%" PetscInt_FMT ",%" PetscInt_FMT ")", m, n, i, j, mi, ni);
16930189643fSJed Brown       }
16940189643fSJed Brown     }
169576bd3646SJed Brown   }
1696a061e289SJed Brown 
1697a061e289SJed Brown   /* Set A->assembled if all non-null blocks are currently assembled */
1698a061e289SJed Brown   for (i = 0; i < vs->nr; i++) {
1699a061e289SJed Brown     for (j = 0; j < vs->nc; j++) {
17003ba16761SJacob Faibussowitsch       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(PETSC_SUCCESS);
1701a061e289SJed Brown     }
1702a061e289SJed Brown   }
1703a061e289SJed Brown   A->assembled = PETSC_TRUE;
17043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1705d8588912SDave May }
1706d8588912SDave May 
170745c38901SJed Brown /*@C
170811a5261eSBarry Smith   MatCreateNest - Creates a new `MATNEST` matrix containing several nested submatrices, each stored separately
1709659c6bb0SJed Brown 
171011a5261eSBarry Smith   Collective
1711659c6bb0SJed Brown 
1712d8d19677SJose E. Roman   Input Parameters:
171311a5261eSBarry Smith + comm   - Communicator for the new `MATNEST`
1714659c6bb0SJed Brown . nr     - number of nested row blocks
17152ef1f0ffSBarry Smith . is_row - index sets for each nested row block, or `NULL` to make contiguous
1716659c6bb0SJed Brown . nc     - number of nested column blocks
17172ef1f0ffSBarry Smith . is_col - index sets for each nested column block, or `NULL` to make contiguous
1718*e9d3347aSJose E. Roman - a      - array of nr*nc submatrices, empty submatrices can be passed using `NULL`
1719659c6bb0SJed Brown 
1720659c6bb0SJed Brown   Output Parameter:
1721659c6bb0SJed Brown . B - new matrix
1722659c6bb0SJed Brown 
1723*e9d3347aSJose E. Roman   Note:
1724*e9d3347aSJose E. Roman   In both C and Fortran, `a` must be a row-major order array holding references to the matrices.
1725*e9d3347aSJose E. Roman   For instance, to represent the matrix
1726*e9d3347aSJose E. Roman   $\begin{bmatrix} A_{11} & A_{12} \\ A_{21} & A_{22}\end{bmatrix}$
1727*e9d3347aSJose E. Roman   one should use `Mat a[4]={A11,A12,A21,A22}`.
1728*e9d3347aSJose E. Roman 
1729659c6bb0SJed Brown   Level: advanced
1730659c6bb0SJed Brown 
17311cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatCreate()`, `VecCreateNest()`, `DMCreateMatrix()`, `MatNestSetSubMat()`,
1732db781477SPatrick Sanan           `MatNestGetSubMat()`, `MatNestGetLocalISs()`, `MatNestGetSize()`,
1733db781477SPatrick Sanan           `MatNestGetISs()`, `MatNestSetSubMats()`, `MatNestGetSubMats()`
1734659c6bb0SJed Brown @*/
1735d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateNest(MPI_Comm comm, PetscInt nr, const IS is_row[], PetscInt nc, const IS is_col[], const Mat a[], Mat *B)
1736d71ae5a4SJacob Faibussowitsch {
1737d8588912SDave May   Mat A;
1738d8588912SDave May 
1739d8588912SDave May   PetscFunctionBegin;
1740f4259b30SLisandro Dalcin   *B = NULL;
17419566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &A));
17429566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATNEST));
174391a28eb3SBarry Smith   A->preallocated = PETSC_TRUE;
17449566063dSJacob Faibussowitsch   PetscCall(MatNestSetSubMats(A, nr, is_row, nc, is_col, a));
1745d8588912SDave May   *B = A;
17463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1747d8588912SDave May }
1748659c6bb0SJed Brown 
174966976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_Nest_SeqAIJ_fast(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
1750d71ae5a4SJacob Faibussowitsch {
1751b68353e5Sstefano_zampini   Mat_Nest     *nest = (Mat_Nest *)A->data;
175223875855Sstefano_zampini   Mat          *trans;
1753b68353e5Sstefano_zampini   PetscScalar **avv;
1754b68353e5Sstefano_zampini   PetscScalar  *vv;
1755b68353e5Sstefano_zampini   PetscInt    **aii, **ajj;
1756b68353e5Sstefano_zampini   PetscInt     *ii, *jj, *ci;
1757b68353e5Sstefano_zampini   PetscInt      nr, nc, nnz, i, j;
1758b68353e5Sstefano_zampini   PetscBool     done;
1759b68353e5Sstefano_zampini 
1760b68353e5Sstefano_zampini   PetscFunctionBegin;
17619566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &nr, &nc));
1762b68353e5Sstefano_zampini   if (reuse == MAT_REUSE_MATRIX) {
1763b68353e5Sstefano_zampini     PetscInt rnr;
1764b68353e5Sstefano_zampini 
17659566063dSJacob Faibussowitsch     PetscCall(MatGetRowIJ(*newmat, 0, PETSC_FALSE, PETSC_FALSE, &rnr, (const PetscInt **)&ii, (const PetscInt **)&jj, &done));
176628b400f6SJacob Faibussowitsch     PetscCheck(done, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "MatGetRowIJ");
176708401ef6SPierre Jolivet     PetscCheck(rnr == nr, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Cannot reuse matrix, wrong number of rows");
17689566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJGetArray(*newmat, &vv));
1769b68353e5Sstefano_zampini   }
1770b68353e5Sstefano_zampini   /* extract CSR for nested SeqAIJ matrices */
1771b68353e5Sstefano_zampini   nnz = 0;
17729566063dSJacob Faibussowitsch   PetscCall(PetscCalloc4(nest->nr * nest->nc, &aii, nest->nr * nest->nc, &ajj, nest->nr * nest->nc, &avv, nest->nr * nest->nc, &trans));
1773b68353e5Sstefano_zampini   for (i = 0; i < nest->nr; ++i) {
1774b68353e5Sstefano_zampini     for (j = 0; j < nest->nc; ++j) {
1775b68353e5Sstefano_zampini       Mat B = nest->m[i][j];
1776b68353e5Sstefano_zampini       if (B) {
1777b68353e5Sstefano_zampini         PetscScalar *naa;
1778b68353e5Sstefano_zampini         PetscInt    *nii, *njj, nnr;
177923875855Sstefano_zampini         PetscBool    istrans;
1780b68353e5Sstefano_zampini 
1781013e2dc7SBarry Smith         PetscCall(PetscObjectTypeCompare((PetscObject)B, MATTRANSPOSEVIRTUAL, &istrans));
178223875855Sstefano_zampini         if (istrans) {
178323875855Sstefano_zampini           Mat Bt;
178423875855Sstefano_zampini 
17859566063dSJacob Faibussowitsch           PetscCall(MatTransposeGetMat(B, &Bt));
17869566063dSJacob Faibussowitsch           PetscCall(MatTranspose(Bt, MAT_INITIAL_MATRIX, &trans[i * nest->nc + j]));
178723875855Sstefano_zampini           B = trans[i * nest->nc + j];
1788013e2dc7SBarry Smith         } else {
1789013e2dc7SBarry Smith           PetscCall(PetscObjectTypeCompare((PetscObject)B, MATHERMITIANTRANSPOSEVIRTUAL, &istrans));
1790013e2dc7SBarry Smith           if (istrans) {
1791013e2dc7SBarry Smith             Mat Bt;
1792013e2dc7SBarry Smith 
1793013e2dc7SBarry Smith             PetscCall(MatHermitianTransposeGetMat(B, &Bt));
1794013e2dc7SBarry Smith             PetscCall(MatHermitianTranspose(Bt, MAT_INITIAL_MATRIX, &trans[i * nest->nc + j]));
1795013e2dc7SBarry Smith             B = trans[i * nest->nc + j];
1796013e2dc7SBarry Smith           }
179723875855Sstefano_zampini         }
17989566063dSJacob Faibussowitsch         PetscCall(MatGetRowIJ(B, 0, PETSC_FALSE, PETSC_FALSE, &nnr, (const PetscInt **)&nii, (const PetscInt **)&njj, &done));
179928b400f6SJacob Faibussowitsch         PetscCheck(done, PetscObjectComm((PetscObject)B), PETSC_ERR_PLIB, "MatGetRowIJ");
18009566063dSJacob Faibussowitsch         PetscCall(MatSeqAIJGetArray(B, &naa));
1801b68353e5Sstefano_zampini         nnz += nii[nnr];
1802b68353e5Sstefano_zampini 
1803b68353e5Sstefano_zampini         aii[i * nest->nc + j] = nii;
1804b68353e5Sstefano_zampini         ajj[i * nest->nc + j] = njj;
1805b68353e5Sstefano_zampini         avv[i * nest->nc + j] = naa;
1806b68353e5Sstefano_zampini       }
1807b68353e5Sstefano_zampini     }
1808b68353e5Sstefano_zampini   }
1809b68353e5Sstefano_zampini   if (reuse != MAT_REUSE_MATRIX) {
18109566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nr + 1, &ii));
18119566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nnz, &jj));
18129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nnz, &vv));
1813b68353e5Sstefano_zampini   } else {
181408401ef6SPierre Jolivet     PetscCheck(nnz == ii[nr], PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "Cannot reuse matrix, wrong number of nonzeros");
1815b68353e5Sstefano_zampini   }
1816b68353e5Sstefano_zampini 
1817b68353e5Sstefano_zampini   /* new row pointer */
18189566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(ii, nr + 1));
1819b68353e5Sstefano_zampini   for (i = 0; i < nest->nr; ++i) {
1820b68353e5Sstefano_zampini     PetscInt ncr, rst;
1821b68353e5Sstefano_zampini 
18229566063dSJacob Faibussowitsch     PetscCall(ISStrideGetInfo(nest->isglobal.row[i], &rst, NULL));
18239566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(nest->isglobal.row[i], &ncr));
1824b68353e5Sstefano_zampini     for (j = 0; j < nest->nc; ++j) {
1825b68353e5Sstefano_zampini       if (aii[i * nest->nc + j]) {
1826b68353e5Sstefano_zampini         PetscInt *nii = aii[i * nest->nc + j];
1827b68353e5Sstefano_zampini         PetscInt  ir;
1828b68353e5Sstefano_zampini 
1829b68353e5Sstefano_zampini         for (ir = rst; ir < ncr + rst; ++ir) {
1830b68353e5Sstefano_zampini           ii[ir + 1] += nii[1] - nii[0];
1831b68353e5Sstefano_zampini           nii++;
1832b68353e5Sstefano_zampini         }
1833b68353e5Sstefano_zampini       }
1834b68353e5Sstefano_zampini     }
1835b68353e5Sstefano_zampini   }
1836b68353e5Sstefano_zampini   for (i = 0; i < nr; i++) ii[i + 1] += ii[i];
1837b68353e5Sstefano_zampini 
1838b68353e5Sstefano_zampini   /* construct CSR for the new matrix */
18399566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nr, &ci));
1840b68353e5Sstefano_zampini   for (i = 0; i < nest->nr; ++i) {
1841b68353e5Sstefano_zampini     PetscInt ncr, rst;
1842b68353e5Sstefano_zampini 
18439566063dSJacob Faibussowitsch     PetscCall(ISStrideGetInfo(nest->isglobal.row[i], &rst, NULL));
18449566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(nest->isglobal.row[i], &ncr));
1845b68353e5Sstefano_zampini     for (j = 0; j < nest->nc; ++j) {
1846b68353e5Sstefano_zampini       if (aii[i * nest->nc + j]) {
1847b68353e5Sstefano_zampini         PetscScalar *nvv = avv[i * nest->nc + j];
1848b68353e5Sstefano_zampini         PetscInt    *nii = aii[i * nest->nc + j];
1849b68353e5Sstefano_zampini         PetscInt    *njj = ajj[i * nest->nc + j];
1850b68353e5Sstefano_zampini         PetscInt     ir, cst;
1851b68353e5Sstefano_zampini 
18529566063dSJacob Faibussowitsch         PetscCall(ISStrideGetInfo(nest->isglobal.col[j], &cst, NULL));
1853b68353e5Sstefano_zampini         for (ir = rst; ir < ncr + rst; ++ir) {
1854b68353e5Sstefano_zampini           PetscInt ij, rsize = nii[1] - nii[0], ist = ii[ir] + ci[ir];
1855b68353e5Sstefano_zampini 
1856b68353e5Sstefano_zampini           for (ij = 0; ij < rsize; ij++) {
1857b68353e5Sstefano_zampini             jj[ist + ij] = *njj + cst;
1858b68353e5Sstefano_zampini             vv[ist + ij] = *nvv;
1859b68353e5Sstefano_zampini             njj++;
1860b68353e5Sstefano_zampini             nvv++;
1861b68353e5Sstefano_zampini           }
1862b68353e5Sstefano_zampini           ci[ir] += rsize;
1863b68353e5Sstefano_zampini           nii++;
1864b68353e5Sstefano_zampini         }
1865b68353e5Sstefano_zampini       }
1866b68353e5Sstefano_zampini     }
1867b68353e5Sstefano_zampini   }
18689566063dSJacob Faibussowitsch   PetscCall(PetscFree(ci));
1869b68353e5Sstefano_zampini 
1870b68353e5Sstefano_zampini   /* restore info */
1871b68353e5Sstefano_zampini   for (i = 0; i < nest->nr; ++i) {
1872b68353e5Sstefano_zampini     for (j = 0; j < nest->nc; ++j) {
1873b68353e5Sstefano_zampini       Mat B = nest->m[i][j];
1874b68353e5Sstefano_zampini       if (B) {
1875b68353e5Sstefano_zampini         PetscInt nnr = 0, k = i * nest->nc + j;
187623875855Sstefano_zampini 
187723875855Sstefano_zampini         B = (trans[k] ? trans[k] : B);
18789566063dSJacob Faibussowitsch         PetscCall(MatRestoreRowIJ(B, 0, PETSC_FALSE, PETSC_FALSE, &nnr, (const PetscInt **)&aii[k], (const PetscInt **)&ajj[k], &done));
187928b400f6SJacob Faibussowitsch         PetscCheck(done, PetscObjectComm((PetscObject)B), PETSC_ERR_PLIB, "MatRestoreRowIJ");
18809566063dSJacob Faibussowitsch         PetscCall(MatSeqAIJRestoreArray(B, &avv[k]));
18819566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&trans[k]));
1882b68353e5Sstefano_zampini       }
1883b68353e5Sstefano_zampini     }
1884b68353e5Sstefano_zampini   }
18859566063dSJacob Faibussowitsch   PetscCall(PetscFree4(aii, ajj, avv, trans));
1886b68353e5Sstefano_zampini 
1887b68353e5Sstefano_zampini   /* finalize newmat */
1888b68353e5Sstefano_zampini   if (reuse == MAT_INITIAL_MATRIX) {
18899566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A), nr, nc, ii, jj, vv, newmat));
1890b68353e5Sstefano_zampini   } else if (reuse == MAT_INPLACE_MATRIX) {
1891b68353e5Sstefano_zampini     Mat B;
1892b68353e5Sstefano_zampini 
18939566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A), nr, nc, ii, jj, vv, &B));
18949566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
1895b68353e5Sstefano_zampini   }
18969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*newmat, MAT_FINAL_ASSEMBLY));
18979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*newmat, MAT_FINAL_ASSEMBLY));
1898b68353e5Sstefano_zampini   {
1899b68353e5Sstefano_zampini     Mat_SeqAIJ *a = (Mat_SeqAIJ *)((*newmat)->data);
1900b68353e5Sstefano_zampini     a->free_a     = PETSC_TRUE;
1901b68353e5Sstefano_zampini     a->free_ij    = PETSC_TRUE;
1902b68353e5Sstefano_zampini   }
19033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1904b68353e5Sstefano_zampini }
1905b68353e5Sstefano_zampini 
1906d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatAXPY_Dense_Nest(Mat Y, PetscScalar a, Mat X)
1907d71ae5a4SJacob Faibussowitsch {
1908be705e3aSPierre Jolivet   Mat_Nest *nest = (Mat_Nest *)X->data;
1909be705e3aSPierre Jolivet   PetscInt  i, j, k, rstart;
1910be705e3aSPierre Jolivet   PetscBool flg;
1911be705e3aSPierre Jolivet 
1912be705e3aSPierre Jolivet   PetscFunctionBegin;
1913be705e3aSPierre Jolivet   /* Fill by row */
1914be705e3aSPierre Jolivet   for (j = 0; j < nest->nc; ++j) {
1915be705e3aSPierre Jolivet     /* Using global column indices and ISAllGather() is not scalable. */
1916be705e3aSPierre Jolivet     IS              bNis;
1917be705e3aSPierre Jolivet     PetscInt        bN;
1918be705e3aSPierre Jolivet     const PetscInt *bNindices;
19199566063dSJacob Faibussowitsch     PetscCall(ISAllGather(nest->isglobal.col[j], &bNis));
19209566063dSJacob Faibussowitsch     PetscCall(ISGetSize(bNis, &bN));
19219566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(bNis, &bNindices));
1922be705e3aSPierre Jolivet     for (i = 0; i < nest->nr; ++i) {
1923fd8a7442SPierre Jolivet       Mat             B = nest->m[i][j], D = NULL;
1924be705e3aSPierre Jolivet       PetscInt        bm, br;
1925be705e3aSPierre Jolivet       const PetscInt *bmindices;
1926be705e3aSPierre Jolivet       if (!B) continue;
1927013e2dc7SBarry Smith       PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
1928be705e3aSPierre Jolivet       if (flg) {
1929cac4c232SBarry Smith         PetscTryMethod(B, "MatTransposeGetMat_C", (Mat, Mat *), (B, &D));
1930cac4c232SBarry Smith         PetscTryMethod(B, "MatHermitianTransposeGetMat_C", (Mat, Mat *), (B, &D));
19319566063dSJacob Faibussowitsch         PetscCall(MatConvert(B, ((PetscObject)D)->type_name, MAT_INITIAL_MATRIX, &D));
1932be705e3aSPierre Jolivet         B = D;
1933be705e3aSPierre Jolivet       }
19349566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQSBAIJ, MATMPISBAIJ, ""));
1935be705e3aSPierre Jolivet       if (flg) {
1936fd8a7442SPierre Jolivet         if (D) PetscCall(MatConvert(D, MATBAIJ, MAT_INPLACE_MATRIX, &D));
1937fd8a7442SPierre Jolivet         else PetscCall(MatConvert(B, MATBAIJ, MAT_INITIAL_MATRIX, &D));
1938be705e3aSPierre Jolivet         B = D;
1939be705e3aSPierre Jolivet       }
19409566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(nest->isglobal.row[i], &bm));
19419566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(nest->isglobal.row[i], &bmindices));
19429566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRange(B, &rstart, NULL));
1943be705e3aSPierre Jolivet       for (br = 0; br < bm; ++br) {
1944be705e3aSPierre Jolivet         PetscInt           row = bmindices[br], brncols, *cols;
1945be705e3aSPierre Jolivet         const PetscInt    *brcols;
1946be705e3aSPierre Jolivet         const PetscScalar *brcoldata;
1947be705e3aSPierre Jolivet         PetscScalar       *vals = NULL;
19489566063dSJacob Faibussowitsch         PetscCall(MatGetRow(B, br + rstart, &brncols, &brcols, &brcoldata));
19499566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(brncols, &cols));
1950be705e3aSPierre Jolivet         for (k = 0; k < brncols; k++) cols[k] = bNindices[brcols[k]];
1951be705e3aSPierre Jolivet         /*
1952be705e3aSPierre Jolivet           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1953be705e3aSPierre Jolivet           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1954be705e3aSPierre Jolivet          */
1955be705e3aSPierre Jolivet         if (a != 1.0) {
19569566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(brncols, &vals));
1957be705e3aSPierre Jolivet           for (k = 0; k < brncols; k++) vals[k] = a * brcoldata[k];
19589566063dSJacob Faibussowitsch           PetscCall(MatSetValues(Y, 1, &row, brncols, cols, vals, ADD_VALUES));
19599566063dSJacob Faibussowitsch           PetscCall(PetscFree(vals));
1960be705e3aSPierre Jolivet         } else {
19619566063dSJacob Faibussowitsch           PetscCall(MatSetValues(Y, 1, &row, brncols, cols, brcoldata, ADD_VALUES));
1962be705e3aSPierre Jolivet         }
19639566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(B, br + rstart, &brncols, &brcols, &brcoldata));
19649566063dSJacob Faibussowitsch         PetscCall(PetscFree(cols));
1965be705e3aSPierre Jolivet       }
1966fd8a7442SPierre Jolivet       if (D) PetscCall(MatDestroy(&D));
19679566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(nest->isglobal.row[i], &bmindices));
1968be705e3aSPierre Jolivet     }
19699566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(bNis, &bNindices));
19709566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&bNis));
1971be705e3aSPierre Jolivet   }
19729566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Y, MAT_FINAL_ASSEMBLY));
19739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Y, MAT_FINAL_ASSEMBLY));
19743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1975be705e3aSPierre Jolivet }
1976be705e3aSPierre Jolivet 
197766976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_Nest_AIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
1978d71ae5a4SJacob Faibussowitsch {
1979629c3df2SDmitry Karpeev   Mat_Nest   *nest = (Mat_Nest *)A->data;
1980e30678d3SPierre Jolivet   PetscInt    m, n, M, N, i, j, k, *dnnz, *onnz = NULL, rstart, cstart, cend;
1981b68353e5Sstefano_zampini   PetscMPIInt size;
1982629c3df2SDmitry Karpeev   Mat         C;
1983629c3df2SDmitry Karpeev 
1984629c3df2SDmitry Karpeev   PetscFunctionBegin;
19859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1986b68353e5Sstefano_zampini   if (size == 1) { /* look for a special case with SeqAIJ matrices and strided-1, contiguous, blocks */
1987b68353e5Sstefano_zampini     PetscInt  nf;
1988b68353e5Sstefano_zampini     PetscBool fast;
1989b68353e5Sstefano_zampini 
19909566063dSJacob Faibussowitsch     PetscCall(PetscStrcmp(newtype, MATAIJ, &fast));
199148a46eb9SPierre Jolivet     if (!fast) PetscCall(PetscStrcmp(newtype, MATSEQAIJ, &fast));
1992b68353e5Sstefano_zampini     for (i = 0; i < nest->nr && fast; ++i) {
1993b68353e5Sstefano_zampini       for (j = 0; j < nest->nc && fast; ++j) {
1994b68353e5Sstefano_zampini         Mat B = nest->m[i][j];
1995b68353e5Sstefano_zampini         if (B) {
19969566063dSJacob Faibussowitsch           PetscCall(PetscObjectTypeCompare((PetscObject)B, MATSEQAIJ, &fast));
199723875855Sstefano_zampini           if (!fast) {
199823875855Sstefano_zampini             PetscBool istrans;
199923875855Sstefano_zampini 
2000013e2dc7SBarry Smith             PetscCall(PetscObjectTypeCompare((PetscObject)B, MATTRANSPOSEVIRTUAL, &istrans));
200123875855Sstefano_zampini             if (istrans) {
200223875855Sstefano_zampini               Mat Bt;
200323875855Sstefano_zampini 
20049566063dSJacob Faibussowitsch               PetscCall(MatTransposeGetMat(B, &Bt));
20059566063dSJacob Faibussowitsch               PetscCall(PetscObjectTypeCompare((PetscObject)Bt, MATSEQAIJ, &fast));
2006013e2dc7SBarry Smith             } else {
2007013e2dc7SBarry Smith               PetscCall(PetscObjectTypeCompare((PetscObject)B, MATHERMITIANTRANSPOSEVIRTUAL, &istrans));
2008013e2dc7SBarry Smith               if (istrans) {
2009013e2dc7SBarry Smith                 Mat Bt;
2010013e2dc7SBarry Smith 
2011013e2dc7SBarry Smith                 PetscCall(MatHermitianTransposeGetMat(B, &Bt));
2012013e2dc7SBarry Smith                 PetscCall(PetscObjectTypeCompare((PetscObject)Bt, MATSEQAIJ, &fast));
2013013e2dc7SBarry Smith               }
201423875855Sstefano_zampini             }
2015b68353e5Sstefano_zampini           }
2016b68353e5Sstefano_zampini         }
2017b68353e5Sstefano_zampini       }
2018b68353e5Sstefano_zampini     }
2019b68353e5Sstefano_zampini     for (i = 0, nf = 0; i < nest->nr && fast; ++i) {
20209566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)nest->isglobal.row[i], ISSTRIDE, &fast));
2021b68353e5Sstefano_zampini       if (fast) {
2022b68353e5Sstefano_zampini         PetscInt f, s;
2023b68353e5Sstefano_zampini 
20249566063dSJacob Faibussowitsch         PetscCall(ISStrideGetInfo(nest->isglobal.row[i], &f, &s));
20259371c9d4SSatish Balay         if (f != nf || s != 1) {
20269371c9d4SSatish Balay           fast = PETSC_FALSE;
20279371c9d4SSatish Balay         } else {
20289566063dSJacob Faibussowitsch           PetscCall(ISGetSize(nest->isglobal.row[i], &f));
2029b68353e5Sstefano_zampini           nf += f;
2030b68353e5Sstefano_zampini         }
2031b68353e5Sstefano_zampini       }
2032b68353e5Sstefano_zampini     }
2033b68353e5Sstefano_zampini     for (i = 0, nf = 0; i < nest->nc && fast; ++i) {
20349566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)nest->isglobal.col[i], ISSTRIDE, &fast));
2035b68353e5Sstefano_zampini       if (fast) {
2036b68353e5Sstefano_zampini         PetscInt f, s;
2037b68353e5Sstefano_zampini 
20389566063dSJacob Faibussowitsch         PetscCall(ISStrideGetInfo(nest->isglobal.col[i], &f, &s));
20399371c9d4SSatish Balay         if (f != nf || s != 1) {
20409371c9d4SSatish Balay           fast = PETSC_FALSE;
20419371c9d4SSatish Balay         } else {
20429566063dSJacob Faibussowitsch           PetscCall(ISGetSize(nest->isglobal.col[i], &f));
2043b68353e5Sstefano_zampini           nf += f;
2044b68353e5Sstefano_zampini         }
2045b68353e5Sstefano_zampini       }
2046b68353e5Sstefano_zampini     }
2047b68353e5Sstefano_zampini     if (fast) {
20489566063dSJacob Faibussowitsch       PetscCall(MatConvert_Nest_SeqAIJ_fast(A, newtype, reuse, newmat));
20493ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
2050b68353e5Sstefano_zampini     }
2051b68353e5Sstefano_zampini   }
20529566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
20539566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
20549566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRangeColumn(A, &cstart, &cend));
2055d1487292SPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) C = *newmat;
2056d1487292SPierre Jolivet   else {
20579566063dSJacob Faibussowitsch     PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &C));
20589566063dSJacob Faibussowitsch     PetscCall(MatSetType(C, newtype));
20599566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(C, m, n, M, N));
2060629c3df2SDmitry Karpeev   }
20619566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * m, &dnnz));
2062e30678d3SPierre Jolivet   if (m) {
2063629c3df2SDmitry Karpeev     onnz = dnnz + m;
2064629c3df2SDmitry Karpeev     for (k = 0; k < m; k++) {
2065629c3df2SDmitry Karpeev       dnnz[k] = 0;
2066629c3df2SDmitry Karpeev       onnz[k] = 0;
2067629c3df2SDmitry Karpeev     }
2068e30678d3SPierre Jolivet   }
2069629c3df2SDmitry Karpeev   for (j = 0; j < nest->nc; ++j) {
2070629c3df2SDmitry Karpeev     IS              bNis;
2071629c3df2SDmitry Karpeev     PetscInt        bN;
2072629c3df2SDmitry Karpeev     const PetscInt *bNindices;
2073fd8a7442SPierre Jolivet     PetscBool       flg;
2074629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
20759566063dSJacob Faibussowitsch     PetscCall(ISAllGather(nest->isglobal.col[j], &bNis));
20769566063dSJacob Faibussowitsch     PetscCall(ISGetSize(bNis, &bN));
20779566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(bNis, &bNindices));
2078629c3df2SDmitry Karpeev     for (i = 0; i < nest->nr; ++i) {
2079629c3df2SDmitry Karpeev       PetscSF         bmsf;
2080649b366bSFande Kong       PetscSFNode    *iremote;
2081fd8a7442SPierre Jolivet       Mat             B = nest->m[i][j], D = NULL;
2082649b366bSFande Kong       PetscInt        bm, *sub_dnnz, *sub_onnz, br;
2083629c3df2SDmitry Karpeev       const PetscInt *bmindices;
2084629c3df2SDmitry Karpeev       if (!B) continue;
20859566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(nest->isglobal.row[i], &bm));
20869566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(nest->isglobal.row[i], &bmindices));
20879566063dSJacob Faibussowitsch       PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf));
20889566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(bm, &iremote));
20899566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(bm, &sub_dnnz));
20909566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(bm, &sub_onnz));
2091649b366bSFande Kong       for (k = 0; k < bm; ++k) {
2092649b366bSFande Kong         sub_dnnz[k] = 0;
2093649b366bSFande Kong         sub_onnz[k] = 0;
2094649b366bSFande Kong       }
2095dead4d76SPierre Jolivet       PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATTRANSPOSEVIRTUAL, MATHERMITIANTRANSPOSEVIRTUAL, ""));
2096fd8a7442SPierre Jolivet       if (flg) {
2097fd8a7442SPierre Jolivet         PetscTryMethod(B, "MatTransposeGetMat_C", (Mat, Mat *), (B, &D));
2098fd8a7442SPierre Jolivet         PetscTryMethod(B, "MatHermitianTransposeGetMat_C", (Mat, Mat *), (B, &D));
2099fd8a7442SPierre Jolivet         PetscCall(MatConvert(B, ((PetscObject)D)->type_name, MAT_INITIAL_MATRIX, &D));
2100fd8a7442SPierre Jolivet         B = D;
2101fd8a7442SPierre Jolivet       }
2102fd8a7442SPierre Jolivet       PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &flg, MATSEQSBAIJ, MATMPISBAIJ, ""));
2103fd8a7442SPierre Jolivet       if (flg) {
2104fd8a7442SPierre Jolivet         if (D) PetscCall(MatConvert(D, MATBAIJ, MAT_INPLACE_MATRIX, &D));
2105fd8a7442SPierre Jolivet         else PetscCall(MatConvert(B, MATBAIJ, MAT_INITIAL_MATRIX, &D));
2106fd8a7442SPierre Jolivet         B = D;
2107fd8a7442SPierre Jolivet       }
2108629c3df2SDmitry Karpeev       /*
2109629c3df2SDmitry Karpeev        Locate the owners for all of the locally-owned global row indices for this row block.
2110629c3df2SDmitry Karpeev        These determine the roots of PetscSF used to communicate preallocation data to row owners.
2111629c3df2SDmitry Karpeev        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
2112629c3df2SDmitry Karpeev        */
21139566063dSJacob Faibussowitsch       PetscCall(MatGetOwnershipRange(B, &rstart, NULL));
2114629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
2115131c27b5Sprj-         PetscInt        row = bmindices[br], brncols, col;
2116629c3df2SDmitry Karpeev         const PetscInt *brcols;
2117a4b3d3acSMatthew G Knepley         PetscInt        rowrel   = 0; /* row's relative index on its owner rank */
2118131c27b5Sprj-         PetscMPIInt     rowowner = 0;
21199566063dSJacob Faibussowitsch         PetscCall(PetscLayoutFindOwnerIndex(A->rmap, row, &rowowner, &rowrel));
2120649b366bSFande Kong         /* how many roots  */
21219371c9d4SSatish Balay         iremote[br].rank  = rowowner;
21229371c9d4SSatish Balay         iremote[br].index = rowrel; /* edge from bmdnnz to dnnz */
2123649b366bSFande Kong         /* get nonzero pattern */
21249566063dSJacob Faibussowitsch         PetscCall(MatGetRow(B, br + rstart, &brncols, &brcols, NULL));
2125629c3df2SDmitry Karpeev         for (k = 0; k < brncols; k++) {
2126629c3df2SDmitry Karpeev           col = bNindices[brcols[k]];
2127649b366bSFande Kong           if (col >= A->cmap->range[rowowner] && col < A->cmap->range[rowowner + 1]) {
2128649b366bSFande Kong             sub_dnnz[br]++;
2129649b366bSFande Kong           } else {
2130649b366bSFande Kong             sub_onnz[br]++;
2131649b366bSFande Kong           }
2132629c3df2SDmitry Karpeev         }
21339566063dSJacob Faibussowitsch         PetscCall(MatRestoreRow(B, br + rstart, &brncols, &brcols, NULL));
2134629c3df2SDmitry Karpeev       }
2135fd8a7442SPierre Jolivet       if (D) PetscCall(MatDestroy(&D));
21369566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(nest->isglobal.row[i], &bmindices));
2137629c3df2SDmitry Karpeev       /* bsf will have to take care of disposing of bedges. */
21389566063dSJacob Faibussowitsch       PetscCall(PetscSFSetGraph(bmsf, m, bm, NULL, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
21399566063dSJacob Faibussowitsch       PetscCall(PetscSFReduceBegin(bmsf, MPIU_INT, sub_dnnz, dnnz, MPI_SUM));
21409566063dSJacob Faibussowitsch       PetscCall(PetscSFReduceEnd(bmsf, MPIU_INT, sub_dnnz, dnnz, MPI_SUM));
21419566063dSJacob Faibussowitsch       PetscCall(PetscSFReduceBegin(bmsf, MPIU_INT, sub_onnz, onnz, MPI_SUM));
21429566063dSJacob Faibussowitsch       PetscCall(PetscSFReduceEnd(bmsf, MPIU_INT, sub_onnz, onnz, MPI_SUM));
21439566063dSJacob Faibussowitsch       PetscCall(PetscFree(sub_dnnz));
21449566063dSJacob Faibussowitsch       PetscCall(PetscFree(sub_onnz));
21459566063dSJacob Faibussowitsch       PetscCall(PetscSFDestroy(&bmsf));
2146629c3df2SDmitry Karpeev     }
21479566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(bNis, &bNindices));
21489566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&bNis));
214965a4a0a3Sstefano_zampini   }
215065a4a0a3Sstefano_zampini   /* Resize preallocation if overestimated */
215165a4a0a3Sstefano_zampini   for (i = 0; i < m; i++) {
215265a4a0a3Sstefano_zampini     dnnz[i] = PetscMin(dnnz[i], A->cmap->n);
215365a4a0a3Sstefano_zampini     onnz[i] = PetscMin(onnz[i], A->cmap->N - A->cmap->n);
2154629c3df2SDmitry Karpeev   }
21559566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(C, 0, dnnz));
21569566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(C, 0, dnnz, 0, onnz));
21579566063dSJacob Faibussowitsch   PetscCall(PetscFree(dnnz));
21589566063dSJacob Faibussowitsch   PetscCall(MatAXPY_Dense_Nest(C, 1.0, A));
2159d1487292SPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
21609566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &C));
2161d1487292SPierre Jolivet   } else *newmat = C;
21623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2163be705e3aSPierre Jolivet }
2164629c3df2SDmitry Karpeev 
216566976f2fSJacob Faibussowitsch static PetscErrorCode MatConvert_Nest_Dense(Mat A, MatType newtype, MatReuse reuse, Mat *newmat)
2166d71ae5a4SJacob Faibussowitsch {
2167629c3df2SDmitry Karpeev   Mat      B;
2168be705e3aSPierre Jolivet   PetscInt m, n, M, N;
2169be705e3aSPierre Jolivet 
2170be705e3aSPierre Jolivet   PetscFunctionBegin;
21719566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &M, &N));
21729566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
2173be705e3aSPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
2174be705e3aSPierre Jolivet     B = *newmat;
21759566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(B));
2176be705e3aSPierre Jolivet   } else {
21779566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PetscObjectComm((PetscObject)A), m, PETSC_DECIDE, M, N, NULL, &B));
2178629c3df2SDmitry Karpeev   }
21799566063dSJacob Faibussowitsch   PetscCall(MatAXPY_Dense_Nest(B, 1.0, A));
2180be705e3aSPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
21819566063dSJacob Faibussowitsch     PetscCall(MatHeaderReplace(A, &B));
2182be705e3aSPierre Jolivet   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
21833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2184629c3df2SDmitry Karpeev }
2185629c3df2SDmitry Karpeev 
218666976f2fSJacob Faibussowitsch static PetscErrorCode MatHasOperation_Nest(Mat mat, MatOperation op, PetscBool *has)
2187d71ae5a4SJacob Faibussowitsch {
21888b7d3b4bSBarry Smith   Mat_Nest    *bA = (Mat_Nest *)mat->data;
21893c6db4c4SPierre Jolivet   MatOperation opAdd;
21908b7d3b4bSBarry Smith   PetscInt     i, j, nr = bA->nr, nc = bA->nc;
21918b7d3b4bSBarry Smith   PetscBool    flg;
219252c5f739Sprj-   PetscFunctionBegin;
21938b7d3b4bSBarry Smith 
219452c5f739Sprj-   *has = PETSC_FALSE;
21953c6db4c4SPierre Jolivet   if (op == MATOP_MULT || op == MATOP_MULT_ADD || op == MATOP_MULT_TRANSPOSE || op == MATOP_MULT_TRANSPOSE_ADD) {
21963c6db4c4SPierre Jolivet     opAdd = (op == MATOP_MULT || op == MATOP_MULT_ADD ? MATOP_MULT_ADD : MATOP_MULT_TRANSPOSE_ADD);
21978b7d3b4bSBarry Smith     for (j = 0; j < nc; j++) {
21988b7d3b4bSBarry Smith       for (i = 0; i < nr; i++) {
21998b7d3b4bSBarry Smith         if (!bA->m[i][j]) continue;
22009566063dSJacob Faibussowitsch         PetscCall(MatHasOperation(bA->m[i][j], opAdd, &flg));
22013ba16761SJacob Faibussowitsch         if (!flg) PetscFunctionReturn(PETSC_SUCCESS);
22028b7d3b4bSBarry Smith       }
22038b7d3b4bSBarry Smith     }
22048b7d3b4bSBarry Smith   }
22053c6db4c4SPierre Jolivet   if (((void **)mat->ops)[op]) *has = PETSC_TRUE;
22063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22078b7d3b4bSBarry Smith }
22088b7d3b4bSBarry Smith 
2209659c6bb0SJed Brown /*MC
22102ef1f0ffSBarry Smith   MATNEST -  "nest" - Matrix type consisting of nested submatrices, each stored separately.
2211659c6bb0SJed Brown 
2212659c6bb0SJed Brown   Level: intermediate
2213659c6bb0SJed Brown 
2214659c6bb0SJed Brown   Notes:
221511a5261eSBarry Smith   This matrix type permits scalable use of `PCFIELDSPLIT` and avoids the large memory costs of extracting submatrices.
2216659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
221711a5261eSBarry Smith   It is usually used with `DMCOMPOSITE` and `DMCreateMatrix()`
2218659c6bb0SJed Brown 
22198b7d3b4bSBarry Smith   Each of the submatrices lives on the same MPI communicator as the original nest matrix (though they can have zero
22208b7d3b4bSBarry Smith   rows/columns on some processes.) Thus this is not meant for cases where the submatrices live on far fewer processes
22218b7d3b4bSBarry Smith   than the nest matrix.
22228b7d3b4bSBarry Smith 
22231cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATNEST`, `MatCreate()`, `MatType`, `MatCreateNest()`, `MatNestSetSubMat()`, `MatNestGetSubMat()`,
2224db781477SPatrick Sanan           `VecCreateNest()`, `DMCreateMatrix()`, `DMCOMPOSITE`, `MatNestSetVecType()`, `MatNestGetLocalISs()`,
2225db781477SPatrick Sanan           `MatNestGetISs()`, `MatNestSetSubMats()`, `MatNestGetSubMats()`
2226659c6bb0SJed Brown M*/
2227d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
2228d71ae5a4SJacob Faibussowitsch {
2229c8883902SJed Brown   Mat_Nest *s;
2230c8883902SJed Brown 
2231c8883902SJed Brown   PetscFunctionBegin;
22324dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&s));
2233c8883902SJed Brown   A->data = (void *)s;
2234e7c19651SJed Brown 
2235e7c19651SJed Brown   s->nr            = -1;
2236e7c19651SJed Brown   s->nc            = -1;
22370298fd71SBarry Smith   s->m             = NULL;
2238e7c19651SJed Brown   s->splitassembly = PETSC_FALSE;
2239c8883902SJed Brown 
22409566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(A->ops, sizeof(*A->ops)));
224126fbe8dcSKarl Rupp 
2242c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
22439194d70fSJed Brown   A->ops->multadd               = MatMultAdd_Nest;
2244c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
22459194d70fSJed Brown   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
2246f8170845SAlex Fikl   A->ops->transpose             = MatTranspose_Nest;
2247c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
2248c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
2249c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
2250c222c20dSDavid Ham   A->ops->copy                  = MatCopy_Nest;
22516e76ffeaSPierre Jolivet   A->ops->axpy                  = MatAXPY_Nest;
2252c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
22537dae84e0SHong Zhang   A->ops->createsubmatrix       = MatCreateSubMatrix_Nest;
2254c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
2255c8883902SJed Brown   A->ops->view                  = MatView_Nest;
2256f4259b30SLisandro Dalcin   A->ops->getvecs               = NULL; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
2257c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
2258c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
2259429bac76SJed Brown   A->ops->getdiagonal           = MatGetDiagonal_Nest;
2260429bac76SJed Brown   A->ops->diagonalscale         = MatDiagonalScale_Nest;
2261a061e289SJed Brown   A->ops->scale                 = MatScale_Nest;
2262a061e289SJed Brown   A->ops->shift                 = MatShift_Nest;
226313135bc6SAlex Fikl   A->ops->diagonalset           = MatDiagonalSet_Nest;
2264f8170845SAlex Fikl   A->ops->setrandom             = MatSetRandom_Nest;
22658b7d3b4bSBarry Smith   A->ops->hasoperation          = MatHasOperation_Nest;
2266381b8e50SStefano Zampini   A->ops->missingdiagonal       = MatMissingDiagonal_Nest;
2267c8883902SJed Brown 
2268f4259b30SLisandro Dalcin   A->spptr     = NULL;
2269c8883902SJed Brown   A->assembled = PETSC_FALSE;
2270c8883902SJed Brown 
2271c8883902SJed Brown   /* expose Nest api's */
22729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetSubMat_C", MatNestGetSubMat_Nest));
22739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestSetSubMat_C", MatNestSetSubMat_Nest));
22749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetSubMats_C", MatNestGetSubMats_Nest));
22759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetSize_C", MatNestGetSize_Nest));
22769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetISs_C", MatNestGetISs_Nest));
22779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestGetLocalISs_C", MatNestGetLocalISs_Nest));
22789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestSetVecType_C", MatNestSetVecType_Nest));
22799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatNestSetSubMats_C", MatNestSetSubMats_Nest));
22809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_mpiaij_C", MatConvert_Nest_AIJ));
22819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_seqaij_C", MatConvert_Nest_AIJ));
22829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_aij_C", MatConvert_Nest_AIJ));
22839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_is_C", MatConvert_Nest_IS));
22849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_mpidense_C", MatConvert_Nest_Dense));
22859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConvert_nest_seqdense_C", MatConvert_Nest_Dense));
22869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_nest_seqdense_C", MatProductSetFromOptions_Nest_Dense));
22879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_nest_mpidense_C", MatProductSetFromOptions_Nest_Dense));
22889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_nest_dense_C", MatProductSetFromOptions_Nest_Dense));
2289c8883902SJed Brown 
22909566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATNEST));
22913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2292c8883902SJed Brown }
2293