xref: /petsc/src/mat/impls/nest/matnest.c (revision 06a1af2f5b4ece804d199086219e477db115facf)
1d8588912SDave May 
2aaa7dc30SBarry Smith #include <../src/mat/impls/nest/matnestimpl.h> /*I   "petscmat.h"   I*/
3b68353e5Sstefano_zampini #include <../src/mat/impls/aij/seq/aij.h>
40c312b8eSJed Brown #include <petscsf.h>
5d8588912SDave May 
6c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
7*06a1af2fSStefano Zampini static PetscErrorCode MatCreateVecs_Nest(Mat,Vec*,Vec*);
8*06a1af2fSStefano Zampini static PetscErrorCode MatReset_Nest(Mat);
9*06a1af2fSStefano Zampini 
105e3038f0Sstefano_zampini PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat,MatType,MatReuse,Mat*);
11c8883902SJed Brown 
12d8588912SDave May /* private functions */
138188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N)
14d8588912SDave May {
15d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
168188e55aSJed Brown   PetscInt       i,j;
17d8588912SDave May   PetscErrorCode ierr;
18d8588912SDave May 
19d8588912SDave May   PetscFunctionBegin;
208188e55aSJed Brown   *m = *n = *M = *N = 0;
218188e55aSJed Brown   for (i=0; i<bA->nr; i++) {  /* rows */
228188e55aSJed Brown     PetscInt sm,sM;
238188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr);
248188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr);
258188e55aSJed Brown     *m  += sm;
268188e55aSJed Brown     *M  += sM;
27d8588912SDave May   }
288188e55aSJed Brown   for (j=0; j<bA->nc; j++) {  /* cols */
298188e55aSJed Brown     PetscInt sn,sN;
308188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr);
318188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr);
328188e55aSJed Brown     *n  += sn;
338188e55aSJed Brown     *N  += sN;
34d8588912SDave May   }
35d8588912SDave May   PetscFunctionReturn(0);
36d8588912SDave May }
37d8588912SDave May 
38d8588912SDave May /* operations */
39207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
40d8588912SDave May {
41d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
42207556f9SJed Brown   Vec            *bx = bA->right,*by = bA->left;
43207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
44d8588912SDave May   PetscErrorCode ierr;
45d8588912SDave May 
46d8588912SDave May   PetscFunctionBegin;
47207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
48207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
49207556f9SJed Brown   for (i=0; i<nr; i++) {
50d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
51207556f9SJed Brown     for (j=0; j<nc; j++) {
52207556f9SJed Brown       if (!bA->m[i][j]) continue;
53d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
54d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
55d8588912SDave May     }
56d8588912SDave May   }
57207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
58207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
59d8588912SDave May   PetscFunctionReturn(0);
60d8588912SDave May }
61d8588912SDave May 
629194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z)
639194d70fSJed Brown {
649194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
659194d70fSJed Brown   Vec            *bx = bA->right,*bz = bA->left;
669194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
679194d70fSJed Brown   PetscErrorCode ierr;
689194d70fSJed Brown 
699194d70fSJed Brown   PetscFunctionBegin;
709194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
719194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
729194d70fSJed Brown   for (i=0; i<nr; i++) {
739194d70fSJed Brown     if (y != z) {
749194d70fSJed Brown       Vec by;
759194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
769194d70fSJed Brown       ierr = VecCopy(by,bz[i]);CHKERRQ(ierr);
77336d21e7SJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
789194d70fSJed Brown     }
799194d70fSJed Brown     for (j=0; j<nc; j++) {
809194d70fSJed Brown       if (!bA->m[i][j]) continue;
819194d70fSJed Brown       /* y[i] <- y[i] + A[i][j] * x[j] */
829194d70fSJed Brown       ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr);
839194d70fSJed Brown     }
849194d70fSJed Brown   }
859194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
869194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
879194d70fSJed Brown   PetscFunctionReturn(0);
889194d70fSJed Brown }
899194d70fSJed Brown 
90207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
91d8588912SDave May {
92d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
93207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
94207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
95d8588912SDave May   PetscErrorCode ierr;
96d8588912SDave May 
97d8588912SDave May   PetscFunctionBegin;
98609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
99609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
100207556f9SJed Brown   for (j=0; j<nc; j++) {
101609e31cbSJed Brown     ierr = VecZeroEntries(by[j]);CHKERRQ(ierr);
102609e31cbSJed Brown     for (i=0; i<nr; i++) {
1036c75ac25SJed Brown       if (!bA->m[i][j]) continue;
104609e31cbSJed Brown       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
105609e31cbSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr);
106d8588912SDave May     }
107d8588912SDave May   }
108609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
109609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
110d8588912SDave May   PetscFunctionReturn(0);
111d8588912SDave May }
112d8588912SDave May 
1139194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
1149194d70fSJed Brown {
1159194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
1169194d70fSJed Brown   Vec            *bx = bA->left,*bz = bA->right;
1179194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
1189194d70fSJed Brown   PetscErrorCode ierr;
1199194d70fSJed Brown 
1209194d70fSJed Brown   PetscFunctionBegin;
1219194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1229194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1239194d70fSJed Brown   for (j=0; j<nc; j++) {
1249194d70fSJed Brown     if (y != z) {
1259194d70fSJed Brown       Vec by;
1269194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1279194d70fSJed Brown       ierr = VecCopy(by,bz[j]);CHKERRQ(ierr);
1289194d70fSJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1299194d70fSJed Brown     }
1309194d70fSJed Brown     for (i=0; i<nr; i++) {
1316c75ac25SJed Brown       if (!bA->m[i][j]) continue;
1329194d70fSJed Brown       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
1339194d70fSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr);
1349194d70fSJed Brown     }
1359194d70fSJed Brown   }
1369194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1379194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1389194d70fSJed Brown   PetscFunctionReturn(0);
1399194d70fSJed Brown }
1409194d70fSJed Brown 
141f8170845SAlex Fikl static PetscErrorCode MatTranspose_Nest(Mat A,MatReuse reuse,Mat *B)
142f8170845SAlex Fikl {
143f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data, *bC;
144f8170845SAlex Fikl   Mat            C;
145f8170845SAlex Fikl   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
146f8170845SAlex Fikl   PetscErrorCode ierr;
147f8170845SAlex Fikl 
148f8170845SAlex Fikl   PetscFunctionBegin;
149cf37664fSBarry Smith   if (reuse == MAT_INPLACE_MATRIX && nr != nc) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Square nested matrix only for in-place");
150f8170845SAlex Fikl 
151cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
152f8170845SAlex Fikl     Mat *subs;
153f8170845SAlex Fikl     IS  *is_row,*is_col;
154f8170845SAlex Fikl 
155f8170845SAlex Fikl     ierr = PetscCalloc1(nr * nc,&subs);CHKERRQ(ierr);
156f8170845SAlex Fikl     ierr = PetscMalloc2(nr,&is_row,nc,&is_col);CHKERRQ(ierr);
157f8170845SAlex Fikl     ierr = MatNestGetISs(A,is_row,is_col);CHKERRQ(ierr);
158cf37664fSBarry Smith     if (reuse == MAT_INPLACE_MATRIX) {
159ddeb9bd8SAlex Fikl       for (i=0; i<nr; i++) {
160ddeb9bd8SAlex Fikl         for (j=0; j<nc; j++) {
161ddeb9bd8SAlex Fikl           subs[i + nr * j] = bA->m[i][j];
162ddeb9bd8SAlex Fikl         }
163ddeb9bd8SAlex Fikl       }
164ddeb9bd8SAlex Fikl     }
165ddeb9bd8SAlex Fikl 
166f8170845SAlex Fikl     ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nc,is_col,nr,is_row,subs,&C);CHKERRQ(ierr);
167f8170845SAlex Fikl     ierr = PetscFree(subs);CHKERRQ(ierr);
1683d994f23SBarry Smith     ierr = PetscFree2(is_row,is_col);CHKERRQ(ierr);
169f8170845SAlex Fikl   } else {
170f8170845SAlex Fikl     C = *B;
171f8170845SAlex Fikl   }
172f8170845SAlex Fikl 
173f8170845SAlex Fikl   bC = (Mat_Nest*)C->data;
174f8170845SAlex Fikl   for (i=0; i<nr; i++) {
175f8170845SAlex Fikl     for (j=0; j<nc; j++) {
176f8170845SAlex Fikl       if (bA->m[i][j]) {
177f8170845SAlex Fikl         ierr = MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i]));CHKERRQ(ierr);
178f8170845SAlex Fikl       } else {
179f8170845SAlex Fikl         bC->m[j][i] = NULL;
180f8170845SAlex Fikl       }
181f8170845SAlex Fikl     }
182f8170845SAlex Fikl   }
183f8170845SAlex Fikl 
184cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
185f8170845SAlex Fikl     *B = C;
186f8170845SAlex Fikl   } else {
187f8170845SAlex Fikl     ierr = MatHeaderMerge(A, &C);CHKERRQ(ierr);
188f8170845SAlex Fikl   }
189f8170845SAlex Fikl   PetscFunctionReturn(0);
190f8170845SAlex Fikl }
191f8170845SAlex Fikl 
192e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
193e2d7f03fSJed Brown {
194e2d7f03fSJed Brown   PetscErrorCode ierr;
195e2d7f03fSJed Brown   IS             *lst = *list;
196e2d7f03fSJed Brown   PetscInt       i;
197e2d7f03fSJed Brown 
198e2d7f03fSJed Brown   PetscFunctionBegin;
199e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
2006bf464f9SBarry Smith   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
201e2d7f03fSJed Brown   ierr  = PetscFree(lst);CHKERRQ(ierr);
2020298fd71SBarry Smith   *list = NULL;
203e2d7f03fSJed Brown   PetscFunctionReturn(0);
204e2d7f03fSJed Brown }
205e2d7f03fSJed Brown 
206*06a1af2fSStefano Zampini static PetscErrorCode MatReset_Nest(Mat A)
207d8588912SDave May {
208d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
209d8588912SDave May   PetscInt       i,j;
210d8588912SDave May   PetscErrorCode ierr;
211d8588912SDave May 
212d8588912SDave May   PetscFunctionBegin;
213d8588912SDave May   /* release the matrices and the place holders */
214e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
215e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
216e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
217e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
218d8588912SDave May 
219d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
220d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
221*06a1af2fSStefano Zampini   ierr = PetscFree(vs->nnzstate);CHKERRQ(ierr);
222d8588912SDave May 
223207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
224207556f9SJed Brown 
225d8588912SDave May   /* release the matrices and the place holders */
226d8588912SDave May   if (vs->m) {
227d8588912SDave May     for (i=0; i<vs->nr; i++) {
228d8588912SDave May       for (j=0; j<vs->nc; j++) {
2296bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
230d8588912SDave May       }
231d8588912SDave May       ierr = PetscFree(vs->m[i]);CHKERRQ(ierr);
232d8588912SDave May     }
233d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
234d8588912SDave May   }
235*06a1af2fSStefano Zampini 
236*06a1af2fSStefano Zampini   /* restore defaults */
237*06a1af2fSStefano Zampini   vs->nr = 0;
238*06a1af2fSStefano Zampini   vs->nc = 0;
239*06a1af2fSStefano Zampini   vs->splitassembly = PETSC_FALSE;
240*06a1af2fSStefano Zampini   PetscFunctionReturn(0);
241*06a1af2fSStefano Zampini }
242*06a1af2fSStefano Zampini 
243*06a1af2fSStefano Zampini static PetscErrorCode MatDestroy_Nest(Mat A)
244*06a1af2fSStefano Zampini {
245*06a1af2fSStefano Zampini   PetscErrorCode ierr;
246*06a1af2fSStefano Zampini 
247*06a1af2fSStefano Zampini   ierr = MatReset_Nest(A);CHKERRQ(ierr);
248bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
249d8588912SDave May 
250bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);CHKERRQ(ierr);
251bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);CHKERRQ(ierr);
252bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);CHKERRQ(ierr);
253bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);CHKERRQ(ierr);
254bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);CHKERRQ(ierr);
255bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);CHKERRQ(ierr);
256bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);CHKERRQ(ierr);
257bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);CHKERRQ(ierr);
2580899c546SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_mpiaij_C",0);CHKERRQ(ierr);
2590899c546SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_seqaij_C",0);CHKERRQ(ierr);
2605e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",0);CHKERRQ(ierr);
2615e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",0);CHKERRQ(ierr);
262d8588912SDave May   PetscFunctionReturn(0);
263d8588912SDave May }
264d8588912SDave May 
265207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
266d8588912SDave May {
267d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
268d8588912SDave May   PetscInt       i,j;
269d8588912SDave May   PetscErrorCode ierr;
270*06a1af2fSStefano Zampini   PetscBool      nnzstate = PETSC_FALSE;
271d8588912SDave May 
272d8588912SDave May   PetscFunctionBegin;
273d8588912SDave May   for (i=0; i<vs->nr; i++) {
274d8588912SDave May     for (j=0; j<vs->nc; j++) {
275*06a1af2fSStefano Zampini       PetscObjectState subnnzstate = 0;
276e7c19651SJed Brown       if (vs->m[i][j]) {
277e7c19651SJed Brown         ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr);
278e7c19651SJed Brown         if (!vs->splitassembly) {
279e7c19651SJed Brown           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
280e7c19651SJed 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
281e7c19651SJed Brown            * already performing an assembly, but the result would by more complicated and appears to offer less
282e7c19651SJed Brown            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
283e7c19651SJed Brown            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
284e7c19651SJed Brown            */
285e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
286*06a1af2fSStefano Zampini           ierr = MatGetNonzeroState(vs->m[i][j],&subnnzstate);CHKERRQ(ierr);
287e7c19651SJed Brown         }
288e7c19651SJed Brown       }
289*06a1af2fSStefano Zampini       nnzstate = (PetscBool)(nnzstate || vs->nnzstate[i*vs->nc+j] != subnnzstate);
290*06a1af2fSStefano Zampini       vs->nnzstate[i*vs->nc+j] = subnnzstate;
291d8588912SDave May     }
292d8588912SDave May   }
293*06a1af2fSStefano Zampini   if (nnzstate) A->nonzerostate++;
294d8588912SDave May   PetscFunctionReturn(0);
295d8588912SDave May }
296d8588912SDave May 
297207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
298d8588912SDave May {
299d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
300d8588912SDave May   PetscInt       i,j;
301d8588912SDave May   PetscErrorCode ierr;
302d8588912SDave May 
303d8588912SDave May   PetscFunctionBegin;
304d8588912SDave May   for (i=0; i<vs->nr; i++) {
305d8588912SDave May     for (j=0; j<vs->nc; j++) {
306e7c19651SJed Brown       if (vs->m[i][j]) {
307e7c19651SJed Brown         if (vs->splitassembly) {
308e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
309e7c19651SJed Brown         }
310e7c19651SJed Brown       }
311d8588912SDave May     }
312d8588912SDave May   }
313d8588912SDave May   PetscFunctionReturn(0);
314d8588912SDave May }
315d8588912SDave May 
316f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
317d8588912SDave May {
318207556f9SJed Brown   PetscErrorCode ierr;
319f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
320f349c1fdSJed Brown   PetscInt       j;
321f349c1fdSJed Brown   Mat            sub;
322d8588912SDave May 
323d8588912SDave May   PetscFunctionBegin;
3240298fd71SBarry Smith   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
325f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
3264994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
327f349c1fdSJed Brown   *B = sub;
328f349c1fdSJed Brown   PetscFunctionReturn(0);
329d8588912SDave May }
330d8588912SDave May 
331f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
332f349c1fdSJed Brown {
333207556f9SJed Brown   PetscErrorCode ierr;
334f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
335f349c1fdSJed Brown   PetscInt       i;
336f349c1fdSJed Brown   Mat            sub;
337f349c1fdSJed Brown 
338f349c1fdSJed Brown   PetscFunctionBegin;
3390298fd71SBarry Smith   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
340f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
3414994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
342f349c1fdSJed Brown   *B = sub;
343f349c1fdSJed Brown   PetscFunctionReturn(0);
344d8588912SDave May }
345d8588912SDave May 
346f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
347f349c1fdSJed Brown {
348f349c1fdSJed Brown   PetscErrorCode ierr;
349f349c1fdSJed Brown   PetscInt       i;
350f349c1fdSJed Brown   PetscBool      flg;
351f349c1fdSJed Brown 
352f349c1fdSJed Brown   PetscFunctionBegin;
353f349c1fdSJed Brown   PetscValidPointer(list,3);
354f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
355f349c1fdSJed Brown   PetscValidIntPointer(found,5);
356f349c1fdSJed Brown   *found = -1;
357f349c1fdSJed Brown   for (i=0; i<n; i++) {
358207556f9SJed Brown     if (!list[i]) continue;
359f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
360f349c1fdSJed Brown     if (flg) {
361f349c1fdSJed Brown       *found = i;
362f349c1fdSJed Brown       PetscFunctionReturn(0);
363f349c1fdSJed Brown     }
364f349c1fdSJed Brown   }
365ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set");
366f349c1fdSJed Brown   PetscFunctionReturn(0);
367f349c1fdSJed Brown }
368f349c1fdSJed Brown 
3698188e55aSJed Brown /* Get a block row as a new MatNest */
3708188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
3718188e55aSJed Brown {
3728188e55aSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3738188e55aSJed Brown   char           keyname[256];
3748188e55aSJed Brown   PetscErrorCode ierr;
3758188e55aSJed Brown 
3768188e55aSJed Brown   PetscFunctionBegin;
3770298fd71SBarry Smith   *B   = NULL;
3788caf3d72SBarry Smith   ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr);
3798188e55aSJed Brown   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
3808188e55aSJed Brown   if (*B) PetscFunctionReturn(0);
3818188e55aSJed Brown 
382ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
38326fbe8dcSKarl Rupp 
3848188e55aSJed Brown   (*B)->assembled = A->assembled;
38526fbe8dcSKarl Rupp 
3868188e55aSJed Brown   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
3878188e55aSJed Brown   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
3888188e55aSJed Brown   PetscFunctionReturn(0);
3898188e55aSJed Brown }
3908188e55aSJed Brown 
391f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
392f349c1fdSJed Brown {
393f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3948188e55aSJed Brown   PetscErrorCode ierr;
3956b3a5b13SJed Brown   PetscInt       row,col;
396e072481dSJed Brown   PetscBool      same,isFullCol,isFullColGlobal;
397f349c1fdSJed Brown 
398f349c1fdSJed Brown   PetscFunctionBegin;
3998188e55aSJed Brown   /* Check if full column space. This is a hack */
4008188e55aSJed Brown   isFullCol = PETSC_FALSE;
401251f4c67SDmitry Karpeev   ierr      = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
4028188e55aSJed Brown   if (same) {
40377019fcaSJed Brown     PetscInt n,first,step,i,an,am,afirst,astep;
4048188e55aSJed Brown     ierr      = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
4058188e55aSJed Brown     ierr      = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
40677019fcaSJed Brown     isFullCol = PETSC_TRUE;
40705ce4453SJed Brown     for (i=0,an=A->cmap->rstart; i<vs->nc; i++) {
408*06a1af2fSStefano Zampini       ierr = PetscObjectTypeCompare((PetscObject)is->col[i],ISSTRIDE,&same);CHKERRQ(ierr);
40977019fcaSJed Brown       ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr);
410*06a1af2fSStefano Zampini       if (same) {
411*06a1af2fSStefano Zampini         ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr);
41277019fcaSJed Brown         if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
413*06a1af2fSStefano Zampini       } else isFullCol = PETSC_FALSE;
41477019fcaSJed Brown       an += am;
41577019fcaSJed Brown     }
41605ce4453SJed Brown     if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE;
4178188e55aSJed Brown   }
418b2566f29SBarry Smith   ierr = MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr);
4198188e55aSJed Brown 
420427230ceSLisandro Dalcin   if (isFullColGlobal && vs->nc > 1) {
4218188e55aSJed Brown     PetscInt row;
4228188e55aSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
4238188e55aSJed Brown     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
4248188e55aSJed Brown   } else {
425f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
426f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
427b6480e04SStefano Zampini     if (!vs->m[row][col]) {
428b6480e04SStefano Zampini       PetscInt lr,lc;
429b6480e04SStefano Zampini 
430b6480e04SStefano Zampini       ierr = MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);CHKERRQ(ierr);
431b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.row[row],&lr);CHKERRQ(ierr);
432b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.col[col],&lc);CHKERRQ(ierr);
433b6480e04SStefano Zampini       ierr = MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
434b6480e04SStefano Zampini       ierr = MatSetUp(vs->m[row][col]);CHKERRQ(ierr);
435b6480e04SStefano Zampini       ierr = MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
436b6480e04SStefano Zampini       ierr = MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
437b6480e04SStefano Zampini     }
438f349c1fdSJed Brown     *B = vs->m[row][col];
4398188e55aSJed Brown   }
440f349c1fdSJed Brown   PetscFunctionReturn(0);
441f349c1fdSJed Brown }
442f349c1fdSJed Brown 
443*06a1af2fSStefano Zampini /*
444*06a1af2fSStefano Zampini    TODO: This does not actually returns a submatrix we can modify
445*06a1af2fSStefano Zampini */
4467dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
447f349c1fdSJed Brown {
448f349c1fdSJed Brown   PetscErrorCode ierr;
449f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
450f349c1fdSJed Brown   Mat            sub;
451f349c1fdSJed Brown 
452f349c1fdSJed Brown   PetscFunctionBegin;
453f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
454f349c1fdSJed Brown   switch (reuse) {
455f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
4567874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
457f349c1fdSJed Brown     *B = sub;
458f349c1fdSJed Brown     break;
459f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
460ce94432eSBarry Smith     if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
461f349c1fdSJed Brown     break;
462f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
463f349c1fdSJed Brown     break;
464511c6705SHong Zhang   case MAT_INPLACE_MATRIX:       /* Nothing to do */
465511c6705SHong Zhang     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet");
466511c6705SHong Zhang     break;
467f349c1fdSJed Brown   }
468f349c1fdSJed Brown   PetscFunctionReturn(0);
469f349c1fdSJed Brown }
470f349c1fdSJed Brown 
471f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
472f349c1fdSJed Brown {
473f349c1fdSJed Brown   PetscErrorCode ierr;
474f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
475f349c1fdSJed Brown   Mat            sub;
476f349c1fdSJed Brown 
477f349c1fdSJed Brown   PetscFunctionBegin;
478f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
479f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
480f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
481f349c1fdSJed Brown   *B = sub;
482d8588912SDave May   PetscFunctionReturn(0);
483d8588912SDave May }
484d8588912SDave May 
485207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
486d8588912SDave May {
487d8588912SDave May   PetscErrorCode ierr;
488f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
489f349c1fdSJed Brown   Mat            sub;
490d8588912SDave May 
491d8588912SDave May   PetscFunctionBegin;
492f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
493ce94432eSBarry Smith   if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
494f349c1fdSJed Brown   if (sub) {
495ce94432eSBarry Smith     if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
4966bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
497d8588912SDave May   }
498d8588912SDave May   PetscFunctionReturn(0);
499d8588912SDave May }
500d8588912SDave May 
5017874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
5027874fa86SDave May {
5037874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
5047874fa86SDave May   PetscInt       i;
5057874fa86SDave May   PetscErrorCode ierr;
5067874fa86SDave May 
5077874fa86SDave May   PetscFunctionBegin;
5087874fa86SDave May   for (i=0; i<bA->nr; i++) {
509429bac76SJed Brown     Vec bv;
510429bac76SJed Brown     ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
5117874fa86SDave May     if (bA->m[i][i]) {
512429bac76SJed Brown       ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr);
5137874fa86SDave May     } else {
5145159a857SMatthew G. Knepley       ierr = VecSet(bv,0.0);CHKERRQ(ierr);
5157874fa86SDave May     }
516429bac76SJed Brown     ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
5177874fa86SDave May   }
5187874fa86SDave May   PetscFunctionReturn(0);
5197874fa86SDave May }
5207874fa86SDave May 
5217874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
5227874fa86SDave May {
5237874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
524429bac76SJed Brown   Vec            bl,*br;
5257874fa86SDave May   PetscInt       i,j;
5267874fa86SDave May   PetscErrorCode ierr;
5277874fa86SDave May 
5287874fa86SDave May   PetscFunctionBegin;
5293f800ebeSJed Brown   ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr);
5302e6472ebSElliott Sales de Andrade   if (r) {
531429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5322e6472ebSElliott Sales de Andrade   }
5332e6472ebSElliott Sales de Andrade   bl = NULL;
5347874fa86SDave May   for (i=0; i<bA->nr; i++) {
5352e6472ebSElliott Sales de Andrade     if (l) {
536429bac76SJed Brown       ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5372e6472ebSElliott Sales de Andrade     }
5387874fa86SDave May     for (j=0; j<bA->nc; j++) {
5397874fa86SDave May       if (bA->m[i][j]) {
540429bac76SJed Brown         ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr);
5417874fa86SDave May       }
5427874fa86SDave May     }
5432e6472ebSElliott Sales de Andrade     if (l) {
544a061e289SJed Brown       ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5457874fa86SDave May     }
5462e6472ebSElliott Sales de Andrade   }
5472e6472ebSElliott Sales de Andrade   if (r) {
548429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5492e6472ebSElliott Sales de Andrade   }
550429bac76SJed Brown   ierr = PetscFree(br);CHKERRQ(ierr);
5517874fa86SDave May   PetscFunctionReturn(0);
5527874fa86SDave May }
5537874fa86SDave May 
554a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
555a061e289SJed Brown {
556a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
557a061e289SJed Brown   PetscInt       i,j;
558a061e289SJed Brown   PetscErrorCode ierr;
559a061e289SJed Brown 
560a061e289SJed Brown   PetscFunctionBegin;
561a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
562a061e289SJed Brown     for (j=0; j<bA->nc; j++) {
563a061e289SJed Brown       if (bA->m[i][j]) {
564a061e289SJed Brown         ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr);
565a061e289SJed Brown       }
566a061e289SJed Brown     }
567a061e289SJed Brown   }
568a061e289SJed Brown   PetscFunctionReturn(0);
569a061e289SJed Brown }
570a061e289SJed Brown 
571a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
572a061e289SJed Brown {
573a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
574a061e289SJed Brown   PetscInt       i;
575a061e289SJed Brown   PetscErrorCode ierr;
576*06a1af2fSStefano Zampini   PetscBool      nnzstate = PETSC_FALSE;
577a061e289SJed Brown 
578a061e289SJed Brown   PetscFunctionBegin;
579a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
580*06a1af2fSStefano Zampini     PetscObjectState subnnzstate = 0;
581ce94432eSBarry 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);
582a061e289SJed Brown     ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr);
583*06a1af2fSStefano Zampini     ierr = MatGetNonzeroState(bA->m[i][i],&subnnzstate);CHKERRQ(ierr);
584*06a1af2fSStefano Zampini     nnzstate = (PetscBool)(nnzstate || bA->nnzstate[i*bA->nc+i] != subnnzstate);
585*06a1af2fSStefano Zampini     bA->nnzstate[i*bA->nc+i] = subnnzstate;
586a061e289SJed Brown   }
587*06a1af2fSStefano Zampini   if (nnzstate) A->nonzerostate++;
588a061e289SJed Brown   PetscFunctionReturn(0);
589a061e289SJed Brown }
590a061e289SJed Brown 
59113135bc6SAlex Fikl static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is)
59213135bc6SAlex Fikl {
59313135bc6SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
59413135bc6SAlex Fikl   PetscInt       i;
59513135bc6SAlex Fikl   PetscErrorCode ierr;
596*06a1af2fSStefano Zampini   PetscBool      nnzstate = PETSC_FALSE;
59713135bc6SAlex Fikl 
59813135bc6SAlex Fikl   PetscFunctionBegin;
59913135bc6SAlex Fikl   for (i=0; i<bA->nr; i++) {
600*06a1af2fSStefano Zampini     PetscObjectState subnnzstate = 0;
60113135bc6SAlex Fikl     Vec              bv;
60213135bc6SAlex Fikl     ierr = VecGetSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
60313135bc6SAlex Fikl     if (bA->m[i][i]) {
60413135bc6SAlex Fikl       ierr = MatDiagonalSet(bA->m[i][i],bv,is);CHKERRQ(ierr);
605*06a1af2fSStefano Zampini       ierr = MatGetNonzeroState(bA->m[i][i],&subnnzstate);CHKERRQ(ierr);
60613135bc6SAlex Fikl     }
60713135bc6SAlex Fikl     ierr = VecRestoreSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
608*06a1af2fSStefano Zampini     nnzstate = (PetscBool)(nnzstate || bA->nnzstate[i*bA->nc+i] != subnnzstate);
609*06a1af2fSStefano Zampini     bA->nnzstate[i*bA->nc+i] = subnnzstate;
61013135bc6SAlex Fikl   }
611*06a1af2fSStefano Zampini   if (nnzstate) A->nonzerostate++;
61213135bc6SAlex Fikl   PetscFunctionReturn(0);
61313135bc6SAlex Fikl }
61413135bc6SAlex Fikl 
615f8170845SAlex Fikl static PetscErrorCode MatSetRandom_Nest(Mat A,PetscRandom rctx)
616f8170845SAlex Fikl {
617f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
618f8170845SAlex Fikl   PetscInt       i,j;
619f8170845SAlex Fikl   PetscErrorCode ierr;
620f8170845SAlex Fikl 
621f8170845SAlex Fikl   PetscFunctionBegin;
622f8170845SAlex Fikl   for (i=0; i<bA->nr; i++) {
623f8170845SAlex Fikl     for (j=0; j<bA->nc; j++) {
624f8170845SAlex Fikl       if (bA->m[i][j]) {
625f8170845SAlex Fikl         ierr = MatSetRandom(bA->m[i][j],rctx);CHKERRQ(ierr);
626f8170845SAlex Fikl       }
627f8170845SAlex Fikl     }
628f8170845SAlex Fikl   }
629f8170845SAlex Fikl   PetscFunctionReturn(0);
630f8170845SAlex Fikl }
631f8170845SAlex Fikl 
6322a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left)
633d8588912SDave May {
634d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
635d8588912SDave May   Vec            *L,*R;
636d8588912SDave May   MPI_Comm       comm;
637d8588912SDave May   PetscInt       i,j;
638d8588912SDave May   PetscErrorCode ierr;
639d8588912SDave May 
640d8588912SDave May   PetscFunctionBegin;
641ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
642d8588912SDave May   if (right) {
643d8588912SDave May     /* allocate R */
644854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nc, &R);CHKERRQ(ierr);
645d8588912SDave May     /* Create the right vectors */
646d8588912SDave May     for (j=0; j<bA->nc; j++) {
647d8588912SDave May       for (i=0; i<bA->nr; i++) {
648d8588912SDave May         if (bA->m[i][j]) {
6492a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr);
650d8588912SDave May           break;
651d8588912SDave May         }
652d8588912SDave May       }
6536c4ed002SBarry Smith       if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
654d8588912SDave May     }
655f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
656d8588912SDave May     /* hand back control to the nest vector */
657d8588912SDave May     for (j=0; j<bA->nc; j++) {
6586bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
659d8588912SDave May     }
660d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
661d8588912SDave May   }
662d8588912SDave May 
663d8588912SDave May   if (left) {
664d8588912SDave May     /* allocate L */
665854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nr, &L);CHKERRQ(ierr);
666d8588912SDave May     /* Create the left vectors */
667d8588912SDave May     for (i=0; i<bA->nr; i++) {
668d8588912SDave May       for (j=0; j<bA->nc; j++) {
669d8588912SDave May         if (bA->m[i][j]) {
6702a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr);
671d8588912SDave May           break;
672d8588912SDave May         }
673d8588912SDave May       }
6746c4ed002SBarry Smith       if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
675d8588912SDave May     }
676d8588912SDave May 
677f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
678d8588912SDave May     for (i=0; i<bA->nr; i++) {
6796bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
680d8588912SDave May     }
681d8588912SDave May 
682d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
683d8588912SDave May   }
684d8588912SDave May   PetscFunctionReturn(0);
685d8588912SDave May }
686d8588912SDave May 
687207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
688d8588912SDave May {
689d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
69029e60adbSStefano Zampini   PetscBool      isascii,viewSub = PETSC_FALSE;
691d8588912SDave May   PetscInt       i,j;
692d8588912SDave May   PetscErrorCode ierr;
693d8588912SDave May 
694d8588912SDave May   PetscFunctionBegin;
695251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
696d8588912SDave May   if (isascii) {
697d8588912SDave May 
69829e60adbSStefano Zampini     ierr = PetscOptionsGetBool(((PetscObject)A)->options,((PetscObject)A)->prefix,"-mat_view_nest_sub",&viewSub,NULL);CHKERRQ(ierr);
699d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Matrix object: \n");CHKERRQ(ierr);
700d86155a6SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
701d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "type=nest, rows=%D, cols=%D \n",bA->nr,bA->nc);CHKERRQ(ierr);
702d8588912SDave May 
703d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");CHKERRQ(ierr);
704d8588912SDave May     for (i=0; i<bA->nr; i++) {
705d8588912SDave May       for (j=0; j<bA->nc; j++) {
70619fd82e9SBarry Smith         MatType   type;
707270f95d7SJed Brown         char      name[256] = "",prefix[256] = "";
708d8588912SDave May         PetscInt  NR,NC;
709d8588912SDave May         PetscBool isNest = PETSC_FALSE;
710d8588912SDave May 
711d8588912SDave May         if (!bA->m[i][j]) {
712d86155a6SBarry Smith           CHKERRQ(ierr);PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);CHKERRQ(ierr);
713d8588912SDave May           continue;
714d8588912SDave May         }
715d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
716d8588912SDave May         ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr);
7178caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);}
7188caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);}
719251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
720d8588912SDave May 
721270f95d7SJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr);
722d8588912SDave May 
72329e60adbSStefano Zampini         if (isNest || viewSub) {
724270f95d7SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* push1 */
725d8588912SDave May           ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr);
726270f95d7SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop1 */
727d8588912SDave May         }
728d8588912SDave May       }
729d8588912SDave May     }
730d86155a6SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop0 */
731d8588912SDave May   }
732d8588912SDave May   PetscFunctionReturn(0);
733d8588912SDave May }
734d8588912SDave May 
735207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
736d8588912SDave May {
737d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
738d8588912SDave May   PetscInt       i,j;
739d8588912SDave May   PetscErrorCode ierr;
740d8588912SDave May 
741d8588912SDave May   PetscFunctionBegin;
742d8588912SDave May   for (i=0; i<bA->nr; i++) {
743d8588912SDave May     for (j=0; j<bA->nc; j++) {
744d8588912SDave May       if (!bA->m[i][j]) continue;
745d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
746d8588912SDave May     }
747d8588912SDave May   }
748d8588912SDave May   PetscFunctionReturn(0);
749d8588912SDave May }
750d8588912SDave May 
751c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str)
752c222c20dSDavid Ham {
753c222c20dSDavid Ham   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
754c222c20dSDavid Ham   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
755c222c20dSDavid Ham   PetscErrorCode ierr;
756*06a1af2fSStefano Zampini   PetscBool      nnzstate = PETSC_FALSE;
757c222c20dSDavid Ham 
758c222c20dSDavid Ham   PetscFunctionBegin;
759c222c20dSDavid 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);
760c222c20dSDavid Ham   for (i=0; i<nr; i++) {
761c222c20dSDavid Ham     for (j=0; j<nc; j++) {
762*06a1af2fSStefano Zampini       PetscObjectState subnnzstate = 0;
76346a2b97cSJed Brown       if (bA->m[i][j] && bB->m[i][j]) {
764c222c20dSDavid Ham         ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr);
76546a2b97cSJed 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);
766*06a1af2fSStefano Zampini       ierr = MatGetNonzeroState(bB->m[i][j],&subnnzstate);CHKERRQ(ierr);
767*06a1af2fSStefano Zampini       nnzstate = (PetscBool)(nnzstate || bB->nnzstate[i*nc+j] != subnnzstate);
768*06a1af2fSStefano Zampini       bB->nnzstate[i*nc+j] = subnnzstate;
769c222c20dSDavid Ham     }
770c222c20dSDavid Ham   }
771*06a1af2fSStefano Zampini   if (nnzstate) B->nonzerostate++;
772c222c20dSDavid Ham   PetscFunctionReturn(0);
773c222c20dSDavid Ham }
774c222c20dSDavid Ham 
7756e76ffeaSPierre Jolivet static PetscErrorCode MatAXPY_Nest(Mat Y,PetscScalar a,Mat X,MatStructure str)
7766e76ffeaSPierre Jolivet {
7776e76ffeaSPierre Jolivet   Mat_Nest       *bY = (Mat_Nest*)Y->data,*bX = (Mat_Nest*)X->data;
7786e76ffeaSPierre Jolivet   PetscInt       i,j,nr = bY->nr,nc = bY->nc;
7796e76ffeaSPierre Jolivet   PetscErrorCode ierr;
780*06a1af2fSStefano Zampini   PetscBool      nnzstate = PETSC_FALSE;
7816e76ffeaSPierre Jolivet 
7826e76ffeaSPierre Jolivet   PetscFunctionBegin;
7836e76ffeaSPierre 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);
7846e76ffeaSPierre Jolivet   for (i=0; i<nr; i++) {
7856e76ffeaSPierre Jolivet     for (j=0; j<nc; j++) {
786*06a1af2fSStefano Zampini       PetscObjectState subnnzstate = 0;
7876e76ffeaSPierre Jolivet       if (bY->m[i][j] && bX->m[i][j]) {
7886e76ffeaSPierre Jolivet         ierr = MatAXPY(bY->m[i][j],a,bX->m[i][j],str);CHKERRQ(ierr);
789c066aebcSStefano Zampini       } else if (bX->m[i][j]) {
790c066aebcSStefano Zampini         Mat M;
791c066aebcSStefano Zampini 
792c066aebcSStefano Zampini         if (str != DIFFERENT_NONZERO_PATTERN) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Matrix block does not exist at %D,%D",i,j);
793c066aebcSStefano Zampini         ierr = MatDuplicate(bX->m[i][j],MAT_COPY_VALUES,&M);CHKERRQ(ierr);
794c066aebcSStefano Zampini         ierr = MatNestSetSubMat(Y,i,j,M);CHKERRQ(ierr);
795c066aebcSStefano Zampini         ierr = MatDestroy(&M);CHKERRQ(ierr);
796c066aebcSStefano Zampini       }
797*06a1af2fSStefano Zampini       ierr = MatGetNonzeroState(bY->m[i][j],&subnnzstate);CHKERRQ(ierr);
798*06a1af2fSStefano Zampini       nnzstate = (PetscBool)(nnzstate || bY->nnzstate[i*nc+j] != subnnzstate);
799*06a1af2fSStefano Zampini       bY->nnzstate[i*nc+j] = subnnzstate;
8006e76ffeaSPierre Jolivet     }
8016e76ffeaSPierre Jolivet   }
802*06a1af2fSStefano Zampini   if (nnzstate) Y->nonzerostate++;
8036e76ffeaSPierre Jolivet   PetscFunctionReturn(0);
8046e76ffeaSPierre Jolivet }
8056e76ffeaSPierre Jolivet 
806207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
807d8588912SDave May {
808d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
809841e96a3SJed Brown   Mat            *b;
810841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
811d8588912SDave May   PetscErrorCode ierr;
812d8588912SDave May 
813d8588912SDave May   PetscFunctionBegin;
814785e854fSJed Brown   ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr);
815841e96a3SJed Brown   for (i=0; i<nr; i++) {
816841e96a3SJed Brown     for (j=0; j<nc; j++) {
817841e96a3SJed Brown       if (bA->m[i][j]) {
818841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
819841e96a3SJed Brown       } else {
8200298fd71SBarry Smith         b[i*nc+j] = NULL;
821d8588912SDave May       }
822d8588912SDave May     }
823d8588912SDave May   }
824ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
825841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
826841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
8276bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
828d8588912SDave May   }
829d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
830d8588912SDave May 
831841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
832841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
833d8588912SDave May   PetscFunctionReturn(0);
834d8588912SDave May }
835d8588912SDave May 
836d8588912SDave May /* nest api */
837d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
838d8588912SDave May {
839d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
8405fd66863SKarl Rupp 
841d8588912SDave May   PetscFunctionBegin;
842ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
843ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
844d8588912SDave May   *mat = bA->m[idxm][jdxm];
845d8588912SDave May   PetscFunctionReturn(0);
846d8588912SDave May }
847d8588912SDave May 
8489ba0d327SJed Brown /*@
849d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
850d8588912SDave May 
851d8588912SDave May  Not collective
852d8588912SDave May 
853d8588912SDave May  Input Parameters:
854629881c0SJed Brown +   A  - nest matrix
855d8588912SDave May .   idxm - index of the matrix within the nest matrix
856629881c0SJed Brown -   jdxm - index of the matrix within the nest matrix
857d8588912SDave May 
858d8588912SDave May  Output Parameter:
859d8588912SDave May .   sub - matrix at index idxm,jdxm within the nest matrix
860d8588912SDave May 
861d8588912SDave May  Level: developer
862d8588912SDave May 
863d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats()
864d8588912SDave May @*/
8657087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
866d8588912SDave May {
867699a902aSJed Brown   PetscErrorCode ierr;
868d8588912SDave May 
869d8588912SDave May   PetscFunctionBegin;
870699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
871d8588912SDave May   PetscFunctionReturn(0);
872d8588912SDave May }
873d8588912SDave May 
8740782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
8750782ca92SJed Brown {
8760782ca92SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
8770782ca92SJed Brown   PetscInt       m,n,M,N,mi,ni,Mi,Ni;
8780782ca92SJed Brown   PetscErrorCode ierr;
8790782ca92SJed Brown 
8800782ca92SJed Brown   PetscFunctionBegin;
881ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
882ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
8830782ca92SJed Brown   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
8840782ca92SJed Brown   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
8850782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
8860782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
8870782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
8880782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
889ce94432eSBarry 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);
890ce94432eSBarry 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);
89126fbe8dcSKarl Rupp 
892*06a1af2fSStefano Zampini   /* do not increase object state */
893*06a1af2fSStefano Zampini   if (mat == bA->m[idxm][jdxm]) PetscFunctionReturn(0);
894*06a1af2fSStefano Zampini 
8950782ca92SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
8960782ca92SJed Brown   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
8970782ca92SJed Brown   bA->m[idxm][jdxm] = mat;
898*06a1af2fSStefano Zampini   ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr);
899*06a1af2fSStefano Zampini   ierr = MatGetNonzeroState(mat,&bA->nnzstate[idxm*bA->nc+jdxm]);CHKERRQ(ierr);
900*06a1af2fSStefano Zampini   A->nonzerostate++;
9010782ca92SJed Brown   PetscFunctionReturn(0);
9020782ca92SJed Brown }
9030782ca92SJed Brown 
9049ba0d327SJed Brown /*@
9050782ca92SJed Brown  MatNestSetSubMat - Set a single submatrix in the nest matrix.
9060782ca92SJed Brown 
9070782ca92SJed Brown  Logically collective on the submatrix communicator
9080782ca92SJed Brown 
9090782ca92SJed Brown  Input Parameters:
9100782ca92SJed Brown +   A  - nest matrix
9110782ca92SJed Brown .   idxm - index of the matrix within the nest matrix
9120782ca92SJed Brown .   jdxm - index of the matrix within the nest matrix
9130782ca92SJed Brown -   sub - matrix at index idxm,jdxm within the nest matrix
9140782ca92SJed Brown 
9150782ca92SJed Brown  Notes:
9160782ca92SJed Brown  The new submatrix must have the same size and communicator as that block of the nest.
9170782ca92SJed Brown 
9180782ca92SJed Brown  This increments the reference count of the submatrix.
9190782ca92SJed Brown 
9200782ca92SJed Brown  Level: developer
9210782ca92SJed Brown 
922d5dfb694SBarry Smith .seealso: MatNestSetSubMats(), MatNestGetSubMats()
9230782ca92SJed Brown @*/
9240782ca92SJed Brown PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
9250782ca92SJed Brown {
9260782ca92SJed Brown   PetscErrorCode ierr;
9270782ca92SJed Brown 
9280782ca92SJed Brown   PetscFunctionBegin;
9290782ca92SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
9300782ca92SJed Brown   PetscFunctionReturn(0);
9310782ca92SJed Brown }
9320782ca92SJed Brown 
933d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
934d8588912SDave May {
935d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
9365fd66863SKarl Rupp 
937d8588912SDave May   PetscFunctionBegin;
93826fbe8dcSKarl Rupp   if (M)   *M   = bA->nr;
93926fbe8dcSKarl Rupp   if (N)   *N   = bA->nc;
94026fbe8dcSKarl Rupp   if (mat) *mat = bA->m;
941d8588912SDave May   PetscFunctionReturn(0);
942d8588912SDave May }
943d8588912SDave May 
944d8588912SDave May /*@C
945d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
946d8588912SDave May 
947d8588912SDave May  Not collective
948d8588912SDave May 
949d8588912SDave May  Input Parameters:
950629881c0SJed Brown .   A  - nest matrix
951d8588912SDave May 
952d8588912SDave May  Output Parameter:
953629881c0SJed Brown +   M - number of rows in the nest matrix
954d8588912SDave May .   N - number of cols in the nest matrix
955629881c0SJed Brown -   mat - 2d array of matrices
956d8588912SDave May 
957d8588912SDave May  Notes:
958d8588912SDave May 
959d8588912SDave May  The user should not free the array mat.
960d8588912SDave May 
961351962e3SVincent Le Chenadec  In Fortran, this routine has a calling sequence
962351962e3SVincent Le Chenadec $   call MatNestGetSubMats(A, M, N, mat, ierr)
963351962e3SVincent Le Chenadec  where the space allocated for the optional argument mat is assumed large enough (if provided).
964351962e3SVincent Le Chenadec 
965d8588912SDave May  Level: developer
966d8588912SDave May 
967d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat()
968d8588912SDave May @*/
9697087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
970d8588912SDave May {
971699a902aSJed Brown   PetscErrorCode ierr;
972d8588912SDave May 
973d8588912SDave May   PetscFunctionBegin;
974699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
975d8588912SDave May   PetscFunctionReturn(0);
976d8588912SDave May }
977d8588912SDave May 
9787087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
979d8588912SDave May {
980d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
981d8588912SDave May 
982d8588912SDave May   PetscFunctionBegin;
98326fbe8dcSKarl Rupp   if (M) *M = bA->nr;
98426fbe8dcSKarl Rupp   if (N) *N = bA->nc;
985d8588912SDave May   PetscFunctionReturn(0);
986d8588912SDave May }
987d8588912SDave May 
9889ba0d327SJed Brown /*@
989d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
990d8588912SDave May 
991d8588912SDave May  Not collective
992d8588912SDave May 
993d8588912SDave May  Input Parameters:
994d8588912SDave May .   A  - nest matrix
995d8588912SDave May 
996d8588912SDave May  Output Parameter:
997629881c0SJed Brown +   M - number of rows in the nested mat
998629881c0SJed Brown -   N - number of cols in the nested mat
999d8588912SDave May 
1000d8588912SDave May  Notes:
1001d8588912SDave May 
1002d8588912SDave May  Level: developer
1003d8588912SDave May 
1004d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats()
1005d8588912SDave May @*/
10067087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
1007d8588912SDave May {
1008699a902aSJed Brown   PetscErrorCode ierr;
1009d8588912SDave May 
1010d8588912SDave May   PetscFunctionBegin;
1011699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
1012d8588912SDave May   PetscFunctionReturn(0);
1013d8588912SDave May }
1014d8588912SDave May 
1015f7a08781SBarry Smith static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
1016900e7ff2SJed Brown {
1017900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
1018900e7ff2SJed Brown   PetscInt i;
1019900e7ff2SJed Brown 
1020900e7ff2SJed Brown   PetscFunctionBegin;
1021900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
1022900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
1023900e7ff2SJed Brown   PetscFunctionReturn(0);
1024900e7ff2SJed Brown }
1025900e7ff2SJed Brown 
10263a4d7b9aSSatish Balay /*@C
1027900e7ff2SJed Brown  MatNestGetISs - Returns the index sets partitioning the row and column spaces
1028900e7ff2SJed Brown 
1029900e7ff2SJed Brown  Not collective
1030900e7ff2SJed Brown 
1031900e7ff2SJed Brown  Input Parameters:
1032900e7ff2SJed Brown .   A  - nest matrix
1033900e7ff2SJed Brown 
1034900e7ff2SJed Brown  Output Parameter:
1035900e7ff2SJed Brown +   rows - array of row index sets
1036900e7ff2SJed Brown -   cols - array of column index sets
1037900e7ff2SJed Brown 
1038900e7ff2SJed Brown  Level: advanced
1039900e7ff2SJed Brown 
1040900e7ff2SJed Brown  Notes:
1041900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1042900e7ff2SJed Brown 
1043900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
1044900e7ff2SJed Brown @*/
1045900e7ff2SJed Brown PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
1046900e7ff2SJed Brown {
1047900e7ff2SJed Brown   PetscErrorCode ierr;
1048900e7ff2SJed Brown 
1049900e7ff2SJed Brown   PetscFunctionBegin;
1050900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1051900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1052900e7ff2SJed Brown   PetscFunctionReturn(0);
1053900e7ff2SJed Brown }
1054900e7ff2SJed Brown 
1055f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
1056900e7ff2SJed Brown {
1057900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
1058900e7ff2SJed Brown   PetscInt i;
1059900e7ff2SJed Brown 
1060900e7ff2SJed Brown   PetscFunctionBegin;
1061900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
1062900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
1063900e7ff2SJed Brown   PetscFunctionReturn(0);
1064900e7ff2SJed Brown }
1065900e7ff2SJed Brown 
1066900e7ff2SJed Brown /*@C
1067900e7ff2SJed Brown  MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces
1068900e7ff2SJed Brown 
1069900e7ff2SJed Brown  Not collective
1070900e7ff2SJed Brown 
1071900e7ff2SJed Brown  Input Parameters:
1072900e7ff2SJed Brown .   A  - nest matrix
1073900e7ff2SJed Brown 
1074900e7ff2SJed Brown  Output Parameter:
10750298fd71SBarry Smith +   rows - array of row index sets (or NULL to ignore)
10760298fd71SBarry Smith -   cols - array of column index sets (or NULL to ignore)
1077900e7ff2SJed Brown 
1078900e7ff2SJed Brown  Level: advanced
1079900e7ff2SJed Brown 
1080900e7ff2SJed Brown  Notes:
1081900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1082900e7ff2SJed Brown 
1083900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
1084900e7ff2SJed Brown @*/
1085900e7ff2SJed Brown PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
1086900e7ff2SJed Brown {
1087900e7ff2SJed Brown   PetscErrorCode ierr;
1088900e7ff2SJed Brown 
1089900e7ff2SJed Brown   PetscFunctionBegin;
1090900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1091900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1092900e7ff2SJed Brown   PetscFunctionReturn(0);
1093900e7ff2SJed Brown }
1094900e7ff2SJed Brown 
109519fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
1096207556f9SJed Brown {
1097207556f9SJed Brown   PetscErrorCode ierr;
1098207556f9SJed Brown   PetscBool      flg;
1099207556f9SJed Brown 
1100207556f9SJed Brown   PetscFunctionBegin;
1101207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
1102207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
11032a7a6963SBarry Smith   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
110412b53f24SSatish Balay   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
1105207556f9SJed Brown   PetscFunctionReturn(0);
1106207556f9SJed Brown }
1107207556f9SJed Brown 
1108207556f9SJed Brown /*@C
11092a7a6963SBarry Smith  MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs()
1110207556f9SJed Brown 
1111207556f9SJed Brown  Not collective
1112207556f9SJed Brown 
1113207556f9SJed Brown  Input Parameters:
1114207556f9SJed Brown +  A  - nest matrix
1115207556f9SJed Brown -  vtype - type to use for creating vectors
1116207556f9SJed Brown 
1117207556f9SJed Brown  Notes:
1118207556f9SJed Brown 
1119207556f9SJed Brown  Level: developer
1120207556f9SJed Brown 
11212a7a6963SBarry Smith .seealso: MatCreateVecs()
1122207556f9SJed Brown @*/
112319fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
1124207556f9SJed Brown {
1125207556f9SJed Brown   PetscErrorCode ierr;
1126207556f9SJed Brown 
1127207556f9SJed Brown   PetscFunctionBegin;
112819fd82e9SBarry Smith   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr);
1129207556f9SJed Brown   PetscFunctionReturn(0);
1130207556f9SJed Brown }
1131207556f9SJed Brown 
1132c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1133d8588912SDave May {
1134c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
1135c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
1136d8588912SDave May   PetscErrorCode ierr;
1137*06a1af2fSStefano Zampini   PetscBool      cong;
1138d8588912SDave May 
1139d8588912SDave May   PetscFunctionBegin;
1140*06a1af2fSStefano Zampini   ierr = MatReset_Nest(A);CHKERRQ(ierr);
1141*06a1af2fSStefano Zampini 
1142c8883902SJed Brown   s->nr = nr;
1143c8883902SJed Brown   s->nc = nc;
1144d8588912SDave May 
1145c8883902SJed Brown   /* Create space for submatrices */
1146854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->m);CHKERRQ(ierr);
1147c8883902SJed Brown   for (i=0; i<nr; i++) {
1148854ce69bSBarry Smith     ierr = PetscMalloc1(nc,&s->m[i]);CHKERRQ(ierr);
1149d8588912SDave May   }
1150c8883902SJed Brown   for (i=0; i<nr; i++) {
1151c8883902SJed Brown     for (j=0; j<nc; j++) {
1152c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
1153c8883902SJed Brown       if (a[i*nc+j]) {
1154c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1155d8588912SDave May       }
1156d8588912SDave May     }
1157d8588912SDave May   }
1158d8588912SDave May 
11598188e55aSJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1160d8588912SDave May 
1161854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->row_len);CHKERRQ(ierr);
1162854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&s->col_len);CHKERRQ(ierr);
1163c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
1164c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
1165d8588912SDave May 
1166*06a1af2fSStefano Zampini   ierr = PetscCalloc1(nr*nc,&s->nnzstate);CHKERRQ(ierr);
1167*06a1af2fSStefano Zampini   for (i=0; i<nr; i++) {
1168*06a1af2fSStefano Zampini     for (j=0; j<nc; j++) {
1169*06a1af2fSStefano Zampini       if (s->m[i][j]) {
1170*06a1af2fSStefano Zampini         ierr = MatGetNonzeroState(s->m[i][j],&s->nnzstate[i*nc+j]);CHKERRQ(ierr);
1171*06a1af2fSStefano Zampini       }
1172*06a1af2fSStefano Zampini     }
1173*06a1af2fSStefano Zampini   }
1174*06a1af2fSStefano Zampini 
11758188e55aSJed Brown   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1176d8588912SDave May 
1177c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1178c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1179c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1180c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1181c8883902SJed Brown 
1182c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1183c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1184c8883902SJed Brown 
1185*06a1af2fSStefano Zampini   /* disable operations that are not supported for non-square matrices,
1186*06a1af2fSStefano Zampini      or matrices for which is_row != is_col  */
1187*06a1af2fSStefano Zampini   ierr = MatHasCongruentLayouts(A,&cong);CHKERRQ(ierr);
1188*06a1af2fSStefano Zampini   if (cong && nr != nc) cong = PETSC_FALSE;
1189*06a1af2fSStefano Zampini   if (cong) {
1190*06a1af2fSStefano Zampini     for (i = 0; cong && i < nr; i++) {
1191*06a1af2fSStefano Zampini       ierr = ISEqual(s->isglobal.row[i],s->isglobal.col[i],&cong);CHKERRQ(ierr);
1192*06a1af2fSStefano Zampini     }
1193*06a1af2fSStefano Zampini   }
1194*06a1af2fSStefano Zampini   if (!cong) {
1195*06a1af2fSStefano Zampini     A->ops->getdiagonal = NULL;
1196*06a1af2fSStefano Zampini     A->ops->shift       = NULL;
1197*06a1af2fSStefano Zampini     A->ops->diagonalset = NULL;
1198*06a1af2fSStefano Zampini   }
1199*06a1af2fSStefano Zampini 
12001795a4d1SJed Brown   ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr);
1201*06a1af2fSStefano Zampini   ierr = PetscObjectStateIncrease((PetscObject)A);CHKERRQ(ierr);
1202*06a1af2fSStefano Zampini   A->nonzerostate++;
1203d8588912SDave May   PetscFunctionReturn(0);
1204d8588912SDave May }
1205d8588912SDave May 
1206c8883902SJed Brown /*@
1207c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
1208c8883902SJed Brown 
1209c8883902SJed Brown    Collective on Mat
1210c8883902SJed Brown 
1211c8883902SJed Brown    Input Parameter:
1212ffd6319bSRichard Tran Mills +  A - nested matrix
1213c8883902SJed Brown .  nr - number of nested row blocks
12140298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1215c8883902SJed Brown .  nc - number of nested column blocks
12160298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
12170298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1218c8883902SJed Brown 
1219*06a1af2fSStefano Zampini    Notes: this always resets any submatrix information previously set
1220*06a1af2fSStefano Zampini 
1221c8883902SJed Brown    Level: advanced
1222c8883902SJed Brown 
1223c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
1224c8883902SJed Brown @*/
1225c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1226c8883902SJed Brown {
1227c8883902SJed Brown   PetscErrorCode ierr;
1228*06a1af2fSStefano Zampini   PetscInt       i;
1229c8883902SJed Brown 
1230c8883902SJed Brown   PetscFunctionBegin;
1231c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1232ce94432eSBarry Smith   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1233c8883902SJed Brown   if (nr && is_row) {
1234c8883902SJed Brown     PetscValidPointer(is_row,3);
1235c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1236c8883902SJed Brown   }
1237ce94432eSBarry Smith   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
12381664e352SJed Brown   if (nc && is_col) {
1239c8883902SJed Brown     PetscValidPointer(is_col,5);
12409b30a8f6SBarry Smith     for (i=0; i<nc; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1241c8883902SJed Brown   }
1242*06a1af2fSStefano Zampini   if (nr*nc > 0) PetscValidPointer(a,6);
1243c8883902SJed 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);
1244c8883902SJed Brown   PetscFunctionReturn(0);
1245c8883902SJed Brown }
1246d8588912SDave May 
124745b6f7e9SBarry Smith static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog)
124877019fcaSJed Brown {
124977019fcaSJed Brown   PetscErrorCode ierr;
125077019fcaSJed Brown   PetscBool      flg;
125177019fcaSJed Brown   PetscInt       i,j,m,mi,*ix;
125277019fcaSJed Brown 
125377019fcaSJed Brown   PetscFunctionBegin;
125477019fcaSJed Brown   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
125577019fcaSJed Brown     if (islocal[i]) {
125677019fcaSJed Brown       ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
125777019fcaSJed Brown       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
125877019fcaSJed Brown     } else {
125977019fcaSJed Brown       ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
126077019fcaSJed Brown     }
126177019fcaSJed Brown     m += mi;
126277019fcaSJed Brown   }
126377019fcaSJed Brown   if (flg) {
1264785e854fSJed Brown     ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr);
1265165cd838SBarry Smith     for (i=0,m=0; i<n; i++) {
12660298fd71SBarry Smith       ISLocalToGlobalMapping smap = NULL;
1267e108cb99SStefano Zampini       Mat                    sub = NULL;
1268f6d38dbbSStefano Zampini       PetscSF                sf;
1269f6d38dbbSStefano Zampini       PetscLayout            map;
1270f6d38dbbSStefano Zampini       PetscInt               *ix2;
127177019fcaSJed Brown 
1272165cd838SBarry Smith       if (!colflg) {
127377019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
127477019fcaSJed Brown       } else {
127577019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
127677019fcaSJed Brown       }
1277191fd14bSBarry Smith       if (sub) {
1278191fd14bSBarry Smith         if (!colflg) {
1279191fd14bSBarry Smith           ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);
1280191fd14bSBarry Smith         } else {
1281191fd14bSBarry Smith           ierr = MatGetLocalToGlobalMapping(sub,NULL,&smap);CHKERRQ(ierr);
1282191fd14bSBarry Smith         }
1283191fd14bSBarry Smith       }
128477019fcaSJed Brown       if (islocal[i]) {
128577019fcaSJed Brown         ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
128677019fcaSJed Brown       } else {
128777019fcaSJed Brown         ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
128877019fcaSJed Brown       }
128977019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = j;
129077019fcaSJed Brown       if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);}
1291165cd838SBarry Smith 
129277019fcaSJed Brown       /*
129377019fcaSJed Brown         Now we need to extract the monolithic global indices that correspond to the given split global indices.
129477019fcaSJed 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.
129577019fcaSJed Brown        */
1296f6d38dbbSStefano Zampini       ierr = PetscMalloc1(mi,&ix2);CHKERRQ(ierr);
1297f6d38dbbSStefano Zampini       ierr = PetscSFCreate(((PetscObject)isglobal[i])->comm,&sf);CHKERRQ(ierr);
1298f6d38dbbSStefano Zampini       ierr = PetscLayoutCreate(((PetscObject)isglobal[i])->comm,&map);CHKERRQ(ierr);
1299f6d38dbbSStefano Zampini       ierr = PetscLayoutSetLocalSize(map,mi);CHKERRQ(ierr);
1300f6d38dbbSStefano Zampini       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
1301f6d38dbbSStefano Zampini       ierr = PetscSFSetGraphLayout(sf,map,mi,NULL,PETSC_USE_POINTER,ix+m);CHKERRQ(ierr);
1302f6d38dbbSStefano Zampini       ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
1303f6d38dbbSStefano Zampini       for (j=0; j<mi; j++) ix2[j] = ix[m+j];
1304f6d38dbbSStefano Zampini       ierr = PetscSFBcastBegin(sf,MPIU_INT,ix2,ix + m);CHKERRQ(ierr);
1305f6d38dbbSStefano Zampini       ierr = PetscSFBcastEnd(sf,MPIU_INT,ix2,ix + m);CHKERRQ(ierr);
1306f6d38dbbSStefano Zampini       ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1307f6d38dbbSStefano Zampini       ierr = PetscFree(ix2);CHKERRQ(ierr);
130877019fcaSJed Brown       m   += mi;
130977019fcaSJed Brown     }
1310f0413b6fSBarry Smith     ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
131177019fcaSJed Brown   } else {
13120298fd71SBarry Smith     *ltog = NULL;
131377019fcaSJed Brown   }
131477019fcaSJed Brown   PetscFunctionReturn(0);
131577019fcaSJed Brown }
131677019fcaSJed Brown 
131777019fcaSJed Brown 
1318d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1319d8588912SDave May /*
1320d8588912SDave May   nprocessors = NP
1321d8588912SDave May   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1322d8588912SDave May        proc 0: => (g_0,h_0,)
1323d8588912SDave May        proc 1: => (g_1,h_1,)
1324d8588912SDave May        ...
1325d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1326d8588912SDave May 
1327d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1328d8588912SDave May     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))
1329d8588912SDave May 
1330d8588912SDave May             proc 0:
1331d8588912SDave May     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1332d8588912SDave May             proc 1:
1333d8588912SDave May     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)
1334d8588912SDave May 
1335d8588912SDave May             proc NP-1:
1336d8588912SDave May     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1337d8588912SDave May */
1338841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1339d8588912SDave May {
1340e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
13418188e55aSJed Brown   PetscInt       i,j,offset,n,nsum,bs;
1342d8588912SDave May   PetscErrorCode ierr;
13430298fd71SBarry Smith   Mat            sub = NULL;
1344d8588912SDave May 
1345d8588912SDave May   PetscFunctionBegin;
1346854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&vs->isglobal.row);CHKERRQ(ierr);
1347854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&vs->isglobal.col);CHKERRQ(ierr);
1348d8588912SDave May   if (is_row) { /* valid IS is passed in */
1349d8588912SDave May     /* refs on is[] are incremeneted */
1350e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1351d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
135226fbe8dcSKarl Rupp 
1353e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1354d8588912SDave May     }
13552ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
13568188e55aSJed Brown     nsum = 0;
13578188e55aSJed Brown     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
13588188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1359ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
13600298fd71SBarry Smith       ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
1361ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
13628188e55aSJed Brown       nsum += n;
13638188e55aSJed Brown     }
1364ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
136530bc264bSJed Brown     offset -= nsum;
1366e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1367f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
13680298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
13692ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1370ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1371e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
13722ae74bdbSJed Brown       offset += n;
1373d8588912SDave May     }
1374d8588912SDave May   }
1375d8588912SDave May 
1376d8588912SDave May   if (is_col) { /* valid IS is passed in */
1377d8588912SDave May     /* refs on is[] are incremeneted */
1378e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1379d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
138026fbe8dcSKarl Rupp 
1381e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1382d8588912SDave May     }
13832ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
13842ae74bdbSJed Brown     offset = A->cmap->rstart;
13858188e55aSJed Brown     nsum   = 0;
13868188e55aSJed Brown     for (j=0; j<vs->nc; j++) {
13878188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1388ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
13890298fd71SBarry Smith       ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
1390ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
13918188e55aSJed Brown       nsum += n;
13928188e55aSJed Brown     }
1393ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
139430bc264bSJed Brown     offset -= nsum;
1395e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1396f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
13970298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
13982ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1399ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1400e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
14012ae74bdbSJed Brown       offset += n;
1402d8588912SDave May     }
1403d8588912SDave May   }
1404e2d7f03fSJed Brown 
1405e2d7f03fSJed Brown   /* Set up the local ISs */
1406785e854fSJed Brown   ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
1407785e854fSJed Brown   ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
1408e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1409e2d7f03fSJed Brown     IS                     isloc;
14100298fd71SBarry Smith     ISLocalToGlobalMapping rmap = NULL;
1411e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1412e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
14130298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);}
1414207556f9SJed Brown     if (rmap) {
1415e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1416e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1417e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1418e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1419207556f9SJed Brown     } else {
1420207556f9SJed Brown       nlocal = 0;
14210298fd71SBarry Smith       isloc  = NULL;
1422207556f9SJed Brown     }
1423e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1424e2d7f03fSJed Brown     offset            += nlocal;
1425e2d7f03fSJed Brown   }
14268188e55aSJed Brown   for (i=0,offset=0; i<vs->nc; i++) {
1427e2d7f03fSJed Brown     IS                     isloc;
14280298fd71SBarry Smith     ISLocalToGlobalMapping cmap = NULL;
1429e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1430e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
14310298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);}
1432207556f9SJed Brown     if (cmap) {
1433e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1434e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1435e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1436e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1437207556f9SJed Brown     } else {
1438207556f9SJed Brown       nlocal = 0;
14390298fd71SBarry Smith       isloc  = NULL;
1440207556f9SJed Brown     }
1441e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1442e2d7f03fSJed Brown     offset            += nlocal;
1443e2d7f03fSJed Brown   }
14440189643fSJed Brown 
144577019fcaSJed Brown   /* Set up the aggregate ISLocalToGlobalMapping */
144677019fcaSJed Brown   {
144745b6f7e9SBarry Smith     ISLocalToGlobalMapping rmap,cmap;
144845b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);CHKERRQ(ierr);
144945b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);CHKERRQ(ierr);
145077019fcaSJed Brown     if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);}
145177019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
145277019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
145377019fcaSJed Brown   }
145477019fcaSJed Brown 
14550189643fSJed Brown #if defined(PETSC_USE_DEBUG)
14560189643fSJed Brown   for (i=0; i<vs->nr; i++) {
14570189643fSJed Brown     for (j=0; j<vs->nc; j++) {
14580189643fSJed Brown       PetscInt m,n,M,N,mi,ni,Mi,Ni;
14590189643fSJed Brown       Mat      B = vs->m[i][j];
14600189643fSJed Brown       if (!B) continue;
14610189643fSJed Brown       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
14620189643fSJed Brown       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
14630189643fSJed Brown       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
14640189643fSJed Brown       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
14650189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
14660189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
1467ce94432eSBarry 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);
1468ce94432eSBarry 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);
14690189643fSJed Brown     }
14700189643fSJed Brown   }
14710189643fSJed Brown #endif
1472a061e289SJed Brown 
1473a061e289SJed Brown   /* Set A->assembled if all non-null blocks are currently assembled */
1474a061e289SJed Brown   for (i=0; i<vs->nr; i++) {
1475a061e289SJed Brown     for (j=0; j<vs->nc; j++) {
1476a061e289SJed Brown       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0);
1477a061e289SJed Brown     }
1478a061e289SJed Brown   }
1479a061e289SJed Brown   A->assembled = PETSC_TRUE;
1480d8588912SDave May   PetscFunctionReturn(0);
1481d8588912SDave May }
1482d8588912SDave May 
148345c38901SJed Brown /*@C
1484659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1485659c6bb0SJed Brown 
1486659c6bb0SJed Brown    Collective on Mat
1487659c6bb0SJed Brown 
1488659c6bb0SJed Brown    Input Parameter:
1489659c6bb0SJed Brown +  comm - Communicator for the new Mat
1490659c6bb0SJed Brown .  nr - number of nested row blocks
14910298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1492659c6bb0SJed Brown .  nc - number of nested column blocks
14930298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
14940298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1495659c6bb0SJed Brown 
1496659c6bb0SJed Brown    Output Parameter:
1497659c6bb0SJed Brown .  B - new matrix
1498659c6bb0SJed Brown 
1499659c6bb0SJed Brown    Level: advanced
1500659c6bb0SJed Brown 
1501950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1502659c6bb0SJed Brown @*/
15037087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1504d8588912SDave May {
1505d8588912SDave May   Mat            A;
1506d8588912SDave May   PetscErrorCode ierr;
1507d8588912SDave May 
1508d8588912SDave May   PetscFunctionBegin;
1509c8883902SJed Brown   *B   = 0;
1510d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1511c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
151291a28eb3SBarry Smith   A->preallocated = PETSC_TRUE;
1513c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1514d8588912SDave May   *B   = A;
1515d8588912SDave May   PetscFunctionReturn(0);
1516d8588912SDave May }
1517659c6bb0SJed Brown 
1518b68353e5Sstefano_zampini static PetscErrorCode MatConvert_Nest_SeqAIJ_fast(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1519b68353e5Sstefano_zampini {
1520b68353e5Sstefano_zampini   Mat_Nest       *nest = (Mat_Nest*)A->data;
152123875855Sstefano_zampini   Mat            *trans;
1522b68353e5Sstefano_zampini   PetscScalar    **avv;
1523b68353e5Sstefano_zampini   PetscScalar    *vv;
1524b68353e5Sstefano_zampini   PetscInt       **aii,**ajj;
1525b68353e5Sstefano_zampini   PetscInt       *ii,*jj,*ci;
1526b68353e5Sstefano_zampini   PetscInt       nr,nc,nnz,i,j;
1527b68353e5Sstefano_zampini   PetscBool      done;
1528b68353e5Sstefano_zampini   PetscErrorCode ierr;
1529b68353e5Sstefano_zampini 
1530b68353e5Sstefano_zampini   PetscFunctionBegin;
1531b68353e5Sstefano_zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
1532b68353e5Sstefano_zampini   if (reuse == MAT_REUSE_MATRIX) {
1533b68353e5Sstefano_zampini     PetscInt rnr;
1534b68353e5Sstefano_zampini 
1535b68353e5Sstefano_zampini     ierr = MatGetRowIJ(*newmat,0,PETSC_FALSE,PETSC_FALSE,&rnr,(const PetscInt**)&ii,(const PetscInt**)&jj,&done);CHKERRQ(ierr);
1536b68353e5Sstefano_zampini     if (!done) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatGetRowIJ");
1537b68353e5Sstefano_zampini     if (rnr != nr) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of rows");
1538b68353e5Sstefano_zampini     ierr = MatSeqAIJGetArray(*newmat,&vv);CHKERRQ(ierr);
1539b68353e5Sstefano_zampini   }
1540b68353e5Sstefano_zampini   /* extract CSR for nested SeqAIJ matrices */
1541b68353e5Sstefano_zampini   nnz  = 0;
154223875855Sstefano_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);
1543b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1544b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1545b68353e5Sstefano_zampini       Mat B = nest->m[i][j];
1546b68353e5Sstefano_zampini       if (B) {
1547b68353e5Sstefano_zampini         PetscScalar *naa;
1548b68353e5Sstefano_zampini         PetscInt    *nii,*njj,nnr;
154923875855Sstefano_zampini         PetscBool   istrans;
1550b68353e5Sstefano_zampini 
155123875855Sstefano_zampini         ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);CHKERRQ(ierr);
155223875855Sstefano_zampini         if (istrans) {
155323875855Sstefano_zampini           Mat Bt;
155423875855Sstefano_zampini 
155523875855Sstefano_zampini           ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr);
155623875855Sstefano_zampini           ierr = MatTranspose(Bt,MAT_INITIAL_MATRIX,&trans[i*nest->nc+j]);CHKERRQ(ierr);
155723875855Sstefano_zampini           B    = trans[i*nest->nc+j];
155823875855Sstefano_zampini         }
1559b68353e5Sstefano_zampini         ierr = MatGetRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&nii,(const PetscInt**)&njj,&done);CHKERRQ(ierr);
1560b68353e5Sstefano_zampini         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatGetRowIJ");
1561b68353e5Sstefano_zampini         ierr = MatSeqAIJGetArray(B,&naa);CHKERRQ(ierr);
1562b68353e5Sstefano_zampini         nnz += nii[nnr];
1563b68353e5Sstefano_zampini 
1564b68353e5Sstefano_zampini         aii[i*nest->nc+j] = nii;
1565b68353e5Sstefano_zampini         ajj[i*nest->nc+j] = njj;
1566b68353e5Sstefano_zampini         avv[i*nest->nc+j] = naa;
1567b68353e5Sstefano_zampini       }
1568b68353e5Sstefano_zampini     }
1569b68353e5Sstefano_zampini   }
1570b68353e5Sstefano_zampini   if (reuse != MAT_REUSE_MATRIX) {
1571b68353e5Sstefano_zampini     ierr = PetscMalloc1(nr+1,&ii);CHKERRQ(ierr);
1572b68353e5Sstefano_zampini     ierr = PetscMalloc1(nnz,&jj);CHKERRQ(ierr);
1573b68353e5Sstefano_zampini     ierr = PetscMalloc1(nnz,&vv);CHKERRQ(ierr);
1574b68353e5Sstefano_zampini   } else {
1575b68353e5Sstefano_zampini     if (nnz != ii[nr]) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of nonzeros");
1576b68353e5Sstefano_zampini   }
1577b68353e5Sstefano_zampini 
1578b68353e5Sstefano_zampini   /* new row pointer */
1579580bdb30SBarry Smith   ierr = PetscArrayzero(ii,nr+1);CHKERRQ(ierr);
1580b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1581b68353e5Sstefano_zampini     PetscInt       ncr,rst;
1582b68353e5Sstefano_zampini 
1583b68353e5Sstefano_zampini     ierr = ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);CHKERRQ(ierr);
1584b68353e5Sstefano_zampini     ierr = ISGetLocalSize(nest->isglobal.row[i],&ncr);CHKERRQ(ierr);
1585b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1586b68353e5Sstefano_zampini       if (aii[i*nest->nc+j]) {
1587b68353e5Sstefano_zampini         PetscInt    *nii = aii[i*nest->nc+j];
1588b68353e5Sstefano_zampini         PetscInt    ir;
1589b68353e5Sstefano_zampini 
1590b68353e5Sstefano_zampini         for (ir=rst; ir<ncr+rst; ++ir) {
1591b68353e5Sstefano_zampini           ii[ir+1] += nii[1]-nii[0];
1592b68353e5Sstefano_zampini           nii++;
1593b68353e5Sstefano_zampini         }
1594b68353e5Sstefano_zampini       }
1595b68353e5Sstefano_zampini     }
1596b68353e5Sstefano_zampini   }
1597b68353e5Sstefano_zampini   for (i=0; i<nr; i++) ii[i+1] += ii[i];
1598b68353e5Sstefano_zampini 
1599b68353e5Sstefano_zampini   /* construct CSR for the new matrix */
1600b68353e5Sstefano_zampini   ierr = PetscCalloc1(nr,&ci);CHKERRQ(ierr);
1601b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1602b68353e5Sstefano_zampini     PetscInt       ncr,rst;
1603b68353e5Sstefano_zampini 
1604b68353e5Sstefano_zampini     ierr = ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);CHKERRQ(ierr);
1605b68353e5Sstefano_zampini     ierr = ISGetLocalSize(nest->isglobal.row[i],&ncr);CHKERRQ(ierr);
1606b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1607b68353e5Sstefano_zampini       if (aii[i*nest->nc+j]) {
1608b68353e5Sstefano_zampini         PetscScalar *nvv = avv[i*nest->nc+j];
1609b68353e5Sstefano_zampini         PetscInt    *nii = aii[i*nest->nc+j];
1610b68353e5Sstefano_zampini         PetscInt    *njj = ajj[i*nest->nc+j];
1611b68353e5Sstefano_zampini         PetscInt    ir,cst;
1612b68353e5Sstefano_zampini 
1613b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.col[j],&cst,NULL);CHKERRQ(ierr);
1614b68353e5Sstefano_zampini         for (ir=rst; ir<ncr+rst; ++ir) {
1615b68353e5Sstefano_zampini           PetscInt ij,rsize = nii[1]-nii[0],ist = ii[ir]+ci[ir];
1616b68353e5Sstefano_zampini 
1617b68353e5Sstefano_zampini           for (ij=0;ij<rsize;ij++) {
1618b68353e5Sstefano_zampini             jj[ist+ij] = *njj+cst;
1619b68353e5Sstefano_zampini             vv[ist+ij] = *nvv;
1620b68353e5Sstefano_zampini             njj++;
1621b68353e5Sstefano_zampini             nvv++;
1622b68353e5Sstefano_zampini           }
1623b68353e5Sstefano_zampini           ci[ir] += rsize;
1624b68353e5Sstefano_zampini           nii++;
1625b68353e5Sstefano_zampini         }
1626b68353e5Sstefano_zampini       }
1627b68353e5Sstefano_zampini     }
1628b68353e5Sstefano_zampini   }
1629b68353e5Sstefano_zampini   ierr = PetscFree(ci);CHKERRQ(ierr);
1630b68353e5Sstefano_zampini 
1631b68353e5Sstefano_zampini   /* restore info */
1632b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1633b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1634b68353e5Sstefano_zampini       Mat B = nest->m[i][j];
1635b68353e5Sstefano_zampini       if (B) {
1636b68353e5Sstefano_zampini         PetscInt nnr = 0, k = i*nest->nc+j;
163723875855Sstefano_zampini 
163823875855Sstefano_zampini         B    = (trans[k] ? trans[k] : B);
1639b68353e5Sstefano_zampini         ierr = MatRestoreRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&aii[k],(const PetscInt**)&ajj[k],&done);CHKERRQ(ierr);
1640b68353e5Sstefano_zampini         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatRestoreRowIJ");
1641b68353e5Sstefano_zampini         ierr = MatSeqAIJRestoreArray(B,&avv[k]);CHKERRQ(ierr);
164223875855Sstefano_zampini         ierr = MatDestroy(&trans[k]);CHKERRQ(ierr);
1643b68353e5Sstefano_zampini       }
1644b68353e5Sstefano_zampini     }
1645b68353e5Sstefano_zampini   }
164623875855Sstefano_zampini   ierr = PetscFree4(aii,ajj,avv,trans);CHKERRQ(ierr);
1647b68353e5Sstefano_zampini 
1648b68353e5Sstefano_zampini   /* finalize newmat */
1649b68353e5Sstefano_zampini   if (reuse == MAT_INITIAL_MATRIX) {
1650b68353e5Sstefano_zampini     ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,newmat);CHKERRQ(ierr);
1651b68353e5Sstefano_zampini   } else if (reuse == MAT_INPLACE_MATRIX) {
1652b68353e5Sstefano_zampini     Mat B;
1653b68353e5Sstefano_zampini 
1654b68353e5Sstefano_zampini     ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,&B);CHKERRQ(ierr);
1655b68353e5Sstefano_zampini     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
1656b68353e5Sstefano_zampini   }
1657b68353e5Sstefano_zampini   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1658b68353e5Sstefano_zampini   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1659b68353e5Sstefano_zampini   {
1660b68353e5Sstefano_zampini     Mat_SeqAIJ *a = (Mat_SeqAIJ*)((*newmat)->data);
1661b68353e5Sstefano_zampini     a->free_a     = PETSC_TRUE;
1662b68353e5Sstefano_zampini     a->free_ij    = PETSC_TRUE;
1663b68353e5Sstefano_zampini   }
1664b68353e5Sstefano_zampini   PetscFunctionReturn(0);
1665b68353e5Sstefano_zampini }
1666b68353e5Sstefano_zampini 
1667cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1668629c3df2SDmitry Karpeev {
1669629c3df2SDmitry Karpeev   PetscErrorCode ierr;
1670629c3df2SDmitry Karpeev   Mat_Nest       *nest = (Mat_Nest*)A->data;
167183b1a929SMark Adams   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz,rstart;
1672649b366bSFande Kong   PetscInt       cstart,cend;
1673b68353e5Sstefano_zampini   PetscMPIInt    size;
1674629c3df2SDmitry Karpeev   Mat            C;
1675629c3df2SDmitry Karpeev 
1676629c3df2SDmitry Karpeev   PetscFunctionBegin;
1677b68353e5Sstefano_zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1678b68353e5Sstefano_zampini   if (size == 1) { /* look for a special case with SeqAIJ matrices and strided-1, contiguous, blocks */
1679b68353e5Sstefano_zampini     PetscInt  nf;
1680b68353e5Sstefano_zampini     PetscBool fast;
1681b68353e5Sstefano_zampini 
1682b68353e5Sstefano_zampini     ierr = PetscStrcmp(newtype,MATAIJ,&fast);CHKERRQ(ierr);
1683b68353e5Sstefano_zampini     if (!fast) {
1684b68353e5Sstefano_zampini       ierr = PetscStrcmp(newtype,MATSEQAIJ,&fast);CHKERRQ(ierr);
1685b68353e5Sstefano_zampini     }
1686b68353e5Sstefano_zampini     for (i=0; i<nest->nr && fast; ++i) {
1687b68353e5Sstefano_zampini       for (j=0; j<nest->nc && fast; ++j) {
1688b68353e5Sstefano_zampini         Mat B = nest->m[i][j];
1689b68353e5Sstefano_zampini         if (B) {
1690b68353e5Sstefano_zampini           ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQAIJ,&fast);CHKERRQ(ierr);
169123875855Sstefano_zampini           if (!fast) {
169223875855Sstefano_zampini             PetscBool istrans;
169323875855Sstefano_zampini 
169423875855Sstefano_zampini             ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);CHKERRQ(ierr);
169523875855Sstefano_zampini             if (istrans) {
169623875855Sstefano_zampini               Mat Bt;
169723875855Sstefano_zampini 
169823875855Sstefano_zampini               ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr);
169923875855Sstefano_zampini               ierr = PetscObjectTypeCompare((PetscObject)Bt,MATSEQAIJ,&fast);CHKERRQ(ierr);
170023875855Sstefano_zampini             }
1701b68353e5Sstefano_zampini           }
1702b68353e5Sstefano_zampini         }
1703b68353e5Sstefano_zampini       }
1704b68353e5Sstefano_zampini     }
1705b68353e5Sstefano_zampini     for (i=0, nf=0; i<nest->nr && fast; ++i) {
1706b68353e5Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)nest->isglobal.row[i],ISSTRIDE,&fast);CHKERRQ(ierr);
1707b68353e5Sstefano_zampini       if (fast) {
1708b68353e5Sstefano_zampini         PetscInt f,s;
1709b68353e5Sstefano_zampini 
1710b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.row[i],&f,&s);CHKERRQ(ierr);
1711b68353e5Sstefano_zampini         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1712b68353e5Sstefano_zampini         else {
1713b68353e5Sstefano_zampini           ierr = ISGetSize(nest->isglobal.row[i],&f);CHKERRQ(ierr);
1714b68353e5Sstefano_zampini           nf  += f;
1715b68353e5Sstefano_zampini         }
1716b68353e5Sstefano_zampini       }
1717b68353e5Sstefano_zampini     }
1718b68353e5Sstefano_zampini     for (i=0, nf=0; i<nest->nc && fast; ++i) {
1719b68353e5Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)nest->isglobal.col[i],ISSTRIDE,&fast);CHKERRQ(ierr);
1720b68353e5Sstefano_zampini       if (fast) {
1721b68353e5Sstefano_zampini         PetscInt f,s;
1722b68353e5Sstefano_zampini 
1723b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.col[i],&f,&s);CHKERRQ(ierr);
1724b68353e5Sstefano_zampini         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1725b68353e5Sstefano_zampini         else {
1726b68353e5Sstefano_zampini           ierr = ISGetSize(nest->isglobal.col[i],&f);CHKERRQ(ierr);
1727b68353e5Sstefano_zampini           nf  += f;
1728b68353e5Sstefano_zampini         }
1729b68353e5Sstefano_zampini       }
1730b68353e5Sstefano_zampini     }
1731b68353e5Sstefano_zampini     if (fast) {
1732b68353e5Sstefano_zampini       ierr = MatConvert_Nest_SeqAIJ_fast(A,newtype,reuse,newmat);CHKERRQ(ierr);
1733b68353e5Sstefano_zampini       PetscFunctionReturn(0);
1734b68353e5Sstefano_zampini     }
1735b68353e5Sstefano_zampini   }
1736629c3df2SDmitry Karpeev   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1737629c3df2SDmitry Karpeev   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
1738649b366bSFande Kong   ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr);
1739629c3df2SDmitry Karpeev   switch (reuse) {
1740629c3df2SDmitry Karpeev   case MAT_INITIAL_MATRIX:
1741ce94432eSBarry Smith     ierr    = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
1742629c3df2SDmitry Karpeev     ierr    = MatSetType(C,newtype);CHKERRQ(ierr);
1743629c3df2SDmitry Karpeev     ierr    = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
1744629c3df2SDmitry Karpeev     *newmat = C;
1745629c3df2SDmitry Karpeev     break;
1746629c3df2SDmitry Karpeev   case MAT_REUSE_MATRIX:
1747629c3df2SDmitry Karpeev     C = *newmat;
1748629c3df2SDmitry Karpeev     break;
1749ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1750629c3df2SDmitry Karpeev   }
1751785e854fSJed Brown   ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr);
1752629c3df2SDmitry Karpeev   onnz = dnnz + m;
1753629c3df2SDmitry Karpeev   for (k=0; k<m; k++) {
1754629c3df2SDmitry Karpeev     dnnz[k] = 0;
1755629c3df2SDmitry Karpeev     onnz[k] = 0;
1756629c3df2SDmitry Karpeev   }
1757629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1758629c3df2SDmitry Karpeev     IS             bNis;
1759629c3df2SDmitry Karpeev     PetscInt       bN;
1760629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1761629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1762629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1763629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr);
1764629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1765629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1766629c3df2SDmitry Karpeev       PetscSF        bmsf;
1767649b366bSFande Kong       PetscSFNode    *iremote;
1768629c3df2SDmitry Karpeev       Mat            B;
1769649b366bSFande Kong       PetscInt       bm, *sub_dnnz,*sub_onnz, br;
1770629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1771629c3df2SDmitry Karpeev       B = nest->m[i][j];
1772629c3df2SDmitry Karpeev       if (!B) continue;
1773629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1774629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1775ce94432eSBarry Smith       ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr);
1776649b366bSFande Kong       ierr = PetscMalloc1(bm,&iremote);CHKERRQ(ierr);
1777649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_dnnz);CHKERRQ(ierr);
1778649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_onnz);CHKERRQ(ierr);
1779649b366bSFande Kong       for (k = 0; k < bm; ++k){
1780649b366bSFande Kong     	sub_dnnz[k] = 0;
1781649b366bSFande Kong     	sub_onnz[k] = 0;
1782649b366bSFande Kong       }
1783629c3df2SDmitry Karpeev       /*
1784629c3df2SDmitry Karpeev        Locate the owners for all of the locally-owned global row indices for this row block.
1785629c3df2SDmitry Karpeev        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1786629c3df2SDmitry Karpeev        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1787629c3df2SDmitry Karpeev        */
178883b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1789629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1790649b366bSFande Kong         PetscInt       row = bmindices[br], rowowner = 0, brncols, col;
1791629c3df2SDmitry Karpeev         const PetscInt *brcols;
1792a4b3d3acSMatthew G Knepley         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1793629c3df2SDmitry Karpeev         ierr      = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr);
1794649b366bSFande Kong         /* how many roots  */
1795649b366bSFande Kong         iremote[br].rank = rowowner; iremote[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1796649b366bSFande Kong         /* get nonzero pattern */
179783b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1798629c3df2SDmitry Karpeev         for (k=0; k<brncols; k++) {
1799629c3df2SDmitry Karpeev           col  = bNindices[brcols[k]];
1800649b366bSFande Kong           if (col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]) {
1801649b366bSFande Kong             sub_dnnz[br]++;
1802649b366bSFande Kong           } else {
1803649b366bSFande Kong             sub_onnz[br]++;
1804649b366bSFande Kong           }
1805629c3df2SDmitry Karpeev         }
180683b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1807629c3df2SDmitry Karpeev       }
1808629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1809629c3df2SDmitry Karpeev       /* bsf will have to take care of disposing of bedges. */
1810649b366bSFande Kong       ierr = PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1811649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1812649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1813649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1814649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1815649b366bSFande Kong       ierr = PetscFree(sub_dnnz);CHKERRQ(ierr);
1816649b366bSFande Kong       ierr = PetscFree(sub_onnz);CHKERRQ(ierr);
1817629c3df2SDmitry Karpeev       ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr);
1818629c3df2SDmitry Karpeev     }
181922d28d08SBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1820629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
182165a4a0a3Sstefano_zampini   }
182265a4a0a3Sstefano_zampini   /* Resize preallocation if overestimated */
182365a4a0a3Sstefano_zampini   for (i=0;i<m;i++) {
182465a4a0a3Sstefano_zampini     dnnz[i] = PetscMin(dnnz[i],A->cmap->n);
182565a4a0a3Sstefano_zampini     onnz[i] = PetscMin(onnz[i],A->cmap->N - A->cmap->n);
1826629c3df2SDmitry Karpeev   }
1827629c3df2SDmitry Karpeev   ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr);
1828629c3df2SDmitry Karpeev   ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr);
1829629c3df2SDmitry Karpeev   ierr = PetscFree(dnnz);CHKERRQ(ierr);
1830629c3df2SDmitry Karpeev 
1831629c3df2SDmitry Karpeev   /* Fill by row */
1832629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1833629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1834629c3df2SDmitry Karpeev     IS             bNis;
1835629c3df2SDmitry Karpeev     PetscInt       bN;
1836629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1837629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1838629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr);
1839629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1840629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1841629c3df2SDmitry Karpeev       Mat            B;
1842629c3df2SDmitry Karpeev       PetscInt       bm, br;
1843629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1844629c3df2SDmitry Karpeev       B = nest->m[i][j];
1845629c3df2SDmitry Karpeev       if (!B) continue;
1846629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1847629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
184883b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1849629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1850629c3df2SDmitry Karpeev         PetscInt          row = bmindices[br], brncols,  *cols;
1851629c3df2SDmitry Karpeev         const PetscInt    *brcols;
1852629c3df2SDmitry Karpeev         const PetscScalar *brcoldata;
185383b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1854785e854fSJed Brown         ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr);
185526fbe8dcSKarl Rupp         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
1856629c3df2SDmitry Karpeev         /*
1857629c3df2SDmitry Karpeev           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1858629c3df2SDmitry Karpeev           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1859629c3df2SDmitry Karpeev          */
1860a2ea699eSBarry Smith         ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr);
186183b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1862629c3df2SDmitry Karpeev         ierr = PetscFree(cols);CHKERRQ(ierr);
1863629c3df2SDmitry Karpeev       }
1864629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1865629c3df2SDmitry Karpeev     }
1866a2ea699eSBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1867629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1868629c3df2SDmitry Karpeev   }
1869629c3df2SDmitry Karpeev   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1870629c3df2SDmitry Karpeev   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1871629c3df2SDmitry Karpeev   PetscFunctionReturn(0);
1872629c3df2SDmitry Karpeev }
1873629c3df2SDmitry Karpeev 
18748b7d3b4bSBarry Smith PetscErrorCode  MatHasOperation_Nest(Mat mat,MatOperation op,PetscBool  *has)
18758b7d3b4bSBarry Smith {
18768b7d3b4bSBarry Smith   PetscFunctionBegin;
18778b7d3b4bSBarry Smith   *has = PETSC_FALSE;
18788b7d3b4bSBarry Smith   if (op == MATOP_MULT_TRANSPOSE) {
18798b7d3b4bSBarry Smith     Mat_Nest       *bA = (Mat_Nest*)mat->data;
18808b7d3b4bSBarry Smith     PetscInt       i,j,nr = bA->nr,nc = bA->nc;
18818b7d3b4bSBarry Smith     PetscErrorCode ierr;
18828b7d3b4bSBarry Smith     PetscBool      flg;
18838b7d3b4bSBarry Smith 
18848b7d3b4bSBarry Smith     for (j=0; j<nc; j++) {
18858b7d3b4bSBarry Smith       for (i=0; i<nr; i++) {
18868b7d3b4bSBarry Smith         if (!bA->m[i][j]) continue;
18878b7d3b4bSBarry Smith         ierr = MatHasOperation(bA->m[i][j],MATOP_MULT_TRANSPOSE_ADD,&flg);CHKERRQ(ierr);
18888b7d3b4bSBarry Smith         if (!flg) PetscFunctionReturn(0);
18898b7d3b4bSBarry Smith       }
18908b7d3b4bSBarry Smith     }
18918b7d3b4bSBarry Smith   }
18928b7d3b4bSBarry Smith   if (((void**)mat->ops)[op]) *has =  PETSC_TRUE;
18938b7d3b4bSBarry Smith   PetscFunctionReturn(0);
18948b7d3b4bSBarry Smith }
18958b7d3b4bSBarry Smith 
1896659c6bb0SJed Brown /*MC
1897659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1898659c6bb0SJed Brown 
1899659c6bb0SJed Brown   Level: intermediate
1900659c6bb0SJed Brown 
1901659c6bb0SJed Brown   Notes:
1902659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1903659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1904950540a4SJed Brown   It is usually used with DMComposite and DMCreateMatrix()
1905659c6bb0SJed Brown 
19068b7d3b4bSBarry Smith   Each of the submatrices lives on the same MPI communicator as the original nest matrix (though they can have zero
19078b7d3b4bSBarry Smith   rows/columns on some processes.) Thus this is not meant for cases where the submatrices live on far fewer processes
19088b7d3b4bSBarry Smith   than the nest matrix.
19098b7d3b4bSBarry Smith 
1910659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1911659c6bb0SJed Brown M*/
19128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1913c8883902SJed Brown {
1914c8883902SJed Brown   Mat_Nest       *s;
1915c8883902SJed Brown   PetscErrorCode ierr;
1916c8883902SJed Brown 
1917c8883902SJed Brown   PetscFunctionBegin;
1918b00a9115SJed Brown   ierr    = PetscNewLog(A,&s);CHKERRQ(ierr);
1919c8883902SJed Brown   A->data = (void*)s;
1920e7c19651SJed Brown 
1921e7c19651SJed Brown   s->nr            = -1;
1922e7c19651SJed Brown   s->nc            = -1;
19230298fd71SBarry Smith   s->m             = NULL;
1924e7c19651SJed Brown   s->splitassembly = PETSC_FALSE;
1925c8883902SJed Brown 
1926c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
192726fbe8dcSKarl Rupp 
1928c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
19299194d70fSJed Brown   A->ops->multadd               = MatMultAdd_Nest;
1930c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
19319194d70fSJed Brown   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1932f8170845SAlex Fikl   A->ops->transpose             = MatTranspose_Nest;
1933c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1934c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1935c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1936c222c20dSDavid Ham   A->ops->copy                  = MatCopy_Nest;
19376e76ffeaSPierre Jolivet   A->ops->axpy                  = MatAXPY_Nest;
1938c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
19397dae84e0SHong Zhang   A->ops->createsubmatrix       = MatCreateSubMatrix_Nest;
1940c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1941c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1942c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1943c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1944c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1945429bac76SJed Brown   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1946429bac76SJed Brown   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1947a061e289SJed Brown   A->ops->scale                 = MatScale_Nest;
1948a061e289SJed Brown   A->ops->shift                 = MatShift_Nest;
194913135bc6SAlex Fikl   A->ops->diagonalset           = MatDiagonalSet_Nest;
1950f8170845SAlex Fikl   A->ops->setrandom             = MatSetRandom_Nest;
19518b7d3b4bSBarry Smith   A->ops->hasoperation          = MatHasOperation_Nest;
1952c8883902SJed Brown 
1953c8883902SJed Brown   A->spptr        = 0;
1954c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1955c8883902SJed Brown 
1956c8883902SJed Brown   /* expose Nest api's */
1957bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",      MatNestGetSubMat_Nest);CHKERRQ(ierr);
1958bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",      MatNestSetSubMat_Nest);CHKERRQ(ierr);
1959bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",     MatNestGetSubMats_Nest);CHKERRQ(ierr);
1960bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",        MatNestGetSize_Nest);CHKERRQ(ierr);
1961bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",         MatNestGetISs_Nest);CHKERRQ(ierr);
1962bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",    MatNestGetLocalISs_Nest);CHKERRQ(ierr);
1963bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",     MatNestSetVecType_Nest);CHKERRQ(ierr);
1964bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",     MatNestSetSubMats_Nest);CHKERRQ(ierr);
19650899c546SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_mpiaij_C",MatConvert_Nest_AIJ);CHKERRQ(ierr);
19660899c546SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_seqaij_C",MatConvert_Nest_AIJ);CHKERRQ(ierr);
196783b1a929SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",   MatConvert_Nest_AIJ);CHKERRQ(ierr);
19685e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",    MatConvert_Nest_IS);CHKERRQ(ierr);
1969c8883902SJed Brown 
1970c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1971c8883902SJed Brown   PetscFunctionReturn(0);
1972c8883902SJed Brown }
1973