xref: /petsc/src/mat/impls/nest/matnest.c (revision bb97c47c1ebd420044fb986bcb47505b5e22f105)
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 */
128188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N)
13d8588912SDave May {
14d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
158188e55aSJed Brown   PetscInt       i,j;
16d8588912SDave May   PetscErrorCode ierr;
17d8588912SDave May 
18d8588912SDave May   PetscFunctionBegin;
198188e55aSJed Brown   *m = *n = *M = *N = 0;
208188e55aSJed Brown   for (i=0; i<bA->nr; i++) {  /* rows */
218188e55aSJed Brown     PetscInt sm,sM;
228188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr);
238188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr);
248188e55aSJed Brown     *m  += sm;
258188e55aSJed Brown     *M  += sM;
26d8588912SDave May   }
278188e55aSJed Brown   for (j=0; j<bA->nc; j++) {  /* cols */
288188e55aSJed Brown     PetscInt sn,sN;
298188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr);
308188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr);
318188e55aSJed Brown     *n  += sn;
328188e55aSJed Brown     *N  += sN;
33d8588912SDave May   }
34d8588912SDave May   PetscFunctionReturn(0);
35d8588912SDave May }
36d8588912SDave May 
37d8588912SDave May /* operations */
38207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
39d8588912SDave May {
40d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
41207556f9SJed Brown   Vec            *bx = bA->right,*by = bA->left;
42207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
43d8588912SDave May   PetscErrorCode ierr;
44d8588912SDave May 
45d8588912SDave May   PetscFunctionBegin;
46207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
47207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
48207556f9SJed Brown   for (i=0; i<nr; i++) {
49d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
50207556f9SJed Brown     for (j=0; j<nc; j++) {
51207556f9SJed Brown       if (!bA->m[i][j]) continue;
52d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
53d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
54d8588912SDave May     }
55d8588912SDave May   }
56207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
57207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
58d8588912SDave May   PetscFunctionReturn(0);
59d8588912SDave May }
60d8588912SDave May 
619194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z)
629194d70fSJed Brown {
639194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
649194d70fSJed Brown   Vec            *bx = bA->right,*bz = bA->left;
659194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
669194d70fSJed Brown   PetscErrorCode ierr;
679194d70fSJed Brown 
689194d70fSJed Brown   PetscFunctionBegin;
699194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
709194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
719194d70fSJed Brown   for (i=0; i<nr; i++) {
729194d70fSJed Brown     if (y != z) {
739194d70fSJed Brown       Vec by;
749194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
759194d70fSJed Brown       ierr = VecCopy(by,bz[i]);CHKERRQ(ierr);
76336d21e7SJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
779194d70fSJed Brown     }
789194d70fSJed Brown     for (j=0; j<nc; j++) {
799194d70fSJed Brown       if (!bA->m[i][j]) continue;
809194d70fSJed Brown       /* y[i] <- y[i] + A[i][j] * x[j] */
819194d70fSJed Brown       ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr);
829194d70fSJed Brown     }
839194d70fSJed Brown   }
849194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
859194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
869194d70fSJed Brown   PetscFunctionReturn(0);
879194d70fSJed Brown }
889194d70fSJed Brown 
8952c5f739Sprj- typedef struct {
9052c5f739Sprj-   Mat          *workC;    /* array of Mat with specific containers depending on the underlying MatMatMult implementation */
9152c5f739Sprj-   PetscScalar  *tarray;   /* buffer for storing all temporary products A[i][j] B[j] */
9252c5f739Sprj-   PetscInt     *dm,*dn,k; /* displacements and number of submatrices */
9352c5f739Sprj- } Nest_Dense;
9452c5f739Sprj- 
9552c5f739Sprj- PETSC_INTERN PetscErrorCode MatMatMultNumeric_Nest_Dense(Mat A,Mat B,Mat C)
9652c5f739Sprj- {
9752c5f739Sprj-   Mat_Nest          *bA = (Mat_Nest*)A->data;
9852c5f739Sprj-   PetscContainer    container;
9952c5f739Sprj-   Nest_Dense        *contents;
1004222ddf1SHong Zhang   Mat               viewB,viewC,seq,productB,workC;
10152c5f739Sprj-   const PetscScalar *barray;
10252c5f739Sprj-   PetscScalar       *carray;
10352c5f739Sprj-   PetscInt          i,j,M,N,nr = bA->nr,nc = bA->nc,ldb,ldc;
10452c5f739Sprj-   PetscErrorCode    ierr;
10552c5f739Sprj- 
10652c5f739Sprj-   PetscFunctionBegin;
10752c5f739Sprj-   ierr = PetscObjectQuery((PetscObject)C,"workC",(PetscObject*)&container);CHKERRQ(ierr);
10852c5f739Sprj-   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exist");
10952c5f739Sprj-   ierr = PetscContainerGetPointer(container,(void**)&contents);CHKERRQ(ierr);
11052c5f739Sprj-   ierr = MatDenseGetLDA(B,&ldb);CHKERRQ(ierr);
11152c5f739Sprj-   ierr = MatDenseGetLDA(C,&ldc);CHKERRQ(ierr);
11252c5f739Sprj-   ierr = MatGetSize(B,NULL,&N);CHKERRQ(ierr);
11352c5f739Sprj-   ierr = MatZeroEntries(C);CHKERRQ(ierr);
11452c5f739Sprj-   ierr = MatDenseGetArrayRead(B,&barray);CHKERRQ(ierr);
11552c5f739Sprj-   ierr = MatDenseGetArray(C,&carray);CHKERRQ(ierr);
11652c5f739Sprj-   for (i=0; i<nr; i++) {
11752c5f739Sprj-     ierr = ISGetSize(bA->isglobal.row[i],&M);CHKERRQ(ierr);
11852c5f739Sprj-     ierr = MatCreateDense(PetscObjectComm((PetscObject)A),contents->dm[i+1]-contents->dm[i],PETSC_DECIDE,M,N,carray+contents->dm[i],&viewC);CHKERRQ(ierr);
11952c5f739Sprj-     ierr = MatDenseGetLocalMatrix(viewC,&seq);CHKERRQ(ierr);
12052c5f739Sprj-     ierr = MatSeqDenseSetLDA(seq,ldc);CHKERRQ(ierr);
12152c5f739Sprj-     for (j=0; j<nc; j++) {
12252c5f739Sprj-       if (!bA->m[i][j]) continue;
12352c5f739Sprj-       ierr = ISGetSize(bA->isglobal.col[j],&M);CHKERRQ(ierr);
12452c5f739Sprj-       ierr = MatCreateDense(PetscObjectComm((PetscObject)A),contents->dn[j+1]-contents->dn[j],PETSC_DECIDE,M,N,(PetscScalar*)(barray+contents->dn[j]),&viewB);CHKERRQ(ierr);
12552c5f739Sprj-       ierr = MatDenseGetLocalMatrix(viewB,&seq);CHKERRQ(ierr);
12652c5f739Sprj-       ierr = MatSeqDenseSetLDA(seq,ldb);CHKERRQ(ierr);
1274222ddf1SHong Zhang 
1284222ddf1SHong Zhang       /* MatMatMultNumeric(bA->m[i][j],viewB,contents->workC[i*nc + j]); */
1294222ddf1SHong Zhang       workC             = contents->workC[i*nc + j];
1304222ddf1SHong Zhang       productB          = workC->product->B;
1314222ddf1SHong Zhang       workC->product->B = viewB; /* use newly created dense matrix viewB */
1324222ddf1SHong Zhang       ierr = (workC->ops->productnumeric)(workC);CHKERRQ(ierr);
13352c5f739Sprj-       ierr = MatDestroy(&viewB);CHKERRQ(ierr);
1344222ddf1SHong Zhang       workC->product->B = productB; /* resume original B */
1354222ddf1SHong Zhang 
13652c5f739Sprj-       /* C[i] <- workC + C[i] */
13752c5f739Sprj-       ierr = MatAXPY(viewC,1.0,contents->workC[i*nc + j],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
13852c5f739Sprj-     }
13952c5f739Sprj-     ierr = MatDestroy(&viewC);CHKERRQ(ierr);
14052c5f739Sprj-   }
14152c5f739Sprj-   ierr = MatDenseRestoreArray(C,&carray);CHKERRQ(ierr);
14252c5f739Sprj-   ierr = MatDenseRestoreArrayRead(B,&barray);CHKERRQ(ierr);
1434222ddf1SHong Zhang 
1444222ddf1SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1454222ddf1SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14652c5f739Sprj-   PetscFunctionReturn(0);
14752c5f739Sprj- }
14852c5f739Sprj- 
14952c5f739Sprj- PetscErrorCode MatNest_DenseDestroy(void *ctx)
15052c5f739Sprj- {
15152c5f739Sprj-   Nest_Dense     *contents = (Nest_Dense*)ctx;
15252c5f739Sprj-   PetscInt       i;
15352c5f739Sprj-   PetscErrorCode ierr;
15452c5f739Sprj- 
15552c5f739Sprj-   PetscFunctionBegin;
15652c5f739Sprj-   ierr = PetscFree(contents->tarray);CHKERRQ(ierr);
15752c5f739Sprj-   for (i=0; i<contents->k; i++) {
15852c5f739Sprj-     ierr = MatDestroy(contents->workC + i);CHKERRQ(ierr);
15952c5f739Sprj-   }
16052c5f739Sprj-   ierr = PetscFree3(contents->dm,contents->dn,contents->workC);CHKERRQ(ierr);
16152c5f739Sprj-   ierr = PetscFree(contents);CHKERRQ(ierr);
16252c5f739Sprj-   PetscFunctionReturn(0);
16352c5f739Sprj- }
16452c5f739Sprj- 
1654222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Nest_Dense(Mat A,Mat B,PetscReal fill,Mat C)
16652c5f739Sprj- {
16752c5f739Sprj-   Mat_Nest          *bA = (Mat_Nest*)A->data;
1684222ddf1SHong Zhang   Mat               viewB,viewSeq,workC;
16952c5f739Sprj-   const PetscScalar *barray;
17052c5f739Sprj-   PetscInt          i,j,M,N,m,nr = bA->nr,nc = bA->nc,maxm = 0,ldb;
17152c5f739Sprj-   PetscContainer    container;
1724222ddf1SHong Zhang   Nest_Dense        *contents=NULL;
17352c5f739Sprj-   PetscErrorCode    ierr;
17452c5f739Sprj- 
17552c5f739Sprj-   PetscFunctionBegin;
1764222ddf1SHong Zhang   if (!C->assembled) {
17752c5f739Sprj-     ierr = MatGetSize(B,NULL,&N);CHKERRQ(ierr);
17852c5f739Sprj-     ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
17952c5f739Sprj-     ierr = MatGetSize(A,&M,NULL);CHKERRQ(ierr);
1804222ddf1SHong Zhang 
1814222ddf1SHong Zhang     ierr = MatSetSizes(C,m,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1824222ddf1SHong Zhang     ierr = MatSetType(C,MATDENSE);CHKERRQ(ierr);
1834222ddf1SHong Zhang     ierr = MatSeqDenseSetPreallocation(C,NULL);CHKERRQ(ierr);
1844222ddf1SHong Zhang     ierr = MatMPIDenseSetPreallocation(C,NULL);CHKERRQ(ierr);
18552c5f739Sprj-   }
18652c5f739Sprj- 
18752c5f739Sprj-   ierr = PetscNew(&contents);CHKERRQ(ierr);
18852c5f739Sprj-   ierr = PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);CHKERRQ(ierr);
18952c5f739Sprj-   ierr = PetscContainerSetPointer(container,contents);CHKERRQ(ierr);
19052c5f739Sprj-   ierr = PetscContainerSetUserDestroy(container,MatNest_DenseDestroy);CHKERRQ(ierr);
1914222ddf1SHong Zhang   ierr = PetscObjectCompose((PetscObject)C,"workC",(PetscObject)container);CHKERRQ(ierr);
19252c5f739Sprj-   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
19352c5f739Sprj-   ierr = PetscCalloc3(nr+1,&contents->dm,nc+1,&contents->dn,nr*nc,&contents->workC);CHKERRQ(ierr);
19452c5f739Sprj-   contents->k = nr*nc;
19552c5f739Sprj-   for (i=0; i<nr; i++) {
19652c5f739Sprj-     ierr = ISGetLocalSize(bA->isglobal.row[i],contents->dm + i+1);CHKERRQ(ierr);
19752c5f739Sprj-     maxm = PetscMax(maxm,contents->dm[i+1]);
19852c5f739Sprj-     contents->dm[i+1] += contents->dm[i];
19952c5f739Sprj-   }
20052c5f739Sprj-   for (i=0; i<nc; i++) {
20152c5f739Sprj-     ierr = ISGetLocalSize(bA->isglobal.col[i],contents->dn + i+1);CHKERRQ(ierr);
20252c5f739Sprj-     contents->dn[i+1] += contents->dn[i];
20352c5f739Sprj-   }
20452c5f739Sprj-   ierr = PetscMalloc1(maxm*N,&contents->tarray);CHKERRQ(ierr);
20552c5f739Sprj-   ierr = MatDenseGetLDA(B,&ldb);CHKERRQ(ierr);
20652c5f739Sprj-   ierr = MatGetSize(B,NULL,&N);CHKERRQ(ierr);
20752c5f739Sprj-   ierr = MatDenseGetArrayRead(B,&barray);CHKERRQ(ierr);
20852c5f739Sprj-   /* loops are permuted compared to MatMatMultNumeric so that viewB is created only once per column of A */
20952c5f739Sprj-   for (j=0; j<nc; j++) {
21052c5f739Sprj-     ierr = ISGetSize(bA->isglobal.col[j],&M);CHKERRQ(ierr);
21152c5f739Sprj-     ierr = MatCreateDense(PetscObjectComm((PetscObject)A),contents->dn[j+1]-contents->dn[j],PETSC_DECIDE,M,N,(PetscScalar*)(barray+contents->dn[j]),&viewB);CHKERRQ(ierr);
21252c5f739Sprj-     ierr = MatDenseGetLocalMatrix(viewB,&viewSeq);CHKERRQ(ierr);
21352c5f739Sprj-     ierr = MatSeqDenseSetLDA(viewSeq,ldb);CHKERRQ(ierr);
21452c5f739Sprj-     for (i=0; i<nr; i++) {
21552c5f739Sprj-       if (!bA->m[i][j]) continue;
21652c5f739Sprj-       /* MatMatMultSymbolic may attach a specific container (depending on MatType of bA->m[i][j]) to workC[i][j] */
2174222ddf1SHong Zhang 
2184222ddf1SHong Zhang       ierr = MatProductCreate(bA->m[i][j],viewB,NULL,&contents->workC[i*nc + j]);CHKERRQ(ierr);
2194222ddf1SHong Zhang       workC = contents->workC[i*nc + j];
2204222ddf1SHong Zhang       ierr = MatProductSetType(workC,MATPRODUCT_AB);CHKERRQ(ierr);
2214222ddf1SHong Zhang       ierr = MatProductSetAlgorithm(workC,"default");CHKERRQ(ierr);
2224222ddf1SHong Zhang       ierr = MatProductSetFill(workC,fill);CHKERRQ(ierr);
2234222ddf1SHong Zhang       ierr = MatProductSetFromOptions(workC);CHKERRQ(ierr);
2244222ddf1SHong Zhang       ierr = MatProductSymbolic(workC);CHKERRQ(ierr);
2254222ddf1SHong Zhang 
2264222ddf1SHong Zhang       ierr = MatDenseGetLocalMatrix(workC,&viewSeq);CHKERRQ(ierr);
22752c5f739Sprj-       /* free the memory allocated in MatMatMultSymbolic, since tarray will be shared by all Mat */
22852c5f739Sprj-       ierr = MatSeqDenseSetPreallocation(viewSeq,contents->tarray);CHKERRQ(ierr);
22952c5f739Sprj-     }
23052c5f739Sprj-     ierr = MatDestroy(&viewB);CHKERRQ(ierr);
23152c5f739Sprj-   }
23252c5f739Sprj-   ierr = MatDenseRestoreArrayRead(B,&barray);CHKERRQ(ierr);
23352c5f739Sprj- 
2344222ddf1SHong Zhang   C->ops->matmultnumeric = MatMatMultNumeric_Nest_Dense;
23552c5f739Sprj-   PetscFunctionReturn(0);
23652c5f739Sprj- }
23752c5f739Sprj- 
2384222ddf1SHong Zhang /* --------------------------------------------------------- */
2394222ddf1SHong Zhang static PetscErrorCode MatProductSetFromOptions_Nest_Dense_AB(Mat C)
2404222ddf1SHong Zhang {
2414222ddf1SHong Zhang   PetscFunctionBegin;
2424222ddf1SHong Zhang   C->ops->matmultsymbolic = MatMatMultSymbolic_Nest_Dense;
2434222ddf1SHong Zhang   C->ops->productsymbolic = MatProductSymbolic_AB;
2444222ddf1SHong Zhang   PetscFunctionReturn(0);
2454222ddf1SHong Zhang }
2464222ddf1SHong Zhang 
2474222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Nest_Dense(Mat C)
24852c5f739Sprj- {
24952c5f739Sprj-   PetscErrorCode ierr;
2504222ddf1SHong Zhang   Mat_Product    *product = C->product;
25152c5f739Sprj- 
25252c5f739Sprj-   PetscFunctionBegin;
2534222ddf1SHong Zhang   if (product->type == MATPRODUCT_AB) {
2544222ddf1SHong Zhang     ierr = MatProductSetFromOptions_Nest_Dense_AB(C);CHKERRQ(ierr);
2554222ddf1SHong Zhang   } else SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"MatProduct type is not supported");
25652c5f739Sprj-   PetscFunctionReturn(0);
25752c5f739Sprj- }
2584222ddf1SHong Zhang /* --------------------------------------------------------- */
25952c5f739Sprj- 
260207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
261d8588912SDave May {
262d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
263207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
264207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
265d8588912SDave May   PetscErrorCode ierr;
266d8588912SDave May 
267d8588912SDave May   PetscFunctionBegin;
268609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
269609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
270207556f9SJed Brown   for (j=0; j<nc; j++) {
271609e31cbSJed Brown     ierr = VecZeroEntries(by[j]);CHKERRQ(ierr);
272609e31cbSJed Brown     for (i=0; i<nr; i++) {
2736c75ac25SJed Brown       if (!bA->m[i][j]) continue;
274609e31cbSJed Brown       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
275609e31cbSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr);
276d8588912SDave May     }
277d8588912SDave May   }
278609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
279609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
280d8588912SDave May   PetscFunctionReturn(0);
281d8588912SDave May }
282d8588912SDave May 
2839194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
2849194d70fSJed Brown {
2859194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
2869194d70fSJed Brown   Vec            *bx = bA->left,*bz = bA->right;
2879194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
2889194d70fSJed Brown   PetscErrorCode ierr;
2899194d70fSJed Brown 
2909194d70fSJed Brown   PetscFunctionBegin;
2919194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
2929194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
2939194d70fSJed Brown   for (j=0; j<nc; j++) {
2949194d70fSJed Brown     if (y != z) {
2959194d70fSJed Brown       Vec by;
2969194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
2979194d70fSJed Brown       ierr = VecCopy(by,bz[j]);CHKERRQ(ierr);
2989194d70fSJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
2999194d70fSJed Brown     }
3009194d70fSJed Brown     for (i=0; i<nr; i++) {
3016c75ac25SJed Brown       if (!bA->m[i][j]) continue;
3029194d70fSJed Brown       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
3039194d70fSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr);
3049194d70fSJed Brown     }
3059194d70fSJed Brown   }
3069194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
3079194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
3089194d70fSJed Brown   PetscFunctionReturn(0);
3099194d70fSJed Brown }
3109194d70fSJed Brown 
311f8170845SAlex Fikl static PetscErrorCode MatTranspose_Nest(Mat A,MatReuse reuse,Mat *B)
312f8170845SAlex Fikl {
313f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data, *bC;
314f8170845SAlex Fikl   Mat            C;
315f8170845SAlex Fikl   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
316f8170845SAlex Fikl   PetscErrorCode ierr;
317f8170845SAlex Fikl 
318f8170845SAlex Fikl   PetscFunctionBegin;
319cf37664fSBarry Smith   if (reuse == MAT_INPLACE_MATRIX && nr != nc) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Square nested matrix only for in-place");
320f8170845SAlex Fikl 
321cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
322f8170845SAlex Fikl     Mat *subs;
323f8170845SAlex Fikl     IS  *is_row,*is_col;
324f8170845SAlex Fikl 
325f8170845SAlex Fikl     ierr = PetscCalloc1(nr * nc,&subs);CHKERRQ(ierr);
326f8170845SAlex Fikl     ierr = PetscMalloc2(nr,&is_row,nc,&is_col);CHKERRQ(ierr);
327f8170845SAlex Fikl     ierr = MatNestGetISs(A,is_row,is_col);CHKERRQ(ierr);
328cf37664fSBarry Smith     if (reuse == MAT_INPLACE_MATRIX) {
329ddeb9bd8SAlex Fikl       for (i=0; i<nr; i++) {
330ddeb9bd8SAlex Fikl         for (j=0; j<nc; j++) {
331ddeb9bd8SAlex Fikl           subs[i + nr * j] = bA->m[i][j];
332ddeb9bd8SAlex Fikl         }
333ddeb9bd8SAlex Fikl       }
334ddeb9bd8SAlex Fikl     }
335ddeb9bd8SAlex Fikl 
336f8170845SAlex Fikl     ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nc,is_col,nr,is_row,subs,&C);CHKERRQ(ierr);
337f8170845SAlex Fikl     ierr = PetscFree(subs);CHKERRQ(ierr);
3383d994f23SBarry Smith     ierr = PetscFree2(is_row,is_col);CHKERRQ(ierr);
339f8170845SAlex Fikl   } else {
340f8170845SAlex Fikl     C = *B;
341f8170845SAlex Fikl   }
342f8170845SAlex Fikl 
343f8170845SAlex Fikl   bC = (Mat_Nest*)C->data;
344f8170845SAlex Fikl   for (i=0; i<nr; i++) {
345f8170845SAlex Fikl     for (j=0; j<nc; j++) {
346f8170845SAlex Fikl       if (bA->m[i][j]) {
347f8170845SAlex Fikl         ierr = MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i]));CHKERRQ(ierr);
348f8170845SAlex Fikl       } else {
349f8170845SAlex Fikl         bC->m[j][i] = NULL;
350f8170845SAlex Fikl       }
351f8170845SAlex Fikl     }
352f8170845SAlex Fikl   }
353f8170845SAlex Fikl 
354cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
355f8170845SAlex Fikl     *B = C;
356f8170845SAlex Fikl   } else {
357f8170845SAlex Fikl     ierr = MatHeaderMerge(A, &C);CHKERRQ(ierr);
358f8170845SAlex Fikl   }
359f8170845SAlex Fikl   PetscFunctionReturn(0);
360f8170845SAlex Fikl }
361f8170845SAlex Fikl 
362e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
363e2d7f03fSJed Brown {
364e2d7f03fSJed Brown   PetscErrorCode ierr;
365e2d7f03fSJed Brown   IS             *lst = *list;
366e2d7f03fSJed Brown   PetscInt       i;
367e2d7f03fSJed Brown 
368e2d7f03fSJed Brown   PetscFunctionBegin;
369e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
3706bf464f9SBarry Smith   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
371e2d7f03fSJed Brown   ierr  = PetscFree(lst);CHKERRQ(ierr);
3720298fd71SBarry Smith   *list = NULL;
373e2d7f03fSJed Brown   PetscFunctionReturn(0);
374e2d7f03fSJed Brown }
375e2d7f03fSJed Brown 
37606a1af2fSStefano Zampini static PetscErrorCode MatReset_Nest(Mat A)
377d8588912SDave May {
378d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
379d8588912SDave May   PetscInt       i,j;
380d8588912SDave May   PetscErrorCode ierr;
381d8588912SDave May 
382d8588912SDave May   PetscFunctionBegin;
383d8588912SDave May   /* release the matrices and the place holders */
384e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
385e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
386e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
387e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
388d8588912SDave May 
389d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
390d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
39106a1af2fSStefano Zampini   ierr = PetscFree(vs->nnzstate);CHKERRQ(ierr);
392d8588912SDave May 
393207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
394207556f9SJed Brown 
395d8588912SDave May   /* release the matrices and the place holders */
396d8588912SDave May   if (vs->m) {
397d8588912SDave May     for (i=0; i<vs->nr; i++) {
398d8588912SDave May       for (j=0; j<vs->nc; j++) {
3996bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
400d8588912SDave May       }
401d8588912SDave May       ierr = PetscFree(vs->m[i]);CHKERRQ(ierr);
402d8588912SDave May     }
403d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
404d8588912SDave May   }
40506a1af2fSStefano Zampini 
40606a1af2fSStefano Zampini   /* restore defaults */
40706a1af2fSStefano Zampini   vs->nr = 0;
40806a1af2fSStefano Zampini   vs->nc = 0;
40906a1af2fSStefano Zampini   vs->splitassembly = PETSC_FALSE;
41006a1af2fSStefano Zampini   PetscFunctionReturn(0);
41106a1af2fSStefano Zampini }
41206a1af2fSStefano Zampini 
41306a1af2fSStefano Zampini static PetscErrorCode MatDestroy_Nest(Mat A)
41406a1af2fSStefano Zampini {
41506a1af2fSStefano Zampini   PetscErrorCode ierr;
41606a1af2fSStefano Zampini 
41706a1af2fSStefano Zampini   ierr = MatReset_Nest(A);CHKERRQ(ierr);
418bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
419d8588912SDave May 
420bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);CHKERRQ(ierr);
421bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);CHKERRQ(ierr);
422bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);CHKERRQ(ierr);
423bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);CHKERRQ(ierr);
424bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);CHKERRQ(ierr);
425bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);CHKERRQ(ierr);
426bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);CHKERRQ(ierr);
427bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);CHKERRQ(ierr);
4280899c546SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_mpiaij_C",0);CHKERRQ(ierr);
4290899c546SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_seqaij_C",0);CHKERRQ(ierr);
4305e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",0);CHKERRQ(ierr);
4315e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",0);CHKERRQ(ierr);
4324222ddf1SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_seqdense_C",NULL);CHKERRQ(ierr);
4334222ddf1SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_mpidense_C",NULL);CHKERRQ(ierr);
4344222ddf1SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_dense_C",NULL);CHKERRQ(ierr);
435d8588912SDave May   PetscFunctionReturn(0);
436d8588912SDave May }
437d8588912SDave May 
438381b8e50SStefano Zampini static PetscErrorCode MatMissingDiagonal_Nest(Mat mat,PetscBool *missing,PetscInt *dd)
439381b8e50SStefano Zampini {
440381b8e50SStefano Zampini   Mat_Nest       *vs = (Mat_Nest*)mat->data;
441381b8e50SStefano Zampini   PetscInt       i;
442381b8e50SStefano Zampini   PetscErrorCode ierr;
443381b8e50SStefano Zampini 
444381b8e50SStefano Zampini   PetscFunctionBegin;
445381b8e50SStefano Zampini   if (dd) *dd = 0;
446381b8e50SStefano Zampini   if (!vs->nr) {
447381b8e50SStefano Zampini     *missing = PETSC_TRUE;
448381b8e50SStefano Zampini     PetscFunctionReturn(0);
449381b8e50SStefano Zampini   }
450381b8e50SStefano Zampini   *missing = PETSC_FALSE;
451381b8e50SStefano Zampini   for (i = 0; i < vs->nr && !(*missing); i++) {
452381b8e50SStefano Zampini     *missing = PETSC_TRUE;
453381b8e50SStefano Zampini     if (vs->m[i][i]) {
454381b8e50SStefano Zampini       ierr = MatMissingDiagonal(vs->m[i][i],missing,NULL);CHKERRQ(ierr);
455381b8e50SStefano Zampini       if (*missing && dd) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"First missing entry not yet implemented");
456381b8e50SStefano Zampini     }
457381b8e50SStefano Zampini   }
458381b8e50SStefano Zampini   PetscFunctionReturn(0);
459381b8e50SStefano Zampini }
460381b8e50SStefano Zampini 
461207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
462d8588912SDave May {
463d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
464d8588912SDave May   PetscInt       i,j;
465d8588912SDave May   PetscErrorCode ierr;
46606a1af2fSStefano Zampini   PetscBool      nnzstate = PETSC_FALSE;
467d8588912SDave May 
468d8588912SDave May   PetscFunctionBegin;
469d8588912SDave May   for (i=0; i<vs->nr; i++) {
470d8588912SDave May     for (j=0; j<vs->nc; j++) {
47106a1af2fSStefano Zampini       PetscObjectState subnnzstate = 0;
472e7c19651SJed Brown       if (vs->m[i][j]) {
473e7c19651SJed Brown         ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr);
474e7c19651SJed Brown         if (!vs->splitassembly) {
475e7c19651SJed Brown           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
476e7c19651SJed 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
477e7c19651SJed Brown            * already performing an assembly, but the result would by more complicated and appears to offer less
478e7c19651SJed Brown            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
479e7c19651SJed Brown            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
480e7c19651SJed Brown            */
481e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
48206a1af2fSStefano Zampini           ierr = MatGetNonzeroState(vs->m[i][j],&subnnzstate);CHKERRQ(ierr);
483e7c19651SJed Brown         }
484e7c19651SJed Brown       }
48506a1af2fSStefano Zampini       nnzstate = (PetscBool)(nnzstate || vs->nnzstate[i*vs->nc+j] != subnnzstate);
48606a1af2fSStefano Zampini       vs->nnzstate[i*vs->nc+j] = subnnzstate;
487d8588912SDave May     }
488d8588912SDave May   }
48906a1af2fSStefano Zampini   if (nnzstate) A->nonzerostate++;
490d8588912SDave May   PetscFunctionReturn(0);
491d8588912SDave May }
492d8588912SDave May 
493207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
494d8588912SDave May {
495d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
496d8588912SDave May   PetscInt       i,j;
497d8588912SDave May   PetscErrorCode ierr;
498d8588912SDave May 
499d8588912SDave May   PetscFunctionBegin;
500d8588912SDave May   for (i=0; i<vs->nr; i++) {
501d8588912SDave May     for (j=0; j<vs->nc; j++) {
502e7c19651SJed Brown       if (vs->m[i][j]) {
503e7c19651SJed Brown         if (vs->splitassembly) {
504e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
505e7c19651SJed Brown         }
506e7c19651SJed Brown       }
507d8588912SDave May     }
508d8588912SDave May   }
509d8588912SDave May   PetscFunctionReturn(0);
510d8588912SDave May }
511d8588912SDave May 
512f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
513d8588912SDave May {
514207556f9SJed Brown   PetscErrorCode ierr;
515f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
516f349c1fdSJed Brown   PetscInt       j;
517f349c1fdSJed Brown   Mat            sub;
518d8588912SDave May 
519d8588912SDave May   PetscFunctionBegin;
5200298fd71SBarry Smith   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
521f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
5224994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
523f349c1fdSJed Brown   *B = sub;
524f349c1fdSJed Brown   PetscFunctionReturn(0);
525d8588912SDave May }
526d8588912SDave May 
527f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
528f349c1fdSJed Brown {
529207556f9SJed Brown   PetscErrorCode ierr;
530f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
531f349c1fdSJed Brown   PetscInt       i;
532f349c1fdSJed Brown   Mat            sub;
533f349c1fdSJed Brown 
534f349c1fdSJed Brown   PetscFunctionBegin;
5350298fd71SBarry Smith   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
536f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
5374994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
538f349c1fdSJed Brown   *B = sub;
539f349c1fdSJed Brown   PetscFunctionReturn(0);
540d8588912SDave May }
541d8588912SDave May 
542f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
543f349c1fdSJed Brown {
544f349c1fdSJed Brown   PetscErrorCode ierr;
545f349c1fdSJed Brown   PetscInt       i;
546f349c1fdSJed Brown   PetscBool      flg;
547f349c1fdSJed Brown 
548f349c1fdSJed Brown   PetscFunctionBegin;
549f349c1fdSJed Brown   PetscValidPointer(list,3);
550f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
551f349c1fdSJed Brown   PetscValidIntPointer(found,5);
552f349c1fdSJed Brown   *found = -1;
553f349c1fdSJed Brown   for (i=0; i<n; i++) {
554207556f9SJed Brown     if (!list[i]) continue;
555320466b0SStefano Zampini     ierr = ISEqualUnsorted(list[i],is,&flg);CHKERRQ(ierr);
556f349c1fdSJed Brown     if (flg) {
557f349c1fdSJed Brown       *found = i;
558f349c1fdSJed Brown       PetscFunctionReturn(0);
559f349c1fdSJed Brown     }
560f349c1fdSJed Brown   }
561ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set");
562f349c1fdSJed Brown   PetscFunctionReturn(0);
563f349c1fdSJed Brown }
564f349c1fdSJed Brown 
5658188e55aSJed Brown /* Get a block row as a new MatNest */
5668188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
5678188e55aSJed Brown {
5688188e55aSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
5698188e55aSJed Brown   char           keyname[256];
5708188e55aSJed Brown   PetscErrorCode ierr;
5718188e55aSJed Brown 
5728188e55aSJed Brown   PetscFunctionBegin;
5730298fd71SBarry Smith   *B   = NULL;
5748caf3d72SBarry Smith   ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr);
5758188e55aSJed Brown   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
5768188e55aSJed Brown   if (*B) PetscFunctionReturn(0);
5778188e55aSJed Brown 
578ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
57926fbe8dcSKarl Rupp 
5808188e55aSJed Brown   (*B)->assembled = A->assembled;
58126fbe8dcSKarl Rupp 
5828188e55aSJed Brown   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
5838188e55aSJed Brown   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
5848188e55aSJed Brown   PetscFunctionReturn(0);
5858188e55aSJed Brown }
5868188e55aSJed Brown 
587f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
588f349c1fdSJed Brown {
589f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
5908188e55aSJed Brown   PetscErrorCode ierr;
5916b3a5b13SJed Brown   PetscInt       row,col;
592e072481dSJed Brown   PetscBool      same,isFullCol,isFullColGlobal;
593f349c1fdSJed Brown 
594f349c1fdSJed Brown   PetscFunctionBegin;
5958188e55aSJed Brown   /* Check if full column space. This is a hack */
5968188e55aSJed Brown   isFullCol = PETSC_FALSE;
597251f4c67SDmitry Karpeev   ierr      = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
5988188e55aSJed Brown   if (same) {
59977019fcaSJed Brown     PetscInt n,first,step,i,an,am,afirst,astep;
6008188e55aSJed Brown     ierr      = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
6018188e55aSJed Brown     ierr      = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
60277019fcaSJed Brown     isFullCol = PETSC_TRUE;
60305ce4453SJed Brown     for (i=0,an=A->cmap->rstart; i<vs->nc; i++) {
60406a1af2fSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)is->col[i],ISSTRIDE,&same);CHKERRQ(ierr);
60577019fcaSJed Brown       ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr);
60606a1af2fSStefano Zampini       if (same) {
60706a1af2fSStefano Zampini         ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr);
60877019fcaSJed Brown         if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
60906a1af2fSStefano Zampini       } else isFullCol = PETSC_FALSE;
61077019fcaSJed Brown       an += am;
61177019fcaSJed Brown     }
61205ce4453SJed Brown     if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE;
6138188e55aSJed Brown   }
614b2566f29SBarry Smith   ierr = MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr);
6158188e55aSJed Brown 
616427230ceSLisandro Dalcin   if (isFullColGlobal && vs->nc > 1) {
6178188e55aSJed Brown     PetscInt row;
6188188e55aSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
6198188e55aSJed Brown     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
6208188e55aSJed Brown   } else {
621f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
622f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
623b6480e04SStefano Zampini     if (!vs->m[row][col]) {
624b6480e04SStefano Zampini       PetscInt lr,lc;
625b6480e04SStefano Zampini 
626b6480e04SStefano Zampini       ierr = MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);CHKERRQ(ierr);
627b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.row[row],&lr);CHKERRQ(ierr);
628b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.col[col],&lc);CHKERRQ(ierr);
629b6480e04SStefano Zampini       ierr = MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
630fa9f0909SStefano Zampini       ierr = MatSetType(vs->m[row][col],MATAIJ);CHKERRQ(ierr);
631fa9f0909SStefano Zampini       ierr = MatSeqAIJSetPreallocation(vs->m[row][col],0,NULL);CHKERRQ(ierr);
632fa9f0909SStefano Zampini       ierr = MatMPIAIJSetPreallocation(vs->m[row][col],0,NULL,0,NULL);CHKERRQ(ierr);
633b6480e04SStefano Zampini       ierr = MatSetUp(vs->m[row][col]);CHKERRQ(ierr);
634b6480e04SStefano Zampini       ierr = MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
635b6480e04SStefano Zampini       ierr = MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
636b6480e04SStefano Zampini     }
637f349c1fdSJed Brown     *B = vs->m[row][col];
6388188e55aSJed Brown   }
639f349c1fdSJed Brown   PetscFunctionReturn(0);
640f349c1fdSJed Brown }
641f349c1fdSJed Brown 
64206a1af2fSStefano Zampini /*
64306a1af2fSStefano Zampini    TODO: This does not actually returns a submatrix we can modify
64406a1af2fSStefano Zampini */
6457dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
646f349c1fdSJed Brown {
647f349c1fdSJed Brown   PetscErrorCode ierr;
648f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
649f349c1fdSJed Brown   Mat            sub;
650f349c1fdSJed Brown 
651f349c1fdSJed Brown   PetscFunctionBegin;
652f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
653f349c1fdSJed Brown   switch (reuse) {
654f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
6557874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
656f349c1fdSJed Brown     *B = sub;
657f349c1fdSJed Brown     break;
658f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
659ce94432eSBarry Smith     if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
660f349c1fdSJed Brown     break;
661f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
662f349c1fdSJed Brown     break;
663511c6705SHong Zhang   case MAT_INPLACE_MATRIX:       /* Nothing to do */
664511c6705SHong Zhang     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet");
665511c6705SHong Zhang     break;
666f349c1fdSJed Brown   }
667f349c1fdSJed Brown   PetscFunctionReturn(0);
668f349c1fdSJed Brown }
669f349c1fdSJed Brown 
670f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
671f349c1fdSJed Brown {
672f349c1fdSJed Brown   PetscErrorCode ierr;
673f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
674f349c1fdSJed Brown   Mat            sub;
675f349c1fdSJed Brown 
676f349c1fdSJed Brown   PetscFunctionBegin;
677f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
678f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
679f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
680f349c1fdSJed Brown   *B = sub;
681d8588912SDave May   PetscFunctionReturn(0);
682d8588912SDave May }
683d8588912SDave May 
684207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
685d8588912SDave May {
686d8588912SDave May   PetscErrorCode ierr;
687f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
688f349c1fdSJed Brown   Mat            sub;
689d8588912SDave May 
690d8588912SDave May   PetscFunctionBegin;
691f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
692ce94432eSBarry Smith   if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
693f349c1fdSJed Brown   if (sub) {
694ce94432eSBarry Smith     if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
6956bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
696d8588912SDave May   }
697d8588912SDave May   PetscFunctionReturn(0);
698d8588912SDave May }
699d8588912SDave May 
7007874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
7017874fa86SDave May {
7027874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
7037874fa86SDave May   PetscInt       i;
7047874fa86SDave May   PetscErrorCode ierr;
7057874fa86SDave May 
7067874fa86SDave May   PetscFunctionBegin;
7077874fa86SDave May   for (i=0; i<bA->nr; i++) {
708429bac76SJed Brown     Vec bv;
709429bac76SJed Brown     ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
7107874fa86SDave May     if (bA->m[i][i]) {
711429bac76SJed Brown       ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr);
7127874fa86SDave May     } else {
7135159a857SMatthew G. Knepley       ierr = VecSet(bv,0.0);CHKERRQ(ierr);
7147874fa86SDave May     }
715429bac76SJed Brown     ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
7167874fa86SDave May   }
7177874fa86SDave May   PetscFunctionReturn(0);
7187874fa86SDave May }
7197874fa86SDave May 
7207874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
7217874fa86SDave May {
7227874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
723429bac76SJed Brown   Vec            bl,*br;
7247874fa86SDave May   PetscInt       i,j;
7257874fa86SDave May   PetscErrorCode ierr;
7267874fa86SDave May 
7277874fa86SDave May   PetscFunctionBegin;
7283f800ebeSJed Brown   ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr);
7292e6472ebSElliott Sales de Andrade   if (r) {
730429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
7312e6472ebSElliott Sales de Andrade   }
7322e6472ebSElliott Sales de Andrade   bl = NULL;
7337874fa86SDave May   for (i=0; i<bA->nr; i++) {
7342e6472ebSElliott Sales de Andrade     if (l) {
735429bac76SJed Brown       ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
7362e6472ebSElliott Sales de Andrade     }
7377874fa86SDave May     for (j=0; j<bA->nc; j++) {
7387874fa86SDave May       if (bA->m[i][j]) {
739429bac76SJed Brown         ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr);
7407874fa86SDave May       }
7417874fa86SDave May     }
7422e6472ebSElliott Sales de Andrade     if (l) {
743a061e289SJed Brown       ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
7447874fa86SDave May     }
7452e6472ebSElliott Sales de Andrade   }
7462e6472ebSElliott Sales de Andrade   if (r) {
747429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
7482e6472ebSElliott Sales de Andrade   }
749429bac76SJed Brown   ierr = PetscFree(br);CHKERRQ(ierr);
7507874fa86SDave May   PetscFunctionReturn(0);
7517874fa86SDave May }
7527874fa86SDave May 
753a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
754a061e289SJed Brown {
755a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
756a061e289SJed Brown   PetscInt       i,j;
757a061e289SJed Brown   PetscErrorCode ierr;
758a061e289SJed Brown 
759a061e289SJed Brown   PetscFunctionBegin;
760a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
761a061e289SJed Brown     for (j=0; j<bA->nc; j++) {
762a061e289SJed Brown       if (bA->m[i][j]) {
763a061e289SJed Brown         ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr);
764a061e289SJed Brown       }
765a061e289SJed Brown     }
766a061e289SJed Brown   }
767a061e289SJed Brown   PetscFunctionReturn(0);
768a061e289SJed Brown }
769a061e289SJed Brown 
770a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
771a061e289SJed Brown {
772a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
773a061e289SJed Brown   PetscInt       i;
774a061e289SJed Brown   PetscErrorCode ierr;
77506a1af2fSStefano Zampini   PetscBool      nnzstate = PETSC_FALSE;
776a061e289SJed Brown 
777a061e289SJed Brown   PetscFunctionBegin;
778a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
77906a1af2fSStefano Zampini     PetscObjectState subnnzstate = 0;
780ce94432eSBarry Smith     if (!bA->m[i][i]) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for shifting an empty diagonal block, insert a matrix in block (%D,%D)",i,i);
781a061e289SJed Brown     ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr);
78206a1af2fSStefano Zampini     ierr = MatGetNonzeroState(bA->m[i][i],&subnnzstate);CHKERRQ(ierr);
78306a1af2fSStefano Zampini     nnzstate = (PetscBool)(nnzstate || bA->nnzstate[i*bA->nc+i] != subnnzstate);
78406a1af2fSStefano Zampini     bA->nnzstate[i*bA->nc+i] = subnnzstate;
785a061e289SJed Brown   }
78606a1af2fSStefano Zampini   if (nnzstate) A->nonzerostate++;
787a061e289SJed Brown   PetscFunctionReturn(0);
788a061e289SJed Brown }
789a061e289SJed Brown 
79013135bc6SAlex Fikl static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is)
79113135bc6SAlex Fikl {
79213135bc6SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
79313135bc6SAlex Fikl   PetscInt       i;
79413135bc6SAlex Fikl   PetscErrorCode ierr;
79506a1af2fSStefano Zampini   PetscBool      nnzstate = PETSC_FALSE;
79613135bc6SAlex Fikl 
79713135bc6SAlex Fikl   PetscFunctionBegin;
79813135bc6SAlex Fikl   for (i=0; i<bA->nr; i++) {
79906a1af2fSStefano Zampini     PetscObjectState subnnzstate = 0;
80013135bc6SAlex Fikl     Vec              bv;
80113135bc6SAlex Fikl     ierr = VecGetSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
80213135bc6SAlex Fikl     if (bA->m[i][i]) {
80313135bc6SAlex Fikl       ierr = MatDiagonalSet(bA->m[i][i],bv,is);CHKERRQ(ierr);
80406a1af2fSStefano Zampini       ierr = MatGetNonzeroState(bA->m[i][i],&subnnzstate);CHKERRQ(ierr);
80513135bc6SAlex Fikl     }
80613135bc6SAlex Fikl     ierr = VecRestoreSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
80706a1af2fSStefano Zampini     nnzstate = (PetscBool)(nnzstate || bA->nnzstate[i*bA->nc+i] != subnnzstate);
80806a1af2fSStefano Zampini     bA->nnzstate[i*bA->nc+i] = subnnzstate;
80913135bc6SAlex Fikl   }
81006a1af2fSStefano Zampini   if (nnzstate) A->nonzerostate++;
81113135bc6SAlex Fikl   PetscFunctionReturn(0);
81213135bc6SAlex Fikl }
81313135bc6SAlex Fikl 
814f8170845SAlex Fikl static PetscErrorCode MatSetRandom_Nest(Mat A,PetscRandom rctx)
815f8170845SAlex Fikl {
816f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
817f8170845SAlex Fikl   PetscInt       i,j;
818f8170845SAlex Fikl   PetscErrorCode ierr;
819f8170845SAlex Fikl 
820f8170845SAlex Fikl   PetscFunctionBegin;
821f8170845SAlex Fikl   for (i=0; i<bA->nr; i++) {
822f8170845SAlex Fikl     for (j=0; j<bA->nc; j++) {
823f8170845SAlex Fikl       if (bA->m[i][j]) {
824f8170845SAlex Fikl         ierr = MatSetRandom(bA->m[i][j],rctx);CHKERRQ(ierr);
825f8170845SAlex Fikl       }
826f8170845SAlex Fikl     }
827f8170845SAlex Fikl   }
828f8170845SAlex Fikl   PetscFunctionReturn(0);
829f8170845SAlex Fikl }
830f8170845SAlex Fikl 
8312a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left)
832d8588912SDave May {
833d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
834d8588912SDave May   Vec            *L,*R;
835d8588912SDave May   MPI_Comm       comm;
836d8588912SDave May   PetscInt       i,j;
837d8588912SDave May   PetscErrorCode ierr;
838d8588912SDave May 
839d8588912SDave May   PetscFunctionBegin;
840ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
841d8588912SDave May   if (right) {
842d8588912SDave May     /* allocate R */
843854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nc, &R);CHKERRQ(ierr);
844d8588912SDave May     /* Create the right vectors */
845d8588912SDave May     for (j=0; j<bA->nc; j++) {
846d8588912SDave May       for (i=0; i<bA->nr; i++) {
847d8588912SDave May         if (bA->m[i][j]) {
8482a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr);
849d8588912SDave May           break;
850d8588912SDave May         }
851d8588912SDave May       }
8526c4ed002SBarry Smith       if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
853d8588912SDave May     }
854f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
855d8588912SDave May     /* hand back control to the nest vector */
856d8588912SDave May     for (j=0; j<bA->nc; j++) {
8576bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
858d8588912SDave May     }
859d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
860d8588912SDave May   }
861d8588912SDave May 
862d8588912SDave May   if (left) {
863d8588912SDave May     /* allocate L */
864854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nr, &L);CHKERRQ(ierr);
865d8588912SDave May     /* Create the left vectors */
866d8588912SDave May     for (i=0; i<bA->nr; i++) {
867d8588912SDave May       for (j=0; j<bA->nc; j++) {
868d8588912SDave May         if (bA->m[i][j]) {
8692a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr);
870d8588912SDave May           break;
871d8588912SDave May         }
872d8588912SDave May       }
8736c4ed002SBarry Smith       if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
874d8588912SDave May     }
875d8588912SDave May 
876f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
877d8588912SDave May     for (i=0; i<bA->nr; i++) {
8786bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
879d8588912SDave May     }
880d8588912SDave May 
881d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
882d8588912SDave May   }
883d8588912SDave May   PetscFunctionReturn(0);
884d8588912SDave May }
885d8588912SDave May 
886207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
887d8588912SDave May {
888d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
88929e60adbSStefano Zampini   PetscBool      isascii,viewSub = PETSC_FALSE;
890d8588912SDave May   PetscInt       i,j;
891d8588912SDave May   PetscErrorCode ierr;
892d8588912SDave May 
893d8588912SDave May   PetscFunctionBegin;
894251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
895d8588912SDave May   if (isascii) {
896d8588912SDave May 
89729e60adbSStefano Zampini     ierr = PetscOptionsGetBool(((PetscObject)A)->options,((PetscObject)A)->prefix,"-mat_view_nest_sub",&viewSub,NULL);CHKERRQ(ierr);
898d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Matrix object: \n");CHKERRQ(ierr);
899d86155a6SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
900d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "type=nest, rows=%D, cols=%D \n",bA->nr,bA->nc);CHKERRQ(ierr);
901d8588912SDave May 
902d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");CHKERRQ(ierr);
903d8588912SDave May     for (i=0; i<bA->nr; i++) {
904d8588912SDave May       for (j=0; j<bA->nc; j++) {
90519fd82e9SBarry Smith         MatType   type;
906270f95d7SJed Brown         char      name[256] = "",prefix[256] = "";
907d8588912SDave May         PetscInt  NR,NC;
908d8588912SDave May         PetscBool isNest = PETSC_FALSE;
909d8588912SDave May 
910d8588912SDave May         if (!bA->m[i][j]) {
911d86155a6SBarry Smith           CHKERRQ(ierr);PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);CHKERRQ(ierr);
912d8588912SDave May           continue;
913d8588912SDave May         }
914d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
915d8588912SDave May         ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr);
9168caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);}
9178caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);}
918251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
919d8588912SDave May 
920270f95d7SJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr);
921d8588912SDave May 
92229e60adbSStefano Zampini         if (isNest || viewSub) {
923270f95d7SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* push1 */
924d8588912SDave May           ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr);
925270f95d7SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop1 */
926d8588912SDave May         }
927d8588912SDave May       }
928d8588912SDave May     }
929d86155a6SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop0 */
930d8588912SDave May   }
931d8588912SDave May   PetscFunctionReturn(0);
932d8588912SDave May }
933d8588912SDave May 
934207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
935d8588912SDave May {
936d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
937d8588912SDave May   PetscInt       i,j;
938d8588912SDave May   PetscErrorCode ierr;
939d8588912SDave May 
940d8588912SDave May   PetscFunctionBegin;
941d8588912SDave May   for (i=0; i<bA->nr; i++) {
942d8588912SDave May     for (j=0; j<bA->nc; j++) {
943d8588912SDave May       if (!bA->m[i][j]) continue;
944d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
945d8588912SDave May     }
946d8588912SDave May   }
947d8588912SDave May   PetscFunctionReturn(0);
948d8588912SDave May }
949d8588912SDave May 
950c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str)
951c222c20dSDavid Ham {
952c222c20dSDavid Ham   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
953c222c20dSDavid Ham   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
954c222c20dSDavid Ham   PetscErrorCode ierr;
95506a1af2fSStefano Zampini   PetscBool      nnzstate = PETSC_FALSE;
956c222c20dSDavid Ham 
957c222c20dSDavid Ham   PetscFunctionBegin;
958c222c20dSDavid Ham   if (nr != bB->nr || nc != bB->nc) SETERRQ4(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Cannot copy a Mat_Nest of block size (%D,%D) to a Mat_Nest of block size (%D,%D)",bB->nr,bB->nc,nr,nc);
959c222c20dSDavid Ham   for (i=0; i<nr; i++) {
960c222c20dSDavid Ham     for (j=0; j<nc; j++) {
96106a1af2fSStefano Zampini       PetscObjectState subnnzstate = 0;
96246a2b97cSJed Brown       if (bA->m[i][j] && bB->m[i][j]) {
963c222c20dSDavid Ham         ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr);
96446a2b97cSJed Brown       } else if (bA->m[i][j] || bB->m[i][j]) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Matrix block does not exist at %D,%D",i,j);
96506a1af2fSStefano Zampini       ierr = MatGetNonzeroState(bB->m[i][j],&subnnzstate);CHKERRQ(ierr);
96606a1af2fSStefano Zampini       nnzstate = (PetscBool)(nnzstate || bB->nnzstate[i*nc+j] != subnnzstate);
96706a1af2fSStefano Zampini       bB->nnzstate[i*nc+j] = subnnzstate;
968c222c20dSDavid Ham     }
969c222c20dSDavid Ham   }
97006a1af2fSStefano Zampini   if (nnzstate) B->nonzerostate++;
971c222c20dSDavid Ham   PetscFunctionReturn(0);
972c222c20dSDavid Ham }
973c222c20dSDavid Ham 
9746e76ffeaSPierre Jolivet static PetscErrorCode MatAXPY_Nest(Mat Y,PetscScalar a,Mat X,MatStructure str)
9756e76ffeaSPierre Jolivet {
9766e76ffeaSPierre Jolivet   Mat_Nest       *bY = (Mat_Nest*)Y->data,*bX = (Mat_Nest*)X->data;
9776e76ffeaSPierre Jolivet   PetscInt       i,j,nr = bY->nr,nc = bY->nc;
9786e76ffeaSPierre Jolivet   PetscErrorCode ierr;
97906a1af2fSStefano Zampini   PetscBool      nnzstate = PETSC_FALSE;
9806e76ffeaSPierre Jolivet 
9816e76ffeaSPierre Jolivet   PetscFunctionBegin;
9826e76ffeaSPierre Jolivet   if (nr != bX->nr || nc != bX->nc) SETERRQ4(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Cannot AXPY a MatNest of block size (%D,%D) with a MatNest of block size (%D,%D)",bX->nr,bX->nc,nr,nc);
9836e76ffeaSPierre Jolivet   for (i=0; i<nr; i++) {
9846e76ffeaSPierre Jolivet     for (j=0; j<nc; j++) {
98506a1af2fSStefano Zampini       PetscObjectState subnnzstate = 0;
9866e76ffeaSPierre Jolivet       if (bY->m[i][j] && bX->m[i][j]) {
9876e76ffeaSPierre Jolivet         ierr = MatAXPY(bY->m[i][j],a,bX->m[i][j],str);CHKERRQ(ierr);
988c066aebcSStefano Zampini       } else if (bX->m[i][j]) {
989c066aebcSStefano Zampini         Mat M;
990c066aebcSStefano Zampini 
991060bfc19SStefano Zampini         if (str != DIFFERENT_NONZERO_PATTERN) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Matrix block does not exist at %D,%D. Use DIFFERENT_NONZERO_PATTERN",i,j);
992c066aebcSStefano Zampini         ierr = MatDuplicate(bX->m[i][j],MAT_COPY_VALUES,&M);CHKERRQ(ierr);
993c066aebcSStefano Zampini         ierr = MatNestSetSubMat(Y,i,j,M);CHKERRQ(ierr);
994c066aebcSStefano Zampini         ierr = MatDestroy(&M);CHKERRQ(ierr);
995c066aebcSStefano Zampini       }
996060bfc19SStefano Zampini       if (bY->m[i][j]) { ierr = MatGetNonzeroState(bY->m[i][j],&subnnzstate);CHKERRQ(ierr); }
99706a1af2fSStefano Zampini       nnzstate = (PetscBool)(nnzstate || bY->nnzstate[i*nc+j] != subnnzstate);
99806a1af2fSStefano Zampini       bY->nnzstate[i*nc+j] = subnnzstate;
9996e76ffeaSPierre Jolivet     }
10006e76ffeaSPierre Jolivet   }
100106a1af2fSStefano Zampini   if (nnzstate) Y->nonzerostate++;
10026e76ffeaSPierre Jolivet   PetscFunctionReturn(0);
10036e76ffeaSPierre Jolivet }
10046e76ffeaSPierre Jolivet 
1005207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
1006d8588912SDave May {
1007d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
1008841e96a3SJed Brown   Mat            *b;
1009841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
1010d8588912SDave May   PetscErrorCode ierr;
1011d8588912SDave May 
1012d8588912SDave May   PetscFunctionBegin;
1013785e854fSJed Brown   ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr);
1014841e96a3SJed Brown   for (i=0; i<nr; i++) {
1015841e96a3SJed Brown     for (j=0; j<nc; j++) {
1016841e96a3SJed Brown       if (bA->m[i][j]) {
1017841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
1018841e96a3SJed Brown       } else {
10190298fd71SBarry Smith         b[i*nc+j] = NULL;
1020d8588912SDave May       }
1021d8588912SDave May     }
1022d8588912SDave May   }
1023ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
1024841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
1025841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
10266bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
1027d8588912SDave May   }
1028d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
1029d8588912SDave May 
1030841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1031841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1032d8588912SDave May   PetscFunctionReturn(0);
1033d8588912SDave May }
1034d8588912SDave May 
1035d8588912SDave May /* nest api */
1036d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
1037d8588912SDave May {
1038d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
10395fd66863SKarl Rupp 
1040d8588912SDave May   PetscFunctionBegin;
1041ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
1042ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
1043d8588912SDave May   *mat = bA->m[idxm][jdxm];
1044d8588912SDave May   PetscFunctionReturn(0);
1045d8588912SDave May }
1046d8588912SDave May 
10479ba0d327SJed Brown /*@
1048d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
1049d8588912SDave May 
1050d8588912SDave May  Not collective
1051d8588912SDave May 
1052d8588912SDave May  Input Parameters:
1053629881c0SJed Brown +   A  - nest matrix
1054d8588912SDave May .   idxm - index of the matrix within the nest matrix
1055629881c0SJed Brown -   jdxm - index of the matrix within the nest matrix
1056d8588912SDave May 
1057d8588912SDave May  Output Parameter:
1058d8588912SDave May .   sub - matrix at index idxm,jdxm within the nest matrix
1059d8588912SDave May 
1060d8588912SDave May  Level: developer
1061d8588912SDave May 
1062*bb97c47cSPierre Jolivet .seealso: MatNestGetSize(), MatNestGetSubMats(), MatCreateNest(), MATNEST, MatNestSetSubMat(),
106379798668SBarry Smith           MatNestGetLocalISs(), MatNestGetISs()
1064d8588912SDave May @*/
10657087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
1066d8588912SDave May {
1067699a902aSJed Brown   PetscErrorCode ierr;
1068d8588912SDave May 
1069d8588912SDave May   PetscFunctionBegin;
1070699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
1071d8588912SDave May   PetscFunctionReturn(0);
1072d8588912SDave May }
1073d8588912SDave May 
10740782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
10750782ca92SJed Brown {
10760782ca92SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
10770782ca92SJed Brown   PetscInt       m,n,M,N,mi,ni,Mi,Ni;
10780782ca92SJed Brown   PetscErrorCode ierr;
10790782ca92SJed Brown 
10800782ca92SJed Brown   PetscFunctionBegin;
1081ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
1082ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
10830782ca92SJed Brown   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
10840782ca92SJed Brown   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
10850782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
10860782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
10870782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
10880782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
1089ce94432eSBarry Smith   if (M != Mi || N != Ni) SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"Submatrix dimension (%D,%D) incompatible with nest block (%D,%D)",M,N,Mi,Ni);
1090ce94432eSBarry Smith   if (m != mi || n != ni) SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"Submatrix local dimension (%D,%D) incompatible with nest block (%D,%D)",m,n,mi,ni);
109126fbe8dcSKarl Rupp 
109206a1af2fSStefano Zampini   /* do not increase object state */
109306a1af2fSStefano Zampini   if (mat == bA->m[idxm][jdxm]) PetscFunctionReturn(0);
109406a1af2fSStefano Zampini 
10950782ca92SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
10960782ca92SJed Brown   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
10970782ca92SJed Brown   bA->m[idxm][jdxm] = mat;
109806a1af2fSStefano Zampini   ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr);
109906a1af2fSStefano Zampini   ierr = MatGetNonzeroState(mat,&bA->nnzstate[idxm*bA->nc+jdxm]);CHKERRQ(ierr);
110006a1af2fSStefano Zampini   A->nonzerostate++;
11010782ca92SJed Brown   PetscFunctionReturn(0);
11020782ca92SJed Brown }
11030782ca92SJed Brown 
11049ba0d327SJed Brown /*@
11050782ca92SJed Brown  MatNestSetSubMat - Set a single submatrix in the nest matrix.
11060782ca92SJed Brown 
11070782ca92SJed Brown  Logically collective on the submatrix communicator
11080782ca92SJed Brown 
11090782ca92SJed Brown  Input Parameters:
11100782ca92SJed Brown +   A  - nest matrix
11110782ca92SJed Brown .   idxm - index of the matrix within the nest matrix
11120782ca92SJed Brown .   jdxm - index of the matrix within the nest matrix
11130782ca92SJed Brown -   sub - matrix at index idxm,jdxm within the nest matrix
11140782ca92SJed Brown 
11150782ca92SJed Brown  Notes:
11160782ca92SJed Brown  The new submatrix must have the same size and communicator as that block of the nest.
11170782ca92SJed Brown 
11180782ca92SJed Brown  This increments the reference count of the submatrix.
11190782ca92SJed Brown 
11200782ca92SJed Brown  Level: developer
11210782ca92SJed Brown 
1122*bb97c47cSPierre Jolivet .seealso: MatNestSetSubMats(), MatNestGetSubMats(), MatNestGetLocalISs(), MATNEST, MatCreateNest(),
112379798668SBarry Smith           MatNestGetSubMat(), MatNestGetISs(), MatNestGetSize()
11240782ca92SJed Brown @*/
11250782ca92SJed Brown PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
11260782ca92SJed Brown {
11270782ca92SJed Brown   PetscErrorCode ierr;
11280782ca92SJed Brown 
11290782ca92SJed Brown   PetscFunctionBegin;
11300782ca92SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
11310782ca92SJed Brown   PetscFunctionReturn(0);
11320782ca92SJed Brown }
11330782ca92SJed Brown 
1134d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
1135d8588912SDave May {
1136d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
11375fd66863SKarl Rupp 
1138d8588912SDave May   PetscFunctionBegin;
113926fbe8dcSKarl Rupp   if (M)   *M   = bA->nr;
114026fbe8dcSKarl Rupp   if (N)   *N   = bA->nc;
114126fbe8dcSKarl Rupp   if (mat) *mat = bA->m;
1142d8588912SDave May   PetscFunctionReturn(0);
1143d8588912SDave May }
1144d8588912SDave May 
1145d8588912SDave May /*@C
1146d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
1147d8588912SDave May 
1148d8588912SDave May  Not collective
1149d8588912SDave May 
1150d8588912SDave May  Input Parameters:
1151629881c0SJed Brown .   A  - nest matrix
1152d8588912SDave May 
1153d8588912SDave May  Output Parameter:
1154629881c0SJed Brown +   M - number of rows in the nest matrix
1155d8588912SDave May .   N - number of cols in the nest matrix
1156629881c0SJed Brown -   mat - 2d array of matrices
1157d8588912SDave May 
1158d8588912SDave May  Notes:
1159d8588912SDave May 
1160d8588912SDave May  The user should not free the array mat.
1161d8588912SDave May 
1162351962e3SVincent Le Chenadec  In Fortran, this routine has a calling sequence
1163351962e3SVincent Le Chenadec $   call MatNestGetSubMats(A, M, N, mat, ierr)
1164351962e3SVincent Le Chenadec  where the space allocated for the optional argument mat is assumed large enough (if provided).
1165351962e3SVincent Le Chenadec 
1166d8588912SDave May  Level: developer
1167d8588912SDave May 
1168*bb97c47cSPierre Jolivet .seealso: MatNestGetSize(), MatNestGetSubMat(), MatNestGetLocalISs(), MATNEST, MatCreateNest(),
116979798668SBarry Smith           MatNestSetSubMats(), MatNestGetISs(), MatNestSetSubMat()
1170d8588912SDave May @*/
11717087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
1172d8588912SDave May {
1173699a902aSJed Brown   PetscErrorCode ierr;
1174d8588912SDave May 
1175d8588912SDave May   PetscFunctionBegin;
1176699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
1177d8588912SDave May   PetscFunctionReturn(0);
1178d8588912SDave May }
1179d8588912SDave May 
11807087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
1181d8588912SDave May {
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;
1187d8588912SDave May   PetscFunctionReturn(0);
1188d8588912SDave May }
1189d8588912SDave May 
11909ba0d327SJed Brown /*@
1191d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
1192d8588912SDave May 
1193d8588912SDave May  Not collective
1194d8588912SDave May 
1195d8588912SDave May  Input Parameters:
1196d8588912SDave May .   A  - nest matrix
1197d8588912SDave May 
1198d8588912SDave May  Output Parameter:
1199629881c0SJed Brown +   M - number of rows in the nested mat
1200629881c0SJed Brown -   N - number of cols in the nested mat
1201d8588912SDave May 
1202d8588912SDave May  Notes:
1203d8588912SDave May 
1204d8588912SDave May  Level: developer
1205d8588912SDave May 
1206*bb97c47cSPierre Jolivet .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MATNEST, MatCreateNest(), MatNestGetLocalISs(),
120779798668SBarry Smith           MatNestGetISs()
1208d8588912SDave May @*/
12097087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
1210d8588912SDave May {
1211699a902aSJed Brown   PetscErrorCode ierr;
1212d8588912SDave May 
1213d8588912SDave May   PetscFunctionBegin;
1214699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
1215d8588912SDave May   PetscFunctionReturn(0);
1216d8588912SDave May }
1217d8588912SDave May 
1218f7a08781SBarry Smith static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
1219900e7ff2SJed Brown {
1220900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
1221900e7ff2SJed Brown   PetscInt i;
1222900e7ff2SJed Brown 
1223900e7ff2SJed Brown   PetscFunctionBegin;
1224900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
1225900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
1226900e7ff2SJed Brown   PetscFunctionReturn(0);
1227900e7ff2SJed Brown }
1228900e7ff2SJed Brown 
12293a4d7b9aSSatish Balay /*@C
1230900e7ff2SJed Brown  MatNestGetISs - Returns the index sets partitioning the row and column spaces
1231900e7ff2SJed Brown 
1232900e7ff2SJed Brown  Not collective
1233900e7ff2SJed Brown 
1234900e7ff2SJed Brown  Input Parameters:
1235900e7ff2SJed Brown .   A  - nest matrix
1236900e7ff2SJed Brown 
1237900e7ff2SJed Brown  Output Parameter:
1238900e7ff2SJed Brown +   rows - array of row index sets
1239900e7ff2SJed Brown -   cols - array of column index sets
1240900e7ff2SJed Brown 
1241900e7ff2SJed Brown  Level: advanced
1242900e7ff2SJed Brown 
1243900e7ff2SJed Brown  Notes:
1244900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1245900e7ff2SJed Brown 
124679798668SBarry Smith .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs(), MATNEST,
1247*bb97c47cSPierre Jolivet           MatCreateNest(), MatNestGetSubMats(), MatNestSetSubMats()
1248900e7ff2SJed Brown @*/
1249900e7ff2SJed Brown PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
1250900e7ff2SJed Brown {
1251900e7ff2SJed Brown   PetscErrorCode ierr;
1252900e7ff2SJed Brown 
1253900e7ff2SJed Brown   PetscFunctionBegin;
1254900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1255900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1256900e7ff2SJed Brown   PetscFunctionReturn(0);
1257900e7ff2SJed Brown }
1258900e7ff2SJed Brown 
1259f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
1260900e7ff2SJed Brown {
1261900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
1262900e7ff2SJed Brown   PetscInt i;
1263900e7ff2SJed Brown 
1264900e7ff2SJed Brown   PetscFunctionBegin;
1265900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
1266900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
1267900e7ff2SJed Brown   PetscFunctionReturn(0);
1268900e7ff2SJed Brown }
1269900e7ff2SJed Brown 
1270900e7ff2SJed Brown /*@C
1271900e7ff2SJed Brown  MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces
1272900e7ff2SJed Brown 
1273900e7ff2SJed Brown  Not collective
1274900e7ff2SJed Brown 
1275900e7ff2SJed Brown  Input Parameters:
1276900e7ff2SJed Brown .   A  - nest matrix
1277900e7ff2SJed Brown 
1278900e7ff2SJed Brown  Output Parameter:
12790298fd71SBarry Smith +   rows - array of row index sets (or NULL to ignore)
12800298fd71SBarry Smith -   cols - array of column index sets (or NULL to ignore)
1281900e7ff2SJed Brown 
1282900e7ff2SJed Brown  Level: advanced
1283900e7ff2SJed Brown 
1284900e7ff2SJed Brown  Notes:
1285900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1286900e7ff2SJed Brown 
1287*bb97c47cSPierre Jolivet .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs(), MatCreateNest(),
128879798668SBarry Smith           MATNEST, MatNestSetSubMats(), MatNestSetSubMat()
1289900e7ff2SJed Brown @*/
1290900e7ff2SJed Brown PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
1291900e7ff2SJed Brown {
1292900e7ff2SJed Brown   PetscErrorCode ierr;
1293900e7ff2SJed Brown 
1294900e7ff2SJed Brown   PetscFunctionBegin;
1295900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1296900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1297900e7ff2SJed Brown   PetscFunctionReturn(0);
1298900e7ff2SJed Brown }
1299900e7ff2SJed Brown 
130019fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
1301207556f9SJed Brown {
1302207556f9SJed Brown   PetscErrorCode ierr;
1303207556f9SJed Brown   PetscBool      flg;
1304207556f9SJed Brown 
1305207556f9SJed Brown   PetscFunctionBegin;
1306207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
1307207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
13082a7a6963SBarry Smith   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
130912b53f24SSatish Balay   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
1310207556f9SJed Brown   PetscFunctionReturn(0);
1311207556f9SJed Brown }
1312207556f9SJed Brown 
1313207556f9SJed Brown /*@C
13142a7a6963SBarry Smith  MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs()
1315207556f9SJed Brown 
1316207556f9SJed Brown  Not collective
1317207556f9SJed Brown 
1318207556f9SJed Brown  Input Parameters:
1319207556f9SJed Brown +  A  - nest matrix
1320207556f9SJed Brown -  vtype - type to use for creating vectors
1321207556f9SJed Brown 
1322207556f9SJed Brown  Notes:
1323207556f9SJed Brown 
1324207556f9SJed Brown  Level: developer
1325207556f9SJed Brown 
1326*bb97c47cSPierre Jolivet .seealso: MatCreateVecs(), MATNEST, MatCreateNest()
1327207556f9SJed Brown @*/
132819fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
1329207556f9SJed Brown {
1330207556f9SJed Brown   PetscErrorCode ierr;
1331207556f9SJed Brown 
1332207556f9SJed Brown   PetscFunctionBegin;
133319fd82e9SBarry Smith   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr);
1334207556f9SJed Brown   PetscFunctionReturn(0);
1335207556f9SJed Brown }
1336207556f9SJed Brown 
1337c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1338d8588912SDave May {
1339c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
1340c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
1341d8588912SDave May   PetscErrorCode ierr;
134206a1af2fSStefano Zampini   PetscBool      cong;
1343d8588912SDave May 
1344d8588912SDave May   PetscFunctionBegin;
134506a1af2fSStefano Zampini   ierr = MatReset_Nest(A);CHKERRQ(ierr);
134606a1af2fSStefano Zampini 
1347c8883902SJed Brown   s->nr = nr;
1348c8883902SJed Brown   s->nc = nc;
1349d8588912SDave May 
1350c8883902SJed Brown   /* Create space for submatrices */
1351854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->m);CHKERRQ(ierr);
1352c8883902SJed Brown   for (i=0; i<nr; i++) {
1353854ce69bSBarry Smith     ierr = PetscMalloc1(nc,&s->m[i]);CHKERRQ(ierr);
1354d8588912SDave May   }
1355c8883902SJed Brown   for (i=0; i<nr; i++) {
1356c8883902SJed Brown     for (j=0; j<nc; j++) {
1357c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
1358c8883902SJed Brown       if (a[i*nc+j]) {
1359c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1360d8588912SDave May       }
1361d8588912SDave May     }
1362d8588912SDave May   }
1363d8588912SDave May 
13648188e55aSJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1365d8588912SDave May 
1366854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->row_len);CHKERRQ(ierr);
1367854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&s->col_len);CHKERRQ(ierr);
1368c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
1369c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
1370d8588912SDave May 
137106a1af2fSStefano Zampini   ierr = PetscCalloc1(nr*nc,&s->nnzstate);CHKERRQ(ierr);
137206a1af2fSStefano Zampini   for (i=0; i<nr; i++) {
137306a1af2fSStefano Zampini     for (j=0; j<nc; j++) {
137406a1af2fSStefano Zampini       if (s->m[i][j]) {
137506a1af2fSStefano Zampini         ierr = MatGetNonzeroState(s->m[i][j],&s->nnzstate[i*nc+j]);CHKERRQ(ierr);
137606a1af2fSStefano Zampini       }
137706a1af2fSStefano Zampini     }
137806a1af2fSStefano Zampini   }
137906a1af2fSStefano Zampini 
13808188e55aSJed Brown   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1381d8588912SDave May 
1382c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1383c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1384c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1385c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1386c8883902SJed Brown 
1387c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1388c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1389c8883902SJed Brown 
139006a1af2fSStefano Zampini   /* disable operations that are not supported for non-square matrices,
139106a1af2fSStefano Zampini      or matrices for which is_row != is_col  */
139206a1af2fSStefano Zampini   ierr = MatHasCongruentLayouts(A,&cong);CHKERRQ(ierr);
139306a1af2fSStefano Zampini   if (cong && nr != nc) cong = PETSC_FALSE;
139406a1af2fSStefano Zampini   if (cong) {
139506a1af2fSStefano Zampini     for (i = 0; cong && i < nr; i++) {
1396320466b0SStefano Zampini       ierr = ISEqualUnsorted(s->isglobal.row[i],s->isglobal.col[i],&cong);CHKERRQ(ierr);
139706a1af2fSStefano Zampini     }
139806a1af2fSStefano Zampini   }
139906a1af2fSStefano Zampini   if (!cong) {
1400381b8e50SStefano Zampini     A->ops->missingdiagonal = NULL;
140106a1af2fSStefano Zampini     A->ops->getdiagonal     = NULL;
140206a1af2fSStefano Zampini     A->ops->shift           = NULL;
140306a1af2fSStefano Zampini     A->ops->diagonalset     = NULL;
140406a1af2fSStefano Zampini   }
140506a1af2fSStefano Zampini 
14061795a4d1SJed Brown   ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr);
140706a1af2fSStefano Zampini   ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr);
140806a1af2fSStefano Zampini   A->nonzerostate++;
1409d8588912SDave May   PetscFunctionReturn(0);
1410d8588912SDave May }
1411d8588912SDave May 
1412c8883902SJed Brown /*@
1413c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
1414c8883902SJed Brown 
1415c8883902SJed Brown    Collective on Mat
1416c8883902SJed Brown 
1417c8883902SJed Brown    Input Parameter:
1418ffd6319bSRichard Tran Mills +  A - nested matrix
1419c8883902SJed Brown .  nr - number of nested row blocks
14200298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1421c8883902SJed Brown .  nc - number of nested column blocks
14220298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
14230298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1424c8883902SJed Brown 
142506a1af2fSStefano Zampini    Notes: this always resets any submatrix information previously set
142606a1af2fSStefano Zampini 
1427c8883902SJed Brown    Level: advanced
1428c8883902SJed Brown 
142979798668SBarry Smith .seealso: MatCreateNest(), MATNEST, MatNestSetSubMat(), MatNestGetSubMat(), MatNestGetSubMats()
1430c8883902SJed Brown @*/
1431c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1432c8883902SJed Brown {
1433c8883902SJed Brown   PetscErrorCode ierr;
143406a1af2fSStefano Zampini   PetscInt       i;
1435c8883902SJed Brown 
1436c8883902SJed Brown   PetscFunctionBegin;
1437c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1438ce94432eSBarry Smith   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1439c8883902SJed Brown   if (nr && is_row) {
1440c8883902SJed Brown     PetscValidPointer(is_row,3);
1441c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1442c8883902SJed Brown   }
1443ce94432eSBarry Smith   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
14441664e352SJed Brown   if (nc && is_col) {
1445c8883902SJed Brown     PetscValidPointer(is_col,5);
14469b30a8f6SBarry Smith     for (i=0; i<nc; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1447c8883902SJed Brown   }
144806a1af2fSStefano Zampini   if (nr*nc > 0) PetscValidPointer(a,6);
1449c8883902SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));CHKERRQ(ierr);
1450c8883902SJed Brown   PetscFunctionReturn(0);
1451c8883902SJed Brown }
1452d8588912SDave May 
145345b6f7e9SBarry Smith static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog)
145477019fcaSJed Brown {
145577019fcaSJed Brown   PetscErrorCode ierr;
145677019fcaSJed Brown   PetscBool      flg;
145777019fcaSJed Brown   PetscInt       i,j,m,mi,*ix;
145877019fcaSJed Brown 
145977019fcaSJed Brown   PetscFunctionBegin;
1460aea6d515SStefano Zampini   *ltog = NULL;
146177019fcaSJed Brown   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
146277019fcaSJed Brown     if (islocal[i]) {
1463aea6d515SStefano Zampini       ierr = ISGetLocalSize(islocal[i],&mi);CHKERRQ(ierr);
146477019fcaSJed Brown       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
146577019fcaSJed Brown     } else {
1466aea6d515SStefano Zampini       ierr = ISGetLocalSize(isglobal[i],&mi);CHKERRQ(ierr);
146777019fcaSJed Brown     }
146877019fcaSJed Brown     m += mi;
146977019fcaSJed Brown   }
1470aea6d515SStefano Zampini   if (!flg) PetscFunctionReturn(0);
1471aea6d515SStefano Zampini 
1472785e854fSJed Brown   ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr);
1473165cd838SBarry Smith   for (i=0,m=0; i<n; i++) {
14740298fd71SBarry Smith     ISLocalToGlobalMapping smap = NULL;
1475e108cb99SStefano Zampini     Mat                    sub = NULL;
1476f6d38dbbSStefano Zampini     PetscSF                sf;
1477f6d38dbbSStefano Zampini     PetscLayout            map;
1478aea6d515SStefano Zampini     const PetscInt         *ix2;
147977019fcaSJed Brown 
1480165cd838SBarry Smith     if (!colflg) {
148177019fcaSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
148277019fcaSJed Brown     } else {
148377019fcaSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
148477019fcaSJed Brown     }
1485191fd14bSBarry Smith     if (sub) {
1486191fd14bSBarry Smith       if (!colflg) {
1487191fd14bSBarry Smith         ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);
1488191fd14bSBarry Smith       } else {
1489191fd14bSBarry Smith         ierr = MatGetLocalToGlobalMapping(sub,NULL,&smap);CHKERRQ(ierr);
1490191fd14bSBarry Smith       }
1491191fd14bSBarry Smith     }
149277019fcaSJed Brown     /*
149377019fcaSJed Brown        Now we need to extract the monolithic global indices that correspond to the given split global indices.
149477019fcaSJed 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.
149577019fcaSJed Brown     */
1496aea6d515SStefano Zampini     ierr = ISGetIndices(isglobal[i],&ix2);CHKERRQ(ierr);
1497aea6d515SStefano Zampini     if (islocal[i]) {
1498aea6d515SStefano Zampini       PetscInt *ilocal,*iremote;
1499aea6d515SStefano Zampini       PetscInt mil,nleaves;
1500aea6d515SStefano Zampini 
1501aea6d515SStefano Zampini       ierr = ISGetLocalSize(islocal[i],&mi);CHKERRQ(ierr);
1502aea6d515SStefano Zampini       if (!smap) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing local to global map");
1503aea6d515SStefano Zampini       for (j=0; j<mi; j++) ix[m+j] = j;
1504aea6d515SStefano Zampini       ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);
1505aea6d515SStefano Zampini 
1506aea6d515SStefano Zampini       /* PetscSFSetGraphLayout does not like negative indices */
1507aea6d515SStefano Zampini       ierr = PetscMalloc2(mi,&ilocal,mi,&iremote);CHKERRQ(ierr);
1508aea6d515SStefano Zampini       for (j=0, nleaves = 0; j<mi; j++) {
1509aea6d515SStefano Zampini         if (ix[m+j] < 0) continue;
1510aea6d515SStefano Zampini         ilocal[nleaves]  = j;
1511aea6d515SStefano Zampini         iremote[nleaves] = ix[m+j];
1512aea6d515SStefano Zampini         nleaves++;
1513aea6d515SStefano Zampini       }
1514aea6d515SStefano Zampini       ierr = ISGetLocalSize(isglobal[i],&mil);CHKERRQ(ierr);
1515aea6d515SStefano Zampini       ierr = PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);CHKERRQ(ierr);
1516aea6d515SStefano Zampini       ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)A),&map);CHKERRQ(ierr);
1517aea6d515SStefano Zampini       ierr = PetscLayoutSetLocalSize(map,mil);CHKERRQ(ierr);
1518f6d38dbbSStefano Zampini       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
1519aea6d515SStefano Zampini       ierr = PetscSFSetGraphLayout(sf,map,nleaves,ilocal,PETSC_USE_POINTER,iremote);CHKERRQ(ierr);
1520f6d38dbbSStefano Zampini       ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
1521f6d38dbbSStefano Zampini       ierr = PetscSFBcastBegin(sf,MPIU_INT,ix2,ix + m);CHKERRQ(ierr);
1522f6d38dbbSStefano Zampini       ierr = PetscSFBcastEnd(sf,MPIU_INT,ix2,ix + m);CHKERRQ(ierr);
1523f6d38dbbSStefano Zampini       ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1524aea6d515SStefano Zampini       ierr = PetscFree2(ilocal,iremote);CHKERRQ(ierr);
1525aea6d515SStefano Zampini     } else {
1526aea6d515SStefano Zampini       ierr = ISGetLocalSize(isglobal[i],&mi);CHKERRQ(ierr);
1527aea6d515SStefano Zampini       for (j=0; j<mi; j++) ix[m+j] = ix2[i];
1528aea6d515SStefano Zampini     }
1529aea6d515SStefano Zampini     ierr = ISRestoreIndices(isglobal[i],&ix2);CHKERRQ(ierr);
153077019fcaSJed Brown     m   += mi;
153177019fcaSJed Brown   }
1532f0413b6fSBarry Smith   ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
153377019fcaSJed Brown   PetscFunctionReturn(0);
153477019fcaSJed Brown }
153577019fcaSJed Brown 
153677019fcaSJed Brown 
1537d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1538d8588912SDave May /*
1539d8588912SDave May   nprocessors = NP
1540d8588912SDave May   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1541d8588912SDave May        proc 0: => (g_0,h_0,)
1542d8588912SDave May        proc 1: => (g_1,h_1,)
1543d8588912SDave May        ...
1544d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1545d8588912SDave May 
1546d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1547d8588912SDave May     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))
1548d8588912SDave May 
1549d8588912SDave May             proc 0:
1550d8588912SDave May     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1551d8588912SDave May             proc 1:
1552d8588912SDave May     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)
1553d8588912SDave May 
1554d8588912SDave May             proc NP-1:
1555d8588912SDave May     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1556d8588912SDave May */
1557841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1558d8588912SDave May {
1559e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
15608188e55aSJed Brown   PetscInt       i,j,offset,n,nsum,bs;
1561d8588912SDave May   PetscErrorCode ierr;
15620298fd71SBarry Smith   Mat            sub = NULL;
1563d8588912SDave May 
1564d8588912SDave May   PetscFunctionBegin;
1565854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&vs->isglobal.row);CHKERRQ(ierr);
1566854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&vs->isglobal.col);CHKERRQ(ierr);
1567d8588912SDave May   if (is_row) { /* valid IS is passed in */
1568d8588912SDave May     /* refs on is[] are incremeneted */
1569e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1570d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
157126fbe8dcSKarl Rupp 
1572e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1573d8588912SDave May     }
15742ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
15758188e55aSJed Brown     nsum = 0;
15768188e55aSJed Brown     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
15778188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1578ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
15790298fd71SBarry Smith       ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
1580ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
15818188e55aSJed Brown       nsum += n;
15828188e55aSJed Brown     }
1583ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
158430bc264bSJed Brown     offset -= nsum;
1585e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1586f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
15870298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
15882ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1589ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1590e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
15912ae74bdbSJed Brown       offset += n;
1592d8588912SDave May     }
1593d8588912SDave May   }
1594d8588912SDave May 
1595d8588912SDave May   if (is_col) { /* valid IS is passed in */
1596d8588912SDave May     /* refs on is[] are incremeneted */
1597e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1598d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
159926fbe8dcSKarl Rupp 
1600e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1601d8588912SDave May     }
16022ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
16032ae74bdbSJed Brown     offset = A->cmap->rstart;
16048188e55aSJed Brown     nsum   = 0;
16058188e55aSJed Brown     for (j=0; j<vs->nc; j++) {
16068188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1607ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
16080298fd71SBarry Smith       ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
1609ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
16108188e55aSJed Brown       nsum += n;
16118188e55aSJed Brown     }
1612ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
161330bc264bSJed Brown     offset -= nsum;
1614e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1615f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
16160298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
16172ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1618ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1619e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
16202ae74bdbSJed Brown       offset += n;
1621d8588912SDave May     }
1622d8588912SDave May   }
1623e2d7f03fSJed Brown 
1624e2d7f03fSJed Brown   /* Set up the local ISs */
1625785e854fSJed Brown   ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
1626785e854fSJed Brown   ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
1627e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1628e2d7f03fSJed Brown     IS                     isloc;
16290298fd71SBarry Smith     ISLocalToGlobalMapping rmap = NULL;
1630e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1631e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
16320298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);}
1633207556f9SJed Brown     if (rmap) {
1634e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1635e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1636e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1637e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1638207556f9SJed Brown     } else {
1639207556f9SJed Brown       nlocal = 0;
16400298fd71SBarry Smith       isloc  = NULL;
1641207556f9SJed Brown     }
1642e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1643e2d7f03fSJed Brown     offset            += nlocal;
1644e2d7f03fSJed Brown   }
16458188e55aSJed Brown   for (i=0,offset=0; i<vs->nc; i++) {
1646e2d7f03fSJed Brown     IS                     isloc;
16470298fd71SBarry Smith     ISLocalToGlobalMapping cmap = NULL;
1648e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1649e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
16500298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);}
1651207556f9SJed Brown     if (cmap) {
1652e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1653e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1654e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1655e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1656207556f9SJed Brown     } else {
1657207556f9SJed Brown       nlocal = 0;
16580298fd71SBarry Smith       isloc  = NULL;
1659207556f9SJed Brown     }
1660e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1661e2d7f03fSJed Brown     offset            += nlocal;
1662e2d7f03fSJed Brown   }
16630189643fSJed Brown 
166477019fcaSJed Brown   /* Set up the aggregate ISLocalToGlobalMapping */
166577019fcaSJed Brown   {
166645b6f7e9SBarry Smith     ISLocalToGlobalMapping rmap,cmap;
166745b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);CHKERRQ(ierr);
166845b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);CHKERRQ(ierr);
166977019fcaSJed Brown     if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);}
167077019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
167177019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
167277019fcaSJed Brown   }
167377019fcaSJed Brown 
16740189643fSJed Brown #if defined(PETSC_USE_DEBUG)
16750189643fSJed Brown   for (i=0; i<vs->nr; i++) {
16760189643fSJed Brown     for (j=0; j<vs->nc; j++) {
16770189643fSJed Brown       PetscInt m,n,M,N,mi,ni,Mi,Ni;
16780189643fSJed Brown       Mat      B = vs->m[i][j];
16790189643fSJed Brown       if (!B) continue;
16800189643fSJed Brown       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
16810189643fSJed Brown       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
16820189643fSJed Brown       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
16830189643fSJed Brown       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
16840189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
16850189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
1686ce94432eSBarry Smith       if (M != Mi || N != Ni) SETERRQ6(PetscObjectComm((PetscObject)sub),PETSC_ERR_ARG_INCOMP,"Global sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",M,N,i,j,Mi,Ni);
1687ce94432eSBarry Smith       if (m != mi || n != ni) SETERRQ6(PetscObjectComm((PetscObject)sub),PETSC_ERR_ARG_INCOMP,"Local sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",m,n,i,j,mi,ni);
16880189643fSJed Brown     }
16890189643fSJed Brown   }
16900189643fSJed Brown #endif
1691a061e289SJed Brown 
1692a061e289SJed Brown   /* Set A->assembled if all non-null blocks are currently assembled */
1693a061e289SJed Brown   for (i=0; i<vs->nr; i++) {
1694a061e289SJed Brown     for (j=0; j<vs->nc; j++) {
1695a061e289SJed Brown       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0);
1696a061e289SJed Brown     }
1697a061e289SJed Brown   }
1698a061e289SJed Brown   A->assembled = PETSC_TRUE;
1699d8588912SDave May   PetscFunctionReturn(0);
1700d8588912SDave May }
1701d8588912SDave May 
170245c38901SJed Brown /*@C
1703659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1704659c6bb0SJed Brown 
1705659c6bb0SJed Brown    Collective on Mat
1706659c6bb0SJed Brown 
1707659c6bb0SJed Brown    Input Parameter:
1708659c6bb0SJed Brown +  comm - Communicator for the new Mat
1709659c6bb0SJed Brown .  nr - number of nested row blocks
17100298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1711659c6bb0SJed Brown .  nc - number of nested column blocks
17120298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
17130298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1714659c6bb0SJed Brown 
1715659c6bb0SJed Brown    Output Parameter:
1716659c6bb0SJed Brown .  B - new matrix
1717659c6bb0SJed Brown 
1718659c6bb0SJed Brown    Level: advanced
1719659c6bb0SJed Brown 
172079798668SBarry Smith .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST, MatNestSetSubMat(),
172179798668SBarry Smith           MatNestGetSubMat(), MatNestGetLocalISs(), MatNestGetSize(),
172279798668SBarry Smith           MatNestGetISs(), MatNestSetSubMats(), MatNestGetSubMats()
1723659c6bb0SJed Brown @*/
17247087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1725d8588912SDave May {
1726d8588912SDave May   Mat            A;
1727d8588912SDave May   PetscErrorCode ierr;
1728d8588912SDave May 
1729d8588912SDave May   PetscFunctionBegin;
1730c8883902SJed Brown   *B   = 0;
1731d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1732c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
173391a28eb3SBarry Smith   A->preallocated = PETSC_TRUE;
1734c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1735d8588912SDave May   *B   = A;
1736d8588912SDave May   PetscFunctionReturn(0);
1737d8588912SDave May }
1738659c6bb0SJed Brown 
1739b68353e5Sstefano_zampini static PetscErrorCode MatConvert_Nest_SeqAIJ_fast(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1740b68353e5Sstefano_zampini {
1741b68353e5Sstefano_zampini   Mat_Nest       *nest = (Mat_Nest*)A->data;
174223875855Sstefano_zampini   Mat            *trans;
1743b68353e5Sstefano_zampini   PetscScalar    **avv;
1744b68353e5Sstefano_zampini   PetscScalar    *vv;
1745b68353e5Sstefano_zampini   PetscInt       **aii,**ajj;
1746b68353e5Sstefano_zampini   PetscInt       *ii,*jj,*ci;
1747b68353e5Sstefano_zampini   PetscInt       nr,nc,nnz,i,j;
1748b68353e5Sstefano_zampini   PetscBool      done;
1749b68353e5Sstefano_zampini   PetscErrorCode ierr;
1750b68353e5Sstefano_zampini 
1751b68353e5Sstefano_zampini   PetscFunctionBegin;
1752b68353e5Sstefano_zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
1753b68353e5Sstefano_zampini   if (reuse == MAT_REUSE_MATRIX) {
1754b68353e5Sstefano_zampini     PetscInt rnr;
1755b68353e5Sstefano_zampini 
1756b68353e5Sstefano_zampini     ierr = MatGetRowIJ(*newmat,0,PETSC_FALSE,PETSC_FALSE,&rnr,(const PetscInt**)&ii,(const PetscInt**)&jj,&done);CHKERRQ(ierr);
1757b68353e5Sstefano_zampini     if (!done) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatGetRowIJ");
1758b68353e5Sstefano_zampini     if (rnr != nr) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of rows");
1759b68353e5Sstefano_zampini     ierr = MatSeqAIJGetArray(*newmat,&vv);CHKERRQ(ierr);
1760b68353e5Sstefano_zampini   }
1761b68353e5Sstefano_zampini   /* extract CSR for nested SeqAIJ matrices */
1762b68353e5Sstefano_zampini   nnz  = 0;
176323875855Sstefano_zampini   ierr = PetscCalloc4(nest->nr*nest->nc,&aii,nest->nr*nest->nc,&ajj,nest->nr*nest->nc,&avv,nest->nr*nest->nc,&trans);CHKERRQ(ierr);
1764b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1765b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1766b68353e5Sstefano_zampini       Mat B = nest->m[i][j];
1767b68353e5Sstefano_zampini       if (B) {
1768b68353e5Sstefano_zampini         PetscScalar *naa;
1769b68353e5Sstefano_zampini         PetscInt    *nii,*njj,nnr;
177023875855Sstefano_zampini         PetscBool   istrans;
1771b68353e5Sstefano_zampini 
177223875855Sstefano_zampini         ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);CHKERRQ(ierr);
177323875855Sstefano_zampini         if (istrans) {
177423875855Sstefano_zampini           Mat Bt;
177523875855Sstefano_zampini 
177623875855Sstefano_zampini           ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr);
177723875855Sstefano_zampini           ierr = MatTranspose(Bt,MAT_INITIAL_MATRIX,&trans[i*nest->nc+j]);CHKERRQ(ierr);
177823875855Sstefano_zampini           B    = trans[i*nest->nc+j];
177923875855Sstefano_zampini         }
1780b68353e5Sstefano_zampini         ierr = MatGetRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&nii,(const PetscInt**)&njj,&done);CHKERRQ(ierr);
1781b68353e5Sstefano_zampini         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatGetRowIJ");
1782b68353e5Sstefano_zampini         ierr = MatSeqAIJGetArray(B,&naa);CHKERRQ(ierr);
1783b68353e5Sstefano_zampini         nnz += nii[nnr];
1784b68353e5Sstefano_zampini 
1785b68353e5Sstefano_zampini         aii[i*nest->nc+j] = nii;
1786b68353e5Sstefano_zampini         ajj[i*nest->nc+j] = njj;
1787b68353e5Sstefano_zampini         avv[i*nest->nc+j] = naa;
1788b68353e5Sstefano_zampini       }
1789b68353e5Sstefano_zampini     }
1790b68353e5Sstefano_zampini   }
1791b68353e5Sstefano_zampini   if (reuse != MAT_REUSE_MATRIX) {
1792b68353e5Sstefano_zampini     ierr = PetscMalloc1(nr+1,&ii);CHKERRQ(ierr);
1793b68353e5Sstefano_zampini     ierr = PetscMalloc1(nnz,&jj);CHKERRQ(ierr);
1794b68353e5Sstefano_zampini     ierr = PetscMalloc1(nnz,&vv);CHKERRQ(ierr);
1795b68353e5Sstefano_zampini   } else {
1796b68353e5Sstefano_zampini     if (nnz != ii[nr]) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of nonzeros");
1797b68353e5Sstefano_zampini   }
1798b68353e5Sstefano_zampini 
1799b68353e5Sstefano_zampini   /* new row pointer */
1800580bdb30SBarry Smith   ierr = PetscArrayzero(ii,nr+1);CHKERRQ(ierr);
1801b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1802b68353e5Sstefano_zampini     PetscInt       ncr,rst;
1803b68353e5Sstefano_zampini 
1804b68353e5Sstefano_zampini     ierr = ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);CHKERRQ(ierr);
1805b68353e5Sstefano_zampini     ierr = ISGetLocalSize(nest->isglobal.row[i],&ncr);CHKERRQ(ierr);
1806b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1807b68353e5Sstefano_zampini       if (aii[i*nest->nc+j]) {
1808b68353e5Sstefano_zampini         PetscInt    *nii = aii[i*nest->nc+j];
1809b68353e5Sstefano_zampini         PetscInt    ir;
1810b68353e5Sstefano_zampini 
1811b68353e5Sstefano_zampini         for (ir=rst; ir<ncr+rst; ++ir) {
1812b68353e5Sstefano_zampini           ii[ir+1] += nii[1]-nii[0];
1813b68353e5Sstefano_zampini           nii++;
1814b68353e5Sstefano_zampini         }
1815b68353e5Sstefano_zampini       }
1816b68353e5Sstefano_zampini     }
1817b68353e5Sstefano_zampini   }
1818b68353e5Sstefano_zampini   for (i=0; i<nr; i++) ii[i+1] += ii[i];
1819b68353e5Sstefano_zampini 
1820b68353e5Sstefano_zampini   /* construct CSR for the new matrix */
1821b68353e5Sstefano_zampini   ierr = PetscCalloc1(nr,&ci);CHKERRQ(ierr);
1822b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1823b68353e5Sstefano_zampini     PetscInt       ncr,rst;
1824b68353e5Sstefano_zampini 
1825b68353e5Sstefano_zampini     ierr = ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);CHKERRQ(ierr);
1826b68353e5Sstefano_zampini     ierr = ISGetLocalSize(nest->isglobal.row[i],&ncr);CHKERRQ(ierr);
1827b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1828b68353e5Sstefano_zampini       if (aii[i*nest->nc+j]) {
1829b68353e5Sstefano_zampini         PetscScalar *nvv = avv[i*nest->nc+j];
1830b68353e5Sstefano_zampini         PetscInt    *nii = aii[i*nest->nc+j];
1831b68353e5Sstefano_zampini         PetscInt    *njj = ajj[i*nest->nc+j];
1832b68353e5Sstefano_zampini         PetscInt    ir,cst;
1833b68353e5Sstefano_zampini 
1834b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.col[j],&cst,NULL);CHKERRQ(ierr);
1835b68353e5Sstefano_zampini         for (ir=rst; ir<ncr+rst; ++ir) {
1836b68353e5Sstefano_zampini           PetscInt ij,rsize = nii[1]-nii[0],ist = ii[ir]+ci[ir];
1837b68353e5Sstefano_zampini 
1838b68353e5Sstefano_zampini           for (ij=0;ij<rsize;ij++) {
1839b68353e5Sstefano_zampini             jj[ist+ij] = *njj+cst;
1840b68353e5Sstefano_zampini             vv[ist+ij] = *nvv;
1841b68353e5Sstefano_zampini             njj++;
1842b68353e5Sstefano_zampini             nvv++;
1843b68353e5Sstefano_zampini           }
1844b68353e5Sstefano_zampini           ci[ir] += rsize;
1845b68353e5Sstefano_zampini           nii++;
1846b68353e5Sstefano_zampini         }
1847b68353e5Sstefano_zampini       }
1848b68353e5Sstefano_zampini     }
1849b68353e5Sstefano_zampini   }
1850b68353e5Sstefano_zampini   ierr = PetscFree(ci);CHKERRQ(ierr);
1851b68353e5Sstefano_zampini 
1852b68353e5Sstefano_zampini   /* restore info */
1853b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1854b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1855b68353e5Sstefano_zampini       Mat B = nest->m[i][j];
1856b68353e5Sstefano_zampini       if (B) {
1857b68353e5Sstefano_zampini         PetscInt nnr = 0, k = i*nest->nc+j;
185823875855Sstefano_zampini 
185923875855Sstefano_zampini         B    = (trans[k] ? trans[k] : B);
1860b68353e5Sstefano_zampini         ierr = MatRestoreRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&aii[k],(const PetscInt**)&ajj[k],&done);CHKERRQ(ierr);
1861b68353e5Sstefano_zampini         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatRestoreRowIJ");
1862b68353e5Sstefano_zampini         ierr = MatSeqAIJRestoreArray(B,&avv[k]);CHKERRQ(ierr);
186323875855Sstefano_zampini         ierr = MatDestroy(&trans[k]);CHKERRQ(ierr);
1864b68353e5Sstefano_zampini       }
1865b68353e5Sstefano_zampini     }
1866b68353e5Sstefano_zampini   }
186723875855Sstefano_zampini   ierr = PetscFree4(aii,ajj,avv,trans);CHKERRQ(ierr);
1868b68353e5Sstefano_zampini 
1869b68353e5Sstefano_zampini   /* finalize newmat */
1870b68353e5Sstefano_zampini   if (reuse == MAT_INITIAL_MATRIX) {
1871b68353e5Sstefano_zampini     ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,newmat);CHKERRQ(ierr);
1872b68353e5Sstefano_zampini   } else if (reuse == MAT_INPLACE_MATRIX) {
1873b68353e5Sstefano_zampini     Mat B;
1874b68353e5Sstefano_zampini 
1875b68353e5Sstefano_zampini     ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,&B);CHKERRQ(ierr);
1876b68353e5Sstefano_zampini     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
1877b68353e5Sstefano_zampini   }
1878b68353e5Sstefano_zampini   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1879b68353e5Sstefano_zampini   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1880b68353e5Sstefano_zampini   {
1881b68353e5Sstefano_zampini     Mat_SeqAIJ *a = (Mat_SeqAIJ*)((*newmat)->data);
1882b68353e5Sstefano_zampini     a->free_a     = PETSC_TRUE;
1883b68353e5Sstefano_zampini     a->free_ij    = PETSC_TRUE;
1884b68353e5Sstefano_zampini   }
1885b68353e5Sstefano_zampini   PetscFunctionReturn(0);
1886b68353e5Sstefano_zampini }
1887b68353e5Sstefano_zampini 
1888cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1889629c3df2SDmitry Karpeev {
1890629c3df2SDmitry Karpeev   PetscErrorCode ierr;
1891629c3df2SDmitry Karpeev   Mat_Nest       *nest = (Mat_Nest*)A->data;
189283b1a929SMark Adams   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz,rstart;
1893649b366bSFande Kong   PetscInt       cstart,cend;
1894b68353e5Sstefano_zampini   PetscMPIInt    size;
1895629c3df2SDmitry Karpeev   Mat            C;
1896629c3df2SDmitry Karpeev 
1897629c3df2SDmitry Karpeev   PetscFunctionBegin;
1898b68353e5Sstefano_zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1899b68353e5Sstefano_zampini   if (size == 1) { /* look for a special case with SeqAIJ matrices and strided-1, contiguous, blocks */
1900b68353e5Sstefano_zampini     PetscInt  nf;
1901b68353e5Sstefano_zampini     PetscBool fast;
1902b68353e5Sstefano_zampini 
1903b68353e5Sstefano_zampini     ierr = PetscStrcmp(newtype,MATAIJ,&fast);CHKERRQ(ierr);
1904b68353e5Sstefano_zampini     if (!fast) {
1905b68353e5Sstefano_zampini       ierr = PetscStrcmp(newtype,MATSEQAIJ,&fast);CHKERRQ(ierr);
1906b68353e5Sstefano_zampini     }
1907b68353e5Sstefano_zampini     for (i=0; i<nest->nr && fast; ++i) {
1908b68353e5Sstefano_zampini       for (j=0; j<nest->nc && fast; ++j) {
1909b68353e5Sstefano_zampini         Mat B = nest->m[i][j];
1910b68353e5Sstefano_zampini         if (B) {
1911b68353e5Sstefano_zampini           ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQAIJ,&fast);CHKERRQ(ierr);
191223875855Sstefano_zampini           if (!fast) {
191323875855Sstefano_zampini             PetscBool istrans;
191423875855Sstefano_zampini 
191523875855Sstefano_zampini             ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);CHKERRQ(ierr);
191623875855Sstefano_zampini             if (istrans) {
191723875855Sstefano_zampini               Mat Bt;
191823875855Sstefano_zampini 
191923875855Sstefano_zampini               ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr);
192023875855Sstefano_zampini               ierr = PetscObjectTypeCompare((PetscObject)Bt,MATSEQAIJ,&fast);CHKERRQ(ierr);
192123875855Sstefano_zampini             }
1922b68353e5Sstefano_zampini           }
1923b68353e5Sstefano_zampini         }
1924b68353e5Sstefano_zampini       }
1925b68353e5Sstefano_zampini     }
1926b68353e5Sstefano_zampini     for (i=0, nf=0; i<nest->nr && fast; ++i) {
1927b68353e5Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)nest->isglobal.row[i],ISSTRIDE,&fast);CHKERRQ(ierr);
1928b68353e5Sstefano_zampini       if (fast) {
1929b68353e5Sstefano_zampini         PetscInt f,s;
1930b68353e5Sstefano_zampini 
1931b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.row[i],&f,&s);CHKERRQ(ierr);
1932b68353e5Sstefano_zampini         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1933b68353e5Sstefano_zampini         else {
1934b68353e5Sstefano_zampini           ierr = ISGetSize(nest->isglobal.row[i],&f);CHKERRQ(ierr);
1935b68353e5Sstefano_zampini           nf  += f;
1936b68353e5Sstefano_zampini         }
1937b68353e5Sstefano_zampini       }
1938b68353e5Sstefano_zampini     }
1939b68353e5Sstefano_zampini     for (i=0, nf=0; i<nest->nc && fast; ++i) {
1940b68353e5Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)nest->isglobal.col[i],ISSTRIDE,&fast);CHKERRQ(ierr);
1941b68353e5Sstefano_zampini       if (fast) {
1942b68353e5Sstefano_zampini         PetscInt f,s;
1943b68353e5Sstefano_zampini 
1944b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.col[i],&f,&s);CHKERRQ(ierr);
1945b68353e5Sstefano_zampini         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1946b68353e5Sstefano_zampini         else {
1947b68353e5Sstefano_zampini           ierr = ISGetSize(nest->isglobal.col[i],&f);CHKERRQ(ierr);
1948b68353e5Sstefano_zampini           nf  += f;
1949b68353e5Sstefano_zampini         }
1950b68353e5Sstefano_zampini       }
1951b68353e5Sstefano_zampini     }
1952b68353e5Sstefano_zampini     if (fast) {
1953b68353e5Sstefano_zampini       ierr = MatConvert_Nest_SeqAIJ_fast(A,newtype,reuse,newmat);CHKERRQ(ierr);
1954b68353e5Sstefano_zampini       PetscFunctionReturn(0);
1955b68353e5Sstefano_zampini     }
1956b68353e5Sstefano_zampini   }
1957629c3df2SDmitry Karpeev   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1958629c3df2SDmitry Karpeev   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
1959649b366bSFande Kong   ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr);
1960629c3df2SDmitry Karpeev   switch (reuse) {
1961629c3df2SDmitry Karpeev   case MAT_INITIAL_MATRIX:
1962ce94432eSBarry Smith     ierr    = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
1963629c3df2SDmitry Karpeev     ierr    = MatSetType(C,newtype);CHKERRQ(ierr);
1964629c3df2SDmitry Karpeev     ierr    = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
1965629c3df2SDmitry Karpeev     *newmat = C;
1966629c3df2SDmitry Karpeev     break;
1967629c3df2SDmitry Karpeev   case MAT_REUSE_MATRIX:
1968629c3df2SDmitry Karpeev     C = *newmat;
1969629c3df2SDmitry Karpeev     break;
1970ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1971629c3df2SDmitry Karpeev   }
1972785e854fSJed Brown   ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr);
1973629c3df2SDmitry Karpeev   onnz = dnnz + m;
1974629c3df2SDmitry Karpeev   for (k=0; k<m; k++) {
1975629c3df2SDmitry Karpeev     dnnz[k] = 0;
1976629c3df2SDmitry Karpeev     onnz[k] = 0;
1977629c3df2SDmitry Karpeev   }
1978629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1979629c3df2SDmitry Karpeev     IS             bNis;
1980629c3df2SDmitry Karpeev     PetscInt       bN;
1981629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1982629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1983629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1984629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr);
1985629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1986629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1987629c3df2SDmitry Karpeev       PetscSF        bmsf;
1988649b366bSFande Kong       PetscSFNode    *iremote;
1989629c3df2SDmitry Karpeev       Mat            B;
1990649b366bSFande Kong       PetscInt       bm, *sub_dnnz,*sub_onnz, br;
1991629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1992629c3df2SDmitry Karpeev       B = nest->m[i][j];
1993629c3df2SDmitry Karpeev       if (!B) continue;
1994629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1995629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1996ce94432eSBarry Smith       ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr);
1997649b366bSFande Kong       ierr = PetscMalloc1(bm,&iremote);CHKERRQ(ierr);
1998649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_dnnz);CHKERRQ(ierr);
1999649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_onnz);CHKERRQ(ierr);
2000649b366bSFande Kong       for (k = 0; k < bm; ++k){
2001649b366bSFande Kong     	sub_dnnz[k] = 0;
2002649b366bSFande Kong     	sub_onnz[k] = 0;
2003649b366bSFande Kong       }
2004629c3df2SDmitry Karpeev       /*
2005629c3df2SDmitry Karpeev        Locate the owners for all of the locally-owned global row indices for this row block.
2006629c3df2SDmitry Karpeev        These determine the roots of PetscSF used to communicate preallocation data to row owners.
2007629c3df2SDmitry Karpeev        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
2008629c3df2SDmitry Karpeev        */
200983b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
2010629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
2011131c27b5Sprj-         PetscInt       row = bmindices[br], brncols, col;
2012629c3df2SDmitry Karpeev         const PetscInt *brcols;
2013a4b3d3acSMatthew G Knepley         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
2014131c27b5Sprj-         PetscMPIInt    rowowner = 0;
2015629c3df2SDmitry Karpeev         ierr      = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr);
2016649b366bSFande Kong         /* how many roots  */
2017649b366bSFande Kong         iremote[br].rank = rowowner; iremote[br].index = rowrel;           /* edge from bmdnnz to dnnz */
2018649b366bSFande Kong         /* get nonzero pattern */
201983b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
2020629c3df2SDmitry Karpeev         for (k=0; k<brncols; k++) {
2021629c3df2SDmitry Karpeev           col  = bNindices[brcols[k]];
2022649b366bSFande Kong           if (col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]) {
2023649b366bSFande Kong             sub_dnnz[br]++;
2024649b366bSFande Kong           } else {
2025649b366bSFande Kong             sub_onnz[br]++;
2026649b366bSFande Kong           }
2027629c3df2SDmitry Karpeev         }
202883b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
2029629c3df2SDmitry Karpeev       }
2030629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
2031629c3df2SDmitry Karpeev       /* bsf will have to take care of disposing of bedges. */
2032649b366bSFande Kong       ierr = PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
2033649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
2034649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
2035649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
2036649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
2037649b366bSFande Kong       ierr = PetscFree(sub_dnnz);CHKERRQ(ierr);
2038649b366bSFande Kong       ierr = PetscFree(sub_onnz);CHKERRQ(ierr);
2039629c3df2SDmitry Karpeev       ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr);
2040629c3df2SDmitry Karpeev     }
204122d28d08SBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
2042629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
204365a4a0a3Sstefano_zampini   }
204465a4a0a3Sstefano_zampini   /* Resize preallocation if overestimated */
204565a4a0a3Sstefano_zampini   for (i=0;i<m;i++) {
204665a4a0a3Sstefano_zampini     dnnz[i] = PetscMin(dnnz[i],A->cmap->n);
204765a4a0a3Sstefano_zampini     onnz[i] = PetscMin(onnz[i],A->cmap->N - A->cmap->n);
2048629c3df2SDmitry Karpeev   }
2049629c3df2SDmitry Karpeev   ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr);
2050629c3df2SDmitry Karpeev   ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr);
2051629c3df2SDmitry Karpeev   ierr = PetscFree(dnnz);CHKERRQ(ierr);
2052629c3df2SDmitry Karpeev 
2053629c3df2SDmitry Karpeev   /* Fill by row */
2054629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
2055629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
2056629c3df2SDmitry Karpeev     IS             bNis;
2057629c3df2SDmitry Karpeev     PetscInt       bN;
2058629c3df2SDmitry Karpeev     const PetscInt *bNindices;
2059629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
2060629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr);
2061629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
2062629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
2063629c3df2SDmitry Karpeev       Mat            B;
2064629c3df2SDmitry Karpeev       PetscInt       bm, br;
2065629c3df2SDmitry Karpeev       const PetscInt *bmindices;
2066629c3df2SDmitry Karpeev       B = nest->m[i][j];
2067629c3df2SDmitry Karpeev       if (!B) continue;
2068629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
2069629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
207083b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
2071629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
2072629c3df2SDmitry Karpeev         PetscInt          row = bmindices[br], brncols,  *cols;
2073629c3df2SDmitry Karpeev         const PetscInt    *brcols;
2074629c3df2SDmitry Karpeev         const PetscScalar *brcoldata;
207583b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
2076785e854fSJed Brown         ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr);
207726fbe8dcSKarl Rupp         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
2078629c3df2SDmitry Karpeev         /*
2079629c3df2SDmitry Karpeev           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
2080629c3df2SDmitry Karpeev           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
2081629c3df2SDmitry Karpeev          */
2082a2ea699eSBarry Smith         ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr);
208383b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
2084629c3df2SDmitry Karpeev         ierr = PetscFree(cols);CHKERRQ(ierr);
2085629c3df2SDmitry Karpeev       }
2086629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
2087629c3df2SDmitry Karpeev     }
2088a2ea699eSBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
2089629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
2090629c3df2SDmitry Karpeev   }
2091629c3df2SDmitry Karpeev   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2092629c3df2SDmitry Karpeev   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
2093629c3df2SDmitry Karpeev   PetscFunctionReturn(0);
2094629c3df2SDmitry Karpeev }
2095629c3df2SDmitry Karpeev 
20968b7d3b4bSBarry Smith PetscErrorCode MatHasOperation_Nest(Mat mat,MatOperation op,PetscBool *has)
20978b7d3b4bSBarry Smith {
20988b7d3b4bSBarry Smith   Mat_Nest       *bA = (Mat_Nest*)mat->data;
20998b7d3b4bSBarry Smith   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
21008b7d3b4bSBarry Smith   PetscBool      flg;
210152c5f739Sprj-   PetscErrorCode ierr;
210252c5f739Sprj-   PetscFunctionBegin;
21038b7d3b4bSBarry Smith 
210452c5f739Sprj-   *has = PETSC_FALSE;
210552c5f739Sprj-   if (op == MATOP_MULT_TRANSPOSE || op == MATOP_MAT_MULT) {
21068b7d3b4bSBarry Smith     for (j=0; j<nc; j++) {
21078b7d3b4bSBarry Smith       for (i=0; i<nr; i++) {
21088b7d3b4bSBarry Smith         if (!bA->m[i][j]) continue;
210952c5f739Sprj-         ierr = MatHasOperation(bA->m[i][j],op,&flg);CHKERRQ(ierr);
21108b7d3b4bSBarry Smith         if (!flg) PetscFunctionReturn(0);
21118b7d3b4bSBarry Smith       }
21128b7d3b4bSBarry Smith     }
21138b7d3b4bSBarry Smith   }
211452c5f739Sprj-   if (((void**)mat->ops)[op] || (op == MATOP_MAT_MULT && flg)) *has = PETSC_TRUE;
21158b7d3b4bSBarry Smith   PetscFunctionReturn(0);
21168b7d3b4bSBarry Smith }
21178b7d3b4bSBarry Smith 
2118659c6bb0SJed Brown /*MC
2119659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
2120659c6bb0SJed Brown 
2121659c6bb0SJed Brown   Level: intermediate
2122659c6bb0SJed Brown 
2123659c6bb0SJed Brown   Notes:
2124659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
2125659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
2126950540a4SJed Brown   It is usually used with DMComposite and DMCreateMatrix()
2127659c6bb0SJed Brown 
21288b7d3b4bSBarry Smith   Each of the submatrices lives on the same MPI communicator as the original nest matrix (though they can have zero
21298b7d3b4bSBarry Smith   rows/columns on some processes.) Thus this is not meant for cases where the submatrices live on far fewer processes
21308b7d3b4bSBarry Smith   than the nest matrix.
21318b7d3b4bSBarry Smith 
213279798668SBarry Smith .seealso: MatCreate(), MatType, MatCreateNest(), MatNestSetSubMat(), MatNestGetSubMat(),
213379798668SBarry Smith           VecCreateNest(), DMCreateMatrix(), DMCOMPOSITE, MatNestSetVecType(), MatNestGetLocalISs(),
213479798668SBarry Smith           MatNestGetISs(), MatNestSetSubMats(), MatNestGetSubMats()
2135659c6bb0SJed Brown M*/
21368cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
2137c8883902SJed Brown {
2138c8883902SJed Brown   Mat_Nest       *s;
2139c8883902SJed Brown   PetscErrorCode ierr;
2140c8883902SJed Brown 
2141c8883902SJed Brown   PetscFunctionBegin;
2142b00a9115SJed Brown   ierr    = PetscNewLog(A,&s);CHKERRQ(ierr);
2143c8883902SJed Brown   A->data = (void*)s;
2144e7c19651SJed Brown 
2145e7c19651SJed Brown   s->nr            = -1;
2146e7c19651SJed Brown   s->nc            = -1;
21470298fd71SBarry Smith   s->m             = NULL;
2148e7c19651SJed Brown   s->splitassembly = PETSC_FALSE;
2149c8883902SJed Brown 
2150c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
215126fbe8dcSKarl Rupp 
2152c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
21539194d70fSJed Brown   A->ops->multadd               = MatMultAdd_Nest;
2154c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
21559194d70fSJed Brown   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
2156f8170845SAlex Fikl   A->ops->transpose             = MatTranspose_Nest;
2157c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
2158c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
2159c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
2160c222c20dSDavid Ham   A->ops->copy                  = MatCopy_Nest;
21616e76ffeaSPierre Jolivet   A->ops->axpy                  = MatAXPY_Nest;
2162c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
21637dae84e0SHong Zhang   A->ops->createsubmatrix       = MatCreateSubMatrix_Nest;
2164c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
2165c8883902SJed Brown   A->ops->view                  = MatView_Nest;
2166c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
2167c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
2168c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
2169429bac76SJed Brown   A->ops->getdiagonal           = MatGetDiagonal_Nest;
2170429bac76SJed Brown   A->ops->diagonalscale         = MatDiagonalScale_Nest;
2171a061e289SJed Brown   A->ops->scale                 = MatScale_Nest;
2172a061e289SJed Brown   A->ops->shift                 = MatShift_Nest;
217313135bc6SAlex Fikl   A->ops->diagonalset           = MatDiagonalSet_Nest;
2174f8170845SAlex Fikl   A->ops->setrandom             = MatSetRandom_Nest;
21758b7d3b4bSBarry Smith   A->ops->hasoperation          = MatHasOperation_Nest;
2176381b8e50SStefano Zampini   A->ops->missingdiagonal       = MatMissingDiagonal_Nest;
2177c8883902SJed Brown 
2178c8883902SJed Brown   A->spptr        = 0;
2179c8883902SJed Brown   A->assembled    = PETSC_FALSE;
2180c8883902SJed Brown 
2181c8883902SJed Brown   /* expose Nest api's */
2182bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",        MatNestGetSubMat_Nest);CHKERRQ(ierr);
2183bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",        MatNestSetSubMat_Nest);CHKERRQ(ierr);
2184bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",       MatNestGetSubMats_Nest);CHKERRQ(ierr);
2185bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",          MatNestGetSize_Nest);CHKERRQ(ierr);
2186bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",           MatNestGetISs_Nest);CHKERRQ(ierr);
2187bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",      MatNestGetLocalISs_Nest);CHKERRQ(ierr);
2188bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",       MatNestSetVecType_Nest);CHKERRQ(ierr);
2189bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",       MatNestSetSubMats_Nest);CHKERRQ(ierr);
21900899c546SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_mpiaij_C",  MatConvert_Nest_AIJ);CHKERRQ(ierr);
21910899c546SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_seqaij_C",  MatConvert_Nest_AIJ);CHKERRQ(ierr);
219283b1a929SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",     MatConvert_Nest_AIJ);CHKERRQ(ierr);
21935e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",      MatConvert_Nest_IS);CHKERRQ(ierr);
21944222ddf1SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_seqdense_C",MatProductSetFromOptions_Nest_Dense);CHKERRQ(ierr);
21954222ddf1SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_mpidense_C",MatProductSetFromOptions_Nest_Dense);CHKERRQ(ierr);
21964222ddf1SHong Zhang   ierr = PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_dense_C",MatProductSetFromOptions_Nest_Dense);CHKERRQ(ierr);
2197c8883902SJed Brown 
2198c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
2199c8883902SJed Brown   PetscFunctionReturn(0);
2200c8883902SJed Brown }
2201