xref: /petsc/src/mat/impls/is/matis.c (revision 45471136b65808b32afdd814127af7e99c594e2c)
1b4319ba4SBarry Smith /*
2b4319ba4SBarry Smith     Creates a matrix class for using the Neumann-Neumann type preconditioners.
3b4319ba4SBarry Smith     This stores the matrices in globally unassembled form. Each processor
4b4319ba4SBarry Smith     assembles only its local Neumann problem and the parallel matrix vector
5b4319ba4SBarry Smith     product is handled "implicitly".
6b4319ba4SBarry Smith 
7b4319ba4SBarry Smith     Currently this allows for only one subdomain per processor.
8b4319ba4SBarry Smith */
9b4319ba4SBarry Smith 
10c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h>      /*I "petscmat.h" I*/
114f2d7cafSStefano Zampini #include <petsc/private/sfimpl.h>
1228f4e0baSStefano Zampini 
13f26d0771SStefano Zampini #define MATIS_MAX_ENTRIES_INSERTION 2048
14f26d0771SStefano Zampini 
15cf0a3239SStefano Zampini #undef __FUNCT__
16cf0a3239SStefano Zampini #define __FUNCT__ "MatISSetUpSF"
17cf0a3239SStefano Zampini /*@
183d996552SStefano Zampini    MatISSetUpSF - Setup star forest objects used by MatIS.
19cf0a3239SStefano Zampini 
20cf0a3239SStefano Zampini    Collective on MPI_Comm
21cf0a3239SStefano Zampini 
22cf0a3239SStefano Zampini    Input Parameters:
23cf0a3239SStefano Zampini +  A - the matrix
24cf0a3239SStefano Zampini 
25cf0a3239SStefano Zampini    Level: advanced
26cf0a3239SStefano Zampini 
273d996552SStefano Zampini    Notes: This function does not need to be called by the user.
28cf0a3239SStefano Zampini 
29cf0a3239SStefano Zampini .keywords: matrix
30cf0a3239SStefano Zampini 
31cf0a3239SStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatISSetPreallocation(), MatISGetLocalMat()
32cf0a3239SStefano Zampini @*/
33cf0a3239SStefano Zampini PetscErrorCode  MatISSetUpSF(Mat A)
34cf0a3239SStefano Zampini {
35cf0a3239SStefano Zampini   PetscErrorCode ierr;
36cf0a3239SStefano Zampini 
37cf0a3239SStefano Zampini   PetscFunctionBegin;
38cf0a3239SStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
39cf0a3239SStefano Zampini   PetscValidType(A,1);
40cf0a3239SStefano Zampini   ierr = PetscTryMethod(A,"MatISSetUpSF_C",(Mat),(A));CHKERRQ(ierr);
41cf0a3239SStefano Zampini   PetscFunctionReturn(0);
42cf0a3239SStefano Zampini }
43a72627d2SStefano Zampini 
4428f4e0baSStefano Zampini #undef __FUNCT__
453fd1c9e7SStefano Zampini #define __FUNCT__ "MatDiagonalSet_IS"
463fd1c9e7SStefano Zampini PetscErrorCode  MatDiagonalSet_IS(Mat A,Vec D,InsertMode insmode)
473fd1c9e7SStefano Zampini {
483fd1c9e7SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
493fd1c9e7SStefano Zampini   PetscErrorCode ierr;
503fd1c9e7SStefano Zampini 
513fd1c9e7SStefano Zampini   PetscFunctionBegin;
523fd1c9e7SStefano Zampini   if (!D) { /* this code branch is used by MatShift_IS */
533fd1c9e7SStefano Zampini     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
543fd1c9e7SStefano Zampini   } else {
553fd1c9e7SStefano Zampini     ierr = VecScatterBegin(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
563fd1c9e7SStefano Zampini     ierr = VecScatterEnd(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
573fd1c9e7SStefano Zampini   }
583fd1c9e7SStefano Zampini   ierr = VecPointwiseDivide(is->y,is->y,is->counter);CHKERRQ(ierr);
593fd1c9e7SStefano Zampini   ierr = MatDiagonalSet(is->A,is->y,insmode);CHKERRQ(ierr);
603fd1c9e7SStefano Zampini   PetscFunctionReturn(0);
613fd1c9e7SStefano Zampini }
623fd1c9e7SStefano Zampini 
633fd1c9e7SStefano Zampini #undef __FUNCT__
643fd1c9e7SStefano Zampini #define __FUNCT__ "MatShift_IS"
653fd1c9e7SStefano Zampini PetscErrorCode  MatShift_IS(Mat A,PetscScalar a)
663fd1c9e7SStefano Zampini {
673fd1c9e7SStefano Zampini   PetscErrorCode ierr;
683fd1c9e7SStefano Zampini 
693fd1c9e7SStefano Zampini   PetscFunctionBegin;
703fd1c9e7SStefano Zampini   ierr = MatDiagonalSet_IS(A,NULL,ADD_VALUES);CHKERRQ(ierr);
713fd1c9e7SStefano Zampini   PetscFunctionReturn(0);
723fd1c9e7SStefano Zampini }
733fd1c9e7SStefano Zampini 
743fd1c9e7SStefano Zampini #undef __FUNCT__
75f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesLocal_SubMat_IS"
76f26d0771SStefano Zampini static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
77f26d0771SStefano Zampini {
78f26d0771SStefano Zampini   PetscErrorCode ierr;
79f26d0771SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
80f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
81f26d0771SStefano Zampini 
82f26d0771SStefano Zampini   PetscFunctionBegin;
83f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
84f26d0771SStefano Zampini   if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n);
85f26d0771SStefano Zampini #endif
86f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr);
87f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr);
88f26d0771SStefano Zampini   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
89f26d0771SStefano Zampini   PetscFunctionReturn(0);
90f26d0771SStefano Zampini }
91f26d0771SStefano Zampini 
92f26d0771SStefano Zampini #undef __FUNCT__
93f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesBlockedLocal_SubMat_IS"
94f26d0771SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
95f26d0771SStefano Zampini {
96f26d0771SStefano Zampini   PetscErrorCode ierr;
97f26d0771SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
98f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
99f26d0771SStefano Zampini 
100f26d0771SStefano Zampini   PetscFunctionBegin;
101f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
102f26d0771SStefano Zampini   if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column block indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n);
103f26d0771SStefano Zampini #endif
104f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr);
105f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr);
106f26d0771SStefano Zampini   ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
107f26d0771SStefano Zampini   PetscFunctionReturn(0);
108f26d0771SStefano Zampini }
109f26d0771SStefano Zampini 
110f26d0771SStefano Zampini #undef __FUNCT__
111a8116848SStefano Zampini #define __FUNCT__ "PetscLayoutMapLocal_Private"
112f0ae7da4SStefano Zampini static PetscErrorCode PetscLayoutMapLocal_Private(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs)
113a8116848SStefano Zampini {
114a8116848SStefano Zampini   PetscInt      *owners = map->range;
115a8116848SStefano Zampini   PetscInt       n      = map->n;
116a8116848SStefano Zampini   PetscSF        sf;
117a8116848SStefano Zampini   PetscInt      *lidxs,*work = NULL;
118a8116848SStefano Zampini   PetscSFNode   *ridxs;
119a8116848SStefano Zampini   PetscMPIInt    rank;
120a8116848SStefano Zampini   PetscInt       r, p = 0, len = 0;
121a8116848SStefano Zampini   PetscErrorCode ierr;
122a8116848SStefano Zampini 
123a8116848SStefano Zampini   PetscFunctionBegin;
124a8116848SStefano Zampini   /* Create SF where leaves are input idxs and roots are owned idxs (code adapted from MatZeroRowsMapLocal_Private) */
125a8116848SStefano Zampini   ierr = MPI_Comm_rank(map->comm,&rank);CHKERRQ(ierr);
126a8116848SStefano Zampini   ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr);
127a8116848SStefano Zampini   for (r = 0; r < n; ++r) lidxs[r] = -1;
128a8116848SStefano Zampini   ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr);
129a8116848SStefano Zampini   for (r = 0; r < N; ++r) {
130a8116848SStefano Zampini     const PetscInt idx = idxs[r];
131a8116848SStefano Zampini     if (idx < 0 || map->N <= idx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D out of range [0,%D)",idx,map->N);
132a8116848SStefano Zampini     if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */
133a8116848SStefano Zampini       ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr);
134a8116848SStefano Zampini     }
135a8116848SStefano Zampini     ridxs[r].rank = p;
136a8116848SStefano Zampini     ridxs[r].index = idxs[r] - owners[p];
137a8116848SStefano Zampini   }
138a8116848SStefano Zampini   ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr);
139a8116848SStefano Zampini   ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr);
140a8116848SStefano Zampini   ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr);
141a8116848SStefano Zampini   ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr);
142f0ae7da4SStefano Zampini   if (ogidxs) { /* communicate global idxs */
143a8116848SStefano Zampini     PetscInt cum = 0,start,*work2;
144f0ae7da4SStefano Zampini 
145f0ae7da4SStefano Zampini     ierr = PetscMalloc1(n,&work);CHKERRQ(ierr);
146a8116848SStefano Zampini     ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr);
147a8116848SStefano Zampini     for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++;
148a8116848SStefano Zampini     ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRQ(ierr);
149a8116848SStefano Zampini     start -= cum;
150a8116848SStefano Zampini     cum = 0;
151a8116848SStefano Zampini     for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++;
152a8116848SStefano Zampini     ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr);
153a8116848SStefano Zampini     ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr);
154a8116848SStefano Zampini     ierr = PetscFree(work2);CHKERRQ(ierr);
155a8116848SStefano Zampini   }
156a8116848SStefano Zampini   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
157a8116848SStefano Zampini   /* Compress and put in indices */
158a8116848SStefano Zampini   for (r = 0; r < n; ++r)
159a8116848SStefano Zampini     if (lidxs[r] >= 0) {
160a8116848SStefano Zampini       if (work) work[len] = work[r];
161a8116848SStefano Zampini       lidxs[len++] = r;
162a8116848SStefano Zampini     }
163a8116848SStefano Zampini   if (on) *on = len;
164a8116848SStefano Zampini   if (oidxs) *oidxs = lidxs;
165a8116848SStefano Zampini   if (ogidxs) *ogidxs = work;
166a8116848SStefano Zampini   PetscFunctionReturn(0);
167a8116848SStefano Zampini }
168a8116848SStefano Zampini 
169a8116848SStefano Zampini #undef __FUNCT__
170a8116848SStefano Zampini #define __FUNCT__ "MatGetSubMatrix_IS"
171a8116848SStefano Zampini static PetscErrorCode MatGetSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat)
172a8116848SStefano Zampini {
173a8116848SStefano Zampini   Mat               locmat,newlocmat;
174a8116848SStefano Zampini   Mat_IS            *newmatis;
175a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG)
176a8116848SStefano Zampini   Vec               rtest,ltest;
177a8116848SStefano Zampini   const PetscScalar *array;
178a8116848SStefano Zampini #endif
179a8116848SStefano Zampini   const PetscInt    *idxs;
180a8116848SStefano Zampini   PetscInt          i,m,n;
181a8116848SStefano Zampini   PetscErrorCode    ierr;
182a8116848SStefano Zampini 
183a8116848SStefano Zampini   PetscFunctionBegin;
184a8116848SStefano Zampini   if (scall == MAT_REUSE_MATRIX) {
185a8116848SStefano Zampini     PetscBool ismatis;
186a8116848SStefano Zampini 
187a8116848SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr);
188a8116848SStefano Zampini     if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type");
189a8116848SStefano Zampini     newmatis = (Mat_IS*)(*newmat)->data;
190a8116848SStefano Zampini     if (!newmatis->getsub_ris) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS");
191a8116848SStefano Zampini     if (!newmatis->getsub_cis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS");
192a8116848SStefano Zampini   }
193a8116848SStefano Zampini   /* irow and icol may not have duplicate entries */
194a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG)
195a8116848SStefano Zampini   ierr = MatCreateVecs(mat,&ltest,&rtest);CHKERRQ(ierr);
196a8116848SStefano Zampini   ierr = ISGetLocalSize(irow,&n);CHKERRQ(ierr);
197a8116848SStefano Zampini   ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr);
198a8116848SStefano Zampini   for (i=0;i<n;i++) {
199a8116848SStefano Zampini     ierr = VecSetValue(rtest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr);
200a8116848SStefano Zampini   }
201a8116848SStefano Zampini   ierr = VecAssemblyBegin(rtest);CHKERRQ(ierr);
202a8116848SStefano Zampini   ierr = VecAssemblyEnd(rtest);CHKERRQ(ierr);
203a8116848SStefano Zampini   ierr = VecGetLocalSize(rtest,&n);CHKERRQ(ierr);
204a8116848SStefano Zampini   ierr = VecGetOwnershipRange(rtest,&m,NULL);CHKERRQ(ierr);
205a8116848SStefano Zampini   ierr = VecGetArrayRead(rtest,&array);CHKERRQ(ierr);
206fd479f66SMatthew G. Knepley   for (i=0;i<n;i++) if (array[i] != 0. && array[i] != 1.) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Index %D counted %D times! Irow may not have duplicate entries",i+m,(PetscInt)PetscRealPart(array[i]));
207a8116848SStefano Zampini   ierr = VecRestoreArrayRead(rtest,&array);CHKERRQ(ierr);
208a8116848SStefano Zampini   ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr);
209a8116848SStefano Zampini   ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr);
210a8116848SStefano Zampini   ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr);
211a8116848SStefano Zampini   for (i=0;i<n;i++) {
212a8116848SStefano Zampini     ierr = VecSetValue(ltest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr);
213a8116848SStefano Zampini   }
214a8116848SStefano Zampini   ierr = VecAssemblyBegin(ltest);CHKERRQ(ierr);
215a8116848SStefano Zampini   ierr = VecAssemblyEnd(ltest);CHKERRQ(ierr);
216a8116848SStefano Zampini   ierr = VecGetLocalSize(ltest,&n);CHKERRQ(ierr);
217a8116848SStefano Zampini   ierr = VecGetOwnershipRange(ltest,&m,NULL);CHKERRQ(ierr);
218a8116848SStefano Zampini   ierr = VecGetArrayRead(ltest,&array);CHKERRQ(ierr);
219fd479f66SMatthew G. Knepley   for (i=0;i<n;i++) if (array[i] != 0. && array[i] != 1.) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Index %D counted %D times! Icol may not have duplicate entries",i+m,(PetscInt)PetscRealPart(array[i]));
220a8116848SStefano Zampini   ierr = VecRestoreArrayRead(ltest,&array);CHKERRQ(ierr);
221a8116848SStefano Zampini   ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr);
222a8116848SStefano Zampini   ierr = VecDestroy(&rtest);CHKERRQ(ierr);
223a8116848SStefano Zampini   ierr = VecDestroy(&ltest);CHKERRQ(ierr);
224a8116848SStefano Zampini #endif
225a8116848SStefano Zampini   if (scall == MAT_INITIAL_MATRIX) {
226a8116848SStefano Zampini     Mat_IS                 *matis = (Mat_IS*)mat->data;
227a8116848SStefano Zampini     ISLocalToGlobalMapping rl2g;
228a8116848SStefano Zampini     IS                     is;
229a8116848SStefano Zampini     PetscInt               *lidxs,*lgidxs,*newgidxs;
2306dd40735SStefano Zampini     PetscInt               ll,newloc;
231a8116848SStefano Zampini     MPI_Comm               comm;
232a8116848SStefano Zampini 
233a8116848SStefano Zampini     ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
234a8116848SStefano Zampini     ierr = ISGetLocalSize(irow,&m);CHKERRQ(ierr);
235a8116848SStefano Zampini     ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr);
236a8116848SStefano Zampini     ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
237a8116848SStefano Zampini     ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr);
238a8116848SStefano Zampini     ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
239a8116848SStefano Zampini     /* communicate irow to their owners in the layout */
240a8116848SStefano Zampini     ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr);
241f0ae7da4SStefano Zampini     ierr = PetscLayoutMapLocal_Private(mat->rmap,m,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr);
242a8116848SStefano Zampini     ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr);
2433d996552SStefano Zampini     ierr = MatISSetUpSF(mat);CHKERRQ(ierr);
2443d996552SStefano Zampini     ierr = PetscMemzero(matis->sf_rootdata,matis->sf->nroots*sizeof(PetscInt));CHKERRQ(ierr);
245a8116848SStefano Zampini     for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1;
246a8116848SStefano Zampini     ierr = PetscFree(lidxs);CHKERRQ(ierr);
247a8116848SStefano Zampini     ierr = PetscFree(lgidxs);CHKERRQ(ierr);
248a8116848SStefano Zampini     ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
249a8116848SStefano Zampini     ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
2503d996552SStefano Zampini     for (i=0,newloc=0;i<matis->sf->nleaves;i++) if (matis->sf_leafdata[i]) newloc++;
251a8116848SStefano Zampini     ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr);
252a8116848SStefano Zampini     ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr);
2533d996552SStefano Zampini     for (i=0,newloc=0;i<matis->sf->nleaves;i++)
254a8116848SStefano Zampini       if (matis->sf_leafdata[i]) {
255a8116848SStefano Zampini         lidxs[newloc] = i;
256a8116848SStefano Zampini         newgidxs[newloc++] = matis->sf_leafdata[i]-1;
257a8116848SStefano Zampini       }
258a8116848SStefano Zampini     ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
259a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
260a8116848SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
261a8116848SStefano Zampini     /* local is to extract local submatrix */
262a8116848SStefano Zampini     newmatis = (Mat_IS*)(*newmat)->data;
263a8116848SStefano Zampini     ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris);CHKERRQ(ierr);
264a8116848SStefano Zampini     if (mat->congruentlayouts == -1) { /* first time we compare rows and cols layouts */
265a8116848SStefano Zampini       PetscBool cong;
266a8116848SStefano Zampini       ierr = PetscLayoutCompare(mat->rmap,mat->cmap,&cong);CHKERRQ(ierr);
267a8116848SStefano Zampini       if (cong) mat->congruentlayouts = 1;
268a8116848SStefano Zampini       else      mat->congruentlayouts = 0;
269a8116848SStefano Zampini     }
270a8116848SStefano Zampini     if (mat->congruentlayouts && irow == icol) {
271a8116848SStefano Zampini       ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g);CHKERRQ(ierr);
272a8116848SStefano Zampini       ierr = PetscObjectReference((PetscObject)newmatis->getsub_ris);CHKERRQ(ierr);
273a8116848SStefano Zampini       newmatis->getsub_cis = newmatis->getsub_ris;
274a8116848SStefano Zampini     } else {
275a8116848SStefano Zampini       ISLocalToGlobalMapping cl2g;
276a8116848SStefano Zampini 
277a8116848SStefano Zampini       /* communicate icol to their owners in the layout */
278a8116848SStefano Zampini       ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr);
279f0ae7da4SStefano Zampini       ierr = PetscLayoutMapLocal_Private(mat->cmap,n,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr);
280a8116848SStefano Zampini       ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr);
2813d996552SStefano Zampini       ierr = PetscMemzero(matis->csf_rootdata,matis->csf->nroots*sizeof(PetscInt));CHKERRQ(ierr);
282a8116848SStefano Zampini       for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1;
283a8116848SStefano Zampini       ierr = PetscFree(lidxs);CHKERRQ(ierr);
284a8116848SStefano Zampini       ierr = PetscFree(lgidxs);CHKERRQ(ierr);
285a8116848SStefano Zampini       ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr);
286a8116848SStefano Zampini       ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr);
2873d996552SStefano Zampini       for (i=0,newloc=0;i<matis->csf->nleaves;i++) if (matis->csf_leafdata[i]) newloc++;
288a8116848SStefano Zampini       ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr);
289a8116848SStefano Zampini       ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr);
2903d996552SStefano Zampini       for (i=0,newloc=0;i<matis->csf->nleaves;i++)
291a8116848SStefano Zampini         if (matis->csf_leafdata[i]) {
292a8116848SStefano Zampini           lidxs[newloc] = i;
293a8116848SStefano Zampini           newgidxs[newloc++] = matis->csf_leafdata[i]-1;
294a8116848SStefano Zampini         }
295a8116848SStefano Zampini       ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
296a8116848SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
297a8116848SStefano Zampini       ierr = ISDestroy(&is);CHKERRQ(ierr);
298a8116848SStefano Zampini       /* local is to extract local submatrix */
299a8116848SStefano Zampini       ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis);CHKERRQ(ierr);
300a8116848SStefano Zampini       ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr);
301a8116848SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
302a8116848SStefano Zampini     }
303a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
304a8116848SStefano Zampini   } else {
305a8116848SStefano Zampini     ierr = MatISGetLocalMat(*newmat,&newlocmat);CHKERRQ(ierr);
306a8116848SStefano Zampini   }
307a8116848SStefano Zampini   ierr = MatISGetLocalMat(mat,&locmat);CHKERRQ(ierr);
308a8116848SStefano Zampini   newmatis = (Mat_IS*)(*newmat)->data;
309a8116848SStefano Zampini   ierr = MatGetSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat);CHKERRQ(ierr);
310a8116848SStefano Zampini   if (scall == MAT_INITIAL_MATRIX) {
311a8116848SStefano Zampini     ierr = MatISSetLocalMat(*newmat,newlocmat);CHKERRQ(ierr);
312a8116848SStefano Zampini     ierr = MatDestroy(&newlocmat);CHKERRQ(ierr);
313a8116848SStefano Zampini   }
314a8116848SStefano Zampini   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
315a8116848SStefano Zampini   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
316a8116848SStefano Zampini   PetscFunctionReturn(0);
317a8116848SStefano Zampini }
318a8116848SStefano Zampini 
319a8116848SStefano Zampini #undef __FUNCT__
3202b404112SStefano Zampini #define __FUNCT__ "MatCopy_IS"
321a8116848SStefano Zampini static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str)
3222b404112SStefano Zampini {
3232b404112SStefano Zampini   Mat_IS         *a = (Mat_IS*)A->data,*b;
3242b404112SStefano Zampini   PetscBool      ismatis;
3252b404112SStefano Zampini   PetscErrorCode ierr;
3262b404112SStefano Zampini 
3272b404112SStefano Zampini   PetscFunctionBegin;
3282b404112SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr);
3292b404112SStefano Zampini   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented");
3302b404112SStefano Zampini   b = (Mat_IS*)B->data;
3312b404112SStefano Zampini   ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr);
3322b404112SStefano Zampini   PetscFunctionReturn(0);
3332b404112SStefano Zampini }
3342b404112SStefano Zampini 
3352b404112SStefano Zampini #undef __FUNCT__
3366bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS"
337a8116848SStefano Zampini static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool  *missing,PetscInt *d)
3386bd84002SStefano Zampini {
339527b2640SStefano Zampini   Vec               v;
340527b2640SStefano Zampini   const PetscScalar *array;
341527b2640SStefano Zampini   PetscInt          i,n;
3426bd84002SStefano Zampini   PetscErrorCode    ierr;
3436bd84002SStefano Zampini 
3446bd84002SStefano Zampini   PetscFunctionBegin;
345527b2640SStefano Zampini   *missing = PETSC_FALSE;
346527b2640SStefano Zampini   ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr);
347527b2640SStefano Zampini   ierr = MatGetDiagonal(A,v);CHKERRQ(ierr);
348527b2640SStefano Zampini   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
349527b2640SStefano Zampini   ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr);
350527b2640SStefano Zampini   for (i=0;i<n;i++) if (array[i] == 0.) break;
351527b2640SStefano Zampini   ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr);
352527b2640SStefano Zampini   ierr = VecDestroy(&v);CHKERRQ(ierr);
353527b2640SStefano Zampini   if (i != n) *missing = PETSC_TRUE;
354527b2640SStefano Zampini   if (d) {
355527b2640SStefano Zampini     *d = -1;
356527b2640SStefano Zampini     if (*missing) {
357527b2640SStefano Zampini       PetscInt rstart;
358527b2640SStefano Zampini       ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
359527b2640SStefano Zampini       *d = i+rstart;
360527b2640SStefano Zampini     }
361527b2640SStefano Zampini   }
3626bd84002SStefano Zampini   PetscFunctionReturn(0);
3636bd84002SStefano Zampini }
3646bd84002SStefano Zampini 
3656bd84002SStefano Zampini #undef __FUNCT__
366cf0a3239SStefano Zampini #define __FUNCT__ "MatISSetUpSF_IS"
367cf0a3239SStefano Zampini static PetscErrorCode MatISSetUpSF_IS(Mat B)
36828f4e0baSStefano Zampini {
36928f4e0baSStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
37028f4e0baSStefano Zampini   const PetscInt *gidxs;
3714f2d7cafSStefano Zampini   PetscInt       nleaves;
37228f4e0baSStefano Zampini   PetscErrorCode ierr;
37328f4e0baSStefano Zampini 
37428f4e0baSStefano Zampini   PetscFunctionBegin;
3754f2d7cafSStefano Zampini   if (matis->sf) PetscFunctionReturn(0);
37628f4e0baSStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr);
3773bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
3784f2d7cafSStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(B->rmap->mapping,&nleaves);CHKERRQ(ierr);
3794f2d7cafSStefano Zampini   ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
3803bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
3814f2d7cafSStefano Zampini   ierr = PetscMalloc2(matis->sf->nroots,&matis->sf_rootdata,matis->sf->nleaves,&matis->sf_leafdata);CHKERRQ(ierr);
382a8116848SStefano Zampini   if (B->rmap->mapping != B->cmap->mapping) { /* setup SF for columns */
3833d996552SStefano Zampini     ierr = ISLocalToGlobalMappingGetSize(B->cmap->mapping,&nleaves);CHKERRQ(ierr);
384a8116848SStefano Zampini     ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf);CHKERRQ(ierr);
385a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr);
3863d996552SStefano Zampini     ierr = PetscSFSetGraphLayout(matis->csf,B->cmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
387a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr);
3883d996552SStefano Zampini     ierr = PetscMalloc2(matis->csf->nroots,&matis->csf_rootdata,matis->csf->nleaves,&matis->csf_leafdata);CHKERRQ(ierr);
389a8116848SStefano Zampini   } else {
390a8116848SStefano Zampini     matis->csf = matis->sf;
391a8116848SStefano Zampini     matis->csf_leafdata = matis->sf_leafdata;
392a8116848SStefano Zampini     matis->csf_rootdata = matis->sf_rootdata;
393a8116848SStefano Zampini   }
39428f4e0baSStefano Zampini   PetscFunctionReturn(0);
39528f4e0baSStefano Zampini }
3962e1947a5SStefano Zampini 
3972e1947a5SStefano Zampini #undef __FUNCT__
3982e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation"
399eb82efa4SStefano Zampini /*@
400a88811baSStefano Zampini    MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix.
401a88811baSStefano Zampini 
402a88811baSStefano Zampini    Collective on MPI_Comm
403a88811baSStefano Zampini 
404a88811baSStefano Zampini    Input Parameters:
405a88811baSStefano Zampini +  B - the matrix
406a88811baSStefano Zampini .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
407a88811baSStefano Zampini            (same value is used for all local rows)
408a88811baSStefano Zampini .  d_nnz - array containing the number of nonzeros in the various rows of the
409a88811baSStefano Zampini            DIAGONAL portion of the local submatrix (possibly different for each row)
410a88811baSStefano Zampini            or NULL, if d_nz is used to specify the nonzero structure.
411a88811baSStefano Zampini            The size of this array is equal to the number of local rows, i.e 'm'.
412a88811baSStefano Zampini            For matrices that will be factored, you must leave room for (and set)
413a88811baSStefano Zampini            the diagonal entry even if it is zero.
414a88811baSStefano Zampini .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
415a88811baSStefano Zampini            submatrix (same value is used for all local rows).
416a88811baSStefano Zampini -  o_nnz - array containing the number of nonzeros in the various rows of the
417a88811baSStefano Zampini            OFF-DIAGONAL portion of the local submatrix (possibly different for
418a88811baSStefano Zampini            each row) or NULL, if o_nz is used to specify the nonzero
419a88811baSStefano Zampini            structure. The size of this array is equal to the number
420a88811baSStefano Zampini            of local rows, i.e 'm'.
421a88811baSStefano Zampini 
422a88811baSStefano Zampini    If the *_nnz parameter is given then the *_nz parameter is ignored
423a88811baSStefano Zampini 
424a88811baSStefano Zampini    Level: intermediate
425a88811baSStefano Zampini 
426a88811baSStefano Zampini    Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition
427a88811baSStefano Zampini           from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local
428a88811baSStefano Zampini           matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects.
429a88811baSStefano Zampini 
430a88811baSStefano Zampini .keywords: matrix
431a88811baSStefano Zampini 
4323c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS
433a88811baSStefano Zampini @*/
4342e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
4352e1947a5SStefano Zampini {
4362e1947a5SStefano Zampini   PetscErrorCode ierr;
4372e1947a5SStefano Zampini 
4382e1947a5SStefano Zampini   PetscFunctionBegin;
4392e1947a5SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
4402e1947a5SStefano Zampini   PetscValidType(B,1);
4412e1947a5SStefano Zampini   ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
4422e1947a5SStefano Zampini   PetscFunctionReturn(0);
4432e1947a5SStefano Zampini }
4442e1947a5SStefano Zampini 
4452e1947a5SStefano Zampini #undef __FUNCT__
4462e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS"
4477230de76SStefano Zampini static PetscErrorCode  MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
4482e1947a5SStefano Zampini {
4492e1947a5SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
45028f4e0baSStefano Zampini   PetscInt       bs,i,nlocalcols;
4512e1947a5SStefano Zampini   PetscErrorCode ierr;
4522e1947a5SStefano Zampini 
4532e1947a5SStefano Zampini   PetscFunctionBegin;
4546c4ed002SBarry Smith   if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
455cf0a3239SStefano Zampini   ierr = MatISSetUpSF(B);CHKERRQ(ierr);
4564f2d7cafSStefano Zampini 
4574f2d7cafSStefano Zampini   if (!d_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nz;
4584f2d7cafSStefano Zampini   else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nnz[i];
4594f2d7cafSStefano Zampini 
4604f2d7cafSStefano Zampini   if (!o_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nz;
4614f2d7cafSStefano Zampini   else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nnz[i];
4624f2d7cafSStefano Zampini 
46328f4e0baSStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
46428f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr);
46528f4e0baSStefano Zampini   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
46628f4e0baSStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
4674f2d7cafSStefano Zampini 
4684f2d7cafSStefano Zampini   for (i=0;i<matis->sf->nleaves;i++) matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols);
46928f4e0baSStefano Zampini   ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr);
4704f2d7cafSStefano Zampini 
4714f2d7cafSStefano Zampini   for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs;
47228f4e0baSStefano Zampini   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
4734f2d7cafSStefano Zampini 
4744f2d7cafSStefano Zampini   for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i]-i;
47528f4e0baSStefano Zampini   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
4762e1947a5SStefano Zampini   PetscFunctionReturn(0);
4772e1947a5SStefano Zampini }
478b4319ba4SBarry Smith 
479b4319ba4SBarry Smith #undef __FUNCT__
4803927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private"
4813927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
4823927de2eSStefano Zampini {
4833927de2eSStefano Zampini   Mat_IS          *matis = (Mat_IS*)(A->data);
4843927de2eSStefano Zampini   PetscInt        *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership;
485ecf5a873SStefano Zampini   const PetscInt  *global_indices_r,*global_indices_c;
4863927de2eSStefano Zampini   PetscInt        i,j,bs,rows,cols;
4873927de2eSStefano Zampini   PetscInt        lrows,lcols;
4883927de2eSStefano Zampini   PetscInt        local_rows,local_cols;
4893927de2eSStefano Zampini   PetscMPIInt     nsubdomains;
4903927de2eSStefano Zampini   PetscBool       isdense,issbaij;
4913927de2eSStefano Zampini   PetscErrorCode  ierr;
4923927de2eSStefano Zampini 
4933927de2eSStefano Zampini   PetscFunctionBegin;
4943927de2eSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr);
4953927de2eSStefano Zampini   ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr);
4963927de2eSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
4973927de2eSStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
4983927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
4993927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
500ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
501ecf5a873SStefano Zampini   if (A->rmap->mapping != A->cmap->mapping) {
5027230de76SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
503ecf5a873SStefano Zampini   } else {
504ecf5a873SStefano Zampini     global_indices_c = global_indices_r;
505ecf5a873SStefano Zampini   }
506ecf5a873SStefano Zampini 
5073927de2eSStefano Zampini   if (issbaij) {
5083927de2eSStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
5093927de2eSStefano Zampini   }
5103927de2eSStefano Zampini   /*
511ecf5a873SStefano Zampini      An SF reduce is needed to sum up properly on shared rows.
5123927de2eSStefano Zampini      Note that generally preallocation is not exact, since it overestimates nonzeros
5133927de2eSStefano Zampini   */
514cf0a3239SStefano Zampini   ierr = MatISSetUpSF(A);CHKERRQ(ierr);
5153927de2eSStefano Zampini   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
5163927de2eSStefano Zampini   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
5173927de2eSStefano Zampini   /* All processes need to compute entire row ownership */
5183927de2eSStefano Zampini   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
5193927de2eSStefano Zampini   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
5203927de2eSStefano Zampini   for (i=0;i<nsubdomains;i++) {
5213927de2eSStefano Zampini     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
5223927de2eSStefano Zampini       row_ownership[j] = i;
5233927de2eSStefano Zampini     }
5243927de2eSStefano Zampini   }
5257230de76SStefano Zampini   ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
5263927de2eSStefano Zampini 
5273927de2eSStefano Zampini   /*
5283927de2eSStefano Zampini      my_dnz and my_onz contains exact contribution to preallocation from each local mat
5293927de2eSStefano Zampini      then, they will be summed up properly. This way, preallocation is always sufficient
5303927de2eSStefano Zampini   */
5313927de2eSStefano Zampini   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
5323927de2eSStefano Zampini   /* preallocation as a MATAIJ */
5333927de2eSStefano Zampini   if (isdense) { /* special case for dense local matrices */
5343927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
535ecf5a873SStefano Zampini       PetscInt index_row = global_indices_r[i];
5363927de2eSStefano Zampini       for (j=i;j<local_rows;j++) {
5373927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
538ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[j];
5393927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
5403927de2eSStefano Zampini           my_dnz[i] += 1;
5413927de2eSStefano Zampini         } else { /* offdiag block */
5423927de2eSStefano Zampini           my_onz[i] += 1;
5433927de2eSStefano Zampini         }
5443927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
5453927de2eSStefano Zampini         if (i != j) {
5463927de2eSStefano Zampini           owner = row_ownership[index_col];
5473927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
5483927de2eSStefano Zampini             my_dnz[j] += 1;
5493927de2eSStefano Zampini           } else {
5503927de2eSStefano Zampini             my_onz[j] += 1;
5513927de2eSStefano Zampini           }
5523927de2eSStefano Zampini         }
5533927de2eSStefano Zampini       }
5543927de2eSStefano Zampini     }
555bb1015c3SStefano Zampini   } else if (matis->A->ops->getrowij) {
556bb1015c3SStefano Zampini     const PetscInt *ii,*jj,*jptr;
557bb1015c3SStefano Zampini     PetscBool      done;
558bb1015c3SStefano Zampini     ierr = MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr);
559bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__);
560bb1015c3SStefano Zampini     jptr = jj;
561bb1015c3SStefano Zampini     for (i=0;i<local_rows;i++) {
562bb1015c3SStefano Zampini       PetscInt index_row = global_indices_r[i];
563bb1015c3SStefano Zampini       for (j=0;j<ii[i+1]-ii[i];j++,jptr++) {
564bb1015c3SStefano Zampini         PetscInt owner = row_ownership[index_row];
565bb1015c3SStefano Zampini         PetscInt index_col = global_indices_c[*jptr];
566bb1015c3SStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
567bb1015c3SStefano Zampini           my_dnz[i] += 1;
568bb1015c3SStefano Zampini         } else { /* offdiag block */
569bb1015c3SStefano Zampini           my_onz[i] += 1;
570bb1015c3SStefano Zampini         }
571bb1015c3SStefano Zampini         /* same as before, interchanging rows and cols */
572bb1015c3SStefano Zampini         if (issbaij && index_col != index_row) {
573bb1015c3SStefano Zampini           owner = row_ownership[index_col];
574bb1015c3SStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
575bb1015c3SStefano Zampini             my_dnz[*jptr] += 1;
576bb1015c3SStefano Zampini           } else {
577bb1015c3SStefano Zampini             my_onz[*jptr] += 1;
578bb1015c3SStefano Zampini           }
579bb1015c3SStefano Zampini         }
580bb1015c3SStefano Zampini       }
581bb1015c3SStefano Zampini     }
582bb1015c3SStefano Zampini     ierr = MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr);
583bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__);
584bb1015c3SStefano Zampini   } else { /* loop over rows and use MatGetRow */
5853927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
5863927de2eSStefano Zampini       const PetscInt *cols;
587ecf5a873SStefano Zampini       PetscInt       ncols,index_row = global_indices_r[i];
5883927de2eSStefano Zampini       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
5893927de2eSStefano Zampini       for (j=0;j<ncols;j++) {
5903927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
591ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[cols[j]];
5923927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
5933927de2eSStefano Zampini           my_dnz[i] += 1;
5943927de2eSStefano Zampini         } else { /* offdiag block */
5953927de2eSStefano Zampini           my_onz[i] += 1;
5963927de2eSStefano Zampini         }
5973927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
598d9a9e74cSStefano Zampini         if (issbaij && index_col != index_row) {
5993927de2eSStefano Zampini           owner = row_ownership[index_col];
6003927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
601d9a9e74cSStefano Zampini             my_dnz[cols[j]] += 1;
6023927de2eSStefano Zampini           } else {
603d9a9e74cSStefano Zampini             my_onz[cols[j]] += 1;
6043927de2eSStefano Zampini           }
6053927de2eSStefano Zampini         }
6063927de2eSStefano Zampini       }
6073927de2eSStefano Zampini       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
6083927de2eSStefano Zampini     }
6093927de2eSStefano Zampini   }
610ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
611ecf5a873SStefano Zampini   if (global_indices_c != global_indices_r) {
6127230de76SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
613ecf5a873SStefano Zampini   }
6143927de2eSStefano Zampini   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
615ecf5a873SStefano Zampini 
616ecf5a873SStefano Zampini   /* Reduce my_dnz and my_onz */
6173927de2eSStefano Zampini   if (maxreduce) {
6183927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
6193927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
620bb1015c3SStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
6213927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
6223927de2eSStefano Zampini   } else {
6233927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
6243927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
625bb1015c3SStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
6263927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
6273927de2eSStefano Zampini   }
6283927de2eSStefano Zampini   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
6293927de2eSStefano Zampini 
6303927de2eSStefano Zampini   /* Resize preallocation if overestimated */
6313927de2eSStefano Zampini   for (i=0;i<lrows;i++) {
6323927de2eSStefano Zampini     dnz[i] = PetscMin(dnz[i],lcols);
6333927de2eSStefano Zampini     onz[i] = PetscMin(onz[i],cols-lcols);
6343927de2eSStefano Zampini   }
6353927de2eSStefano Zampini   /* set preallocation */
6363927de2eSStefano Zampini   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
6373927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
6383927de2eSStefano Zampini     dnz[i] = dnz[i*bs]/bs;
6393927de2eSStefano Zampini     onz[i] = onz[i*bs]/bs;
6403927de2eSStefano Zampini   }
6413927de2eSStefano Zampini   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
6423927de2eSStefano Zampini   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
6433927de2eSStefano Zampini   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
6443927de2eSStefano Zampini   if (issbaij) {
6453927de2eSStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
6463927de2eSStefano Zampini   }
6473927de2eSStefano Zampini   PetscFunctionReturn(0);
6483927de2eSStefano Zampini }
6493927de2eSStefano Zampini 
6503927de2eSStefano Zampini #undef __FUNCT__
651b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
6527230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
653b7ce53b6SStefano Zampini {
654b7ce53b6SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
655d9a9e74cSStefano Zampini   Mat            local_mat;
656b7ce53b6SStefano Zampini   /* info on mat */
6573cfa4ea4SStefano Zampini   PetscInt       bs,rows,cols,lrows,lcols;
658b7ce53b6SStefano Zampini   PetscInt       local_rows,local_cols;
659b9ed4604SStefano Zampini   PetscBool      isseqdense,isseqsbaij,isseqaij,isseqbaij;
660b9ed4604SStefano Zampini #if defined (PETSC_USE_DEBUG)
661b9ed4604SStefano Zampini   PetscBool      lb[4],bb[4];
662b9ed4604SStefano Zampini #endif
6637c03b4e8SStefano Zampini   PetscMPIInt    nsubdomains;
664b7ce53b6SStefano Zampini   /* values insertion */
665b7ce53b6SStefano Zampini   PetscScalar    *array;
666b7ce53b6SStefano Zampini   /* work */
667b7ce53b6SStefano Zampini   PetscErrorCode ierr;
668b7ce53b6SStefano Zampini 
669b7ce53b6SStefano Zampini   PetscFunctionBegin;
670b7ce53b6SStefano Zampini   /* get info from mat */
6717c03b4e8SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
6727c03b4e8SStefano Zampini   if (nsubdomains == 1) {
6737c03b4e8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
6747c03b4e8SStefano Zampini       ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr);
6757c03b4e8SStefano Zampini     } else {
6767c03b4e8SStefano Zampini       ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
6777c03b4e8SStefano Zampini     }
6787c03b4e8SStefano Zampini     PetscFunctionReturn(0);
6797c03b4e8SStefano Zampini   }
680b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
681b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
6823cfa4ea4SStefano Zampini   ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr);
683b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
684b9ed4604SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isseqdense);CHKERRQ(ierr);
685686e3a49SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
686b9ed4604SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQBAIJ,&isseqbaij);CHKERRQ(ierr);
687b9ed4604SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr);
688b9ed4604SStefano Zampini   if (!isseqdense && !isseqaij && !isseqbaij && !isseqsbaij) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)(matis->A))->type_name);
689b9ed4604SStefano Zampini #if defined (PETSC_USE_DEBUG)
690b9ed4604SStefano Zampini   lb[0] = isseqdense;
691b9ed4604SStefano Zampini   lb[1] = isseqaij;
692b9ed4604SStefano Zampini   lb[2] = isseqbaij;
693b9ed4604SStefano Zampini   lb[3] = isseqsbaij;
694b9ed4604SStefano Zampini   ierr = MPIU_Allreduce(lb,bb,4,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
695b9ed4604SStefano Zampini   if (!bb[0] && !bb[1] && !bb[2] && !bb[3]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local matrices must have the same type");
696b9ed4604SStefano Zampini #endif
697b7ce53b6SStefano Zampini 
698b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
6993927de2eSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
7003cfa4ea4SStefano Zampini     ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr);
7013927de2eSStefano Zampini     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
702b9ed4604SStefano Zampini     if (!isseqsbaij) {
703b9ed4604SStefano Zampini       ierr = MatSetType(*M,MATAIJ);CHKERRQ(ierr);
704b9ed4604SStefano Zampini     } else {
705b9ed4604SStefano Zampini       ierr = MatSetType(*M,MATSBAIJ);CHKERRQ(ierr);
706b9ed4604SStefano Zampini     }
7073927de2eSStefano Zampini     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
708b7ce53b6SStefano Zampini   } else {
7093cfa4ea4SStefano Zampini     PetscInt mbs,mrows,mcols,mlrows,mlcols;
710b7ce53b6SStefano Zampini     /* some checks */
711b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
712b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
7133cfa4ea4SStefano Zampini     ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr);
7146c4ed002SBarry Smith     if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
7156c4ed002SBarry Smith     if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
7166c4ed002SBarry Smith     if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows);
7176c4ed002SBarry Smith     if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols);
7186c4ed002SBarry Smith     if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
719b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
720b7ce53b6SStefano Zampini   }
721d9a9e74cSStefano Zampini 
722b9ed4604SStefano Zampini   if (isseqsbaij) {
723d9a9e74cSStefano Zampini     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr);
724d9a9e74cSStefano Zampini   } else {
725d9a9e74cSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
726d9a9e74cSStefano Zampini     local_mat = matis->A;
727d9a9e74cSStefano Zampini   }
728686e3a49SStefano Zampini 
729b7ce53b6SStefano Zampini   /* Set values */
730ecf5a873SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr);
731b9ed4604SStefano Zampini   if (isseqdense) { /* special case for dense local matrices */
732ecf5a873SStefano Zampini     PetscInt i,*dummy_rows,*dummy_cols;
733ecf5a873SStefano Zampini 
734ecf5a873SStefano Zampini     if (local_rows != local_cols) {
735ecf5a873SStefano Zampini       ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr);
736ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
737ecf5a873SStefano Zampini       for (i=0;i<local_cols;i++) dummy_cols[i] = i;
738ecf5a873SStefano Zampini     } else {
739ecf5a873SStefano Zampini       ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr);
740ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
741ecf5a873SStefano Zampini       dummy_cols = dummy_rows;
742ecf5a873SStefano Zampini     }
743b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
744d9a9e74cSStefano Zampini     ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr);
745ecf5a873SStefano Zampini     ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr);
746d9a9e74cSStefano Zampini     ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr);
747ecf5a873SStefano Zampini     if (dummy_rows != dummy_cols) {
748ecf5a873SStefano Zampini       ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr);
749ecf5a873SStefano Zampini     } else {
750ecf5a873SStefano Zampini       ierr = PetscFree(dummy_rows);CHKERRQ(ierr);
751ecf5a873SStefano Zampini     }
752686e3a49SStefano Zampini   } else if (isseqaij) {
753ecf5a873SStefano Zampini     PetscInt  i,nvtxs,*xadj,*adjncy;
754686e3a49SStefano Zampini     PetscBool done;
755686e3a49SStefano Zampini 
756d9a9e74cSStefano Zampini     ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
757bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__);
758d9a9e74cSStefano Zampini     ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr);
759686e3a49SStefano Zampini     for (i=0;i<nvtxs;i++) {
760ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr);
761686e3a49SStefano Zampini     }
762d9a9e74cSStefano Zampini     ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
763bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__);
764d9a9e74cSStefano Zampini     ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr);
765686e3a49SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
766ecf5a873SStefano Zampini     PetscInt i;
767c0962df8SStefano Zampini 
768686e3a49SStefano Zampini     for (i=0;i<local_rows;i++) {
769686e3a49SStefano Zampini       PetscInt       j;
770ecf5a873SStefano Zampini       const PetscInt *local_indices_cols;
771686e3a49SStefano Zampini 
772ecf5a873SStefano Zampini       ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
773ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr);
774ecf5a873SStefano Zampini       ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
775686e3a49SStefano Zampini     }
776b7ce53b6SStefano Zampini   }
777b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
778d9a9e74cSStefano Zampini   ierr = MatDestroy(&local_mat);CHKERRQ(ierr);
779b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
780b9ed4604SStefano Zampini   if (isseqdense) {
781b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
782b7ce53b6SStefano Zampini   }
783b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
784b7ce53b6SStefano Zampini }
785b7ce53b6SStefano Zampini 
786b7ce53b6SStefano Zampini #undef __FUNCT__
787b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
788b7ce53b6SStefano Zampini /*@
789b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
790b7ce53b6SStefano Zampini 
791b7ce53b6SStefano Zampini   Input Parameter:
792b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
793b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
794b7ce53b6SStefano Zampini 
795b7ce53b6SStefano Zampini   Output Parameter:
796b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
797b7ce53b6SStefano Zampini 
798b7ce53b6SStefano Zampini   Level: developer
799b7ce53b6SStefano Zampini 
800eb82efa4SStefano Zampini   Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested.
801b7ce53b6SStefano Zampini 
802b7ce53b6SStefano Zampini .seealso: MATIS
803b7ce53b6SStefano Zampini @*/
804b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
805b7ce53b6SStefano Zampini {
806b7ce53b6SStefano Zampini   PetscErrorCode ierr;
807b7ce53b6SStefano Zampini 
808b7ce53b6SStefano Zampini   PetscFunctionBegin;
809b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
810b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
811b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
812b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
813b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
814b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
8156c4ed002SBarry Smith     if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix");
816b7ce53b6SStefano Zampini   }
817b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
818b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
819b7ce53b6SStefano Zampini }
820b7ce53b6SStefano Zampini 
821b7ce53b6SStefano Zampini #undef __FUNCT__
822ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
823ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
824ad6194a2SStefano Zampini {
825ad6194a2SStefano Zampini   PetscErrorCode ierr;
826ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
827ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
828ad6194a2SStefano Zampini   Mat            B,localmat;
829ad6194a2SStefano Zampini 
830ad6194a2SStefano Zampini   PetscFunctionBegin;
831ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
832ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
833ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
834e176bc59SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr);
835ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
836ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
837b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
838ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
839ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
840ad6194a2SStefano Zampini   *newmat = B;
841ad6194a2SStefano Zampini   PetscFunctionReturn(0);
842ad6194a2SStefano Zampini }
843ad6194a2SStefano Zampini 
844ad6194a2SStefano Zampini #undef __FUNCT__
84569796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
846a8116848SStefano Zampini static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
84769796d55SStefano Zampini {
84869796d55SStefano Zampini   PetscErrorCode ierr;
84969796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
85069796d55SStefano Zampini   PetscBool      local_sym;
85169796d55SStefano Zampini 
85269796d55SStefano Zampini   PetscFunctionBegin;
85369796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
854b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
85569796d55SStefano Zampini   PetscFunctionReturn(0);
85669796d55SStefano Zampini }
85769796d55SStefano Zampini 
85869796d55SStefano Zampini #undef __FUNCT__
85969796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
860a8116848SStefano Zampini static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
86169796d55SStefano Zampini {
86269796d55SStefano Zampini   PetscErrorCode ierr;
86369796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
86469796d55SStefano Zampini   PetscBool      local_sym;
86569796d55SStefano Zampini 
86669796d55SStefano Zampini   PetscFunctionBegin;
86769796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
868b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
86969796d55SStefano Zampini   PetscFunctionReturn(0);
87069796d55SStefano Zampini }
87169796d55SStefano Zampini 
87269796d55SStefano Zampini #undef __FUNCT__
873*45471136SStefano Zampini #define __FUNCT__ "MatIsStructurallySymmetric_IS"
874*45471136SStefano Zampini static PetscErrorCode MatIsStructurallySymmetric_IS(Mat A,PetscBool  *flg)
875*45471136SStefano Zampini {
876*45471136SStefano Zampini   PetscErrorCode ierr;
877*45471136SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
878*45471136SStefano Zampini   PetscBool      local_sym;
879*45471136SStefano Zampini 
880*45471136SStefano Zampini   PetscFunctionBegin;
881*45471136SStefano Zampini   if (A->rmap->mapping != A->cmap->mapping) {
882*45471136SStefano Zampini     *flg = PETSC_FALSE;
883*45471136SStefano Zampini     PetscFunctionReturn(0);
884*45471136SStefano Zampini   }
885*45471136SStefano Zampini   ierr = MatIsStructurallySymmetric(matis->A,&local_sym);CHKERRQ(ierr);
886*45471136SStefano Zampini   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
887*45471136SStefano Zampini   PetscFunctionReturn(0);
888*45471136SStefano Zampini }
889*45471136SStefano Zampini 
890*45471136SStefano Zampini #undef __FUNCT__
891b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
892a8116848SStefano Zampini static PetscErrorCode MatDestroy_IS(Mat A)
893b4319ba4SBarry Smith {
894dfbe8321SBarry Smith   PetscErrorCode ierr;
895b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
896b4319ba4SBarry Smith 
897b4319ba4SBarry Smith   PetscFunctionBegin;
8986bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
899e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr);
900e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr);
9016bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
9026bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
9033fd1c9e7SStefano Zampini   ierr = VecDestroy(&b->counter);CHKERRQ(ierr);
904a8116848SStefano Zampini   ierr = ISDestroy(&b->getsub_ris);CHKERRQ(ierr);
905a8116848SStefano Zampini   ierr = ISDestroy(&b->getsub_cis);CHKERRQ(ierr);
906a8116848SStefano Zampini   if (b->sf != b->csf) {
907a8116848SStefano Zampini     ierr = PetscSFDestroy(&b->csf);CHKERRQ(ierr);
908a8116848SStefano Zampini     ierr = PetscFree2(b->csf_rootdata,b->csf_leafdata);CHKERRQ(ierr);
909a8116848SStefano Zampini   }
91028f4e0baSStefano Zampini   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
91128f4e0baSStefano Zampini   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
912bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
913dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
914bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
915b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
916b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
9172e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
918cf0a3239SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",NULL);CHKERRQ(ierr);
919b4319ba4SBarry Smith   PetscFunctionReturn(0);
920b4319ba4SBarry Smith }
921b4319ba4SBarry Smith 
922b4319ba4SBarry Smith #undef __FUNCT__
923b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
924a8116848SStefano Zampini static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
925b4319ba4SBarry Smith {
926dfbe8321SBarry Smith   PetscErrorCode ierr;
927b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
928b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
929b4319ba4SBarry Smith 
930b4319ba4SBarry Smith   PetscFunctionBegin;
931b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
932e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
933e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
934b4319ba4SBarry Smith 
935b4319ba4SBarry Smith   /* multiply the local matrix */
936b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
937b4319ba4SBarry Smith 
938b4319ba4SBarry Smith   /* scatter product back into global memory */
9392dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
940e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
941e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
942b4319ba4SBarry Smith   PetscFunctionReturn(0);
943b4319ba4SBarry Smith }
944b4319ba4SBarry Smith 
945b4319ba4SBarry Smith #undef __FUNCT__
9462e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
947a8116848SStefano Zampini static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
9482e74eeadSLisandro Dalcin {
949650997f4SStefano Zampini   Vec            temp_vec;
9502e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9512e74eeadSLisandro Dalcin 
9522e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
953650997f4SStefano Zampini   if (v3 != v2) {
954650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
955650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
956650997f4SStefano Zampini   } else {
957650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
958650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
959650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
960650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
961650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
962650997f4SStefano Zampini   }
9632e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9642e74eeadSLisandro Dalcin }
9652e74eeadSLisandro Dalcin 
9662e74eeadSLisandro Dalcin #undef __FUNCT__
9672e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
968a8116848SStefano Zampini static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x)
9692e74eeadSLisandro Dalcin {
9702e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9712e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9722e74eeadSLisandro Dalcin 
973e176bc59SStefano Zampini   PetscFunctionBegin;
9742e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
975e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
976e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9772e74eeadSLisandro Dalcin 
9782e74eeadSLisandro Dalcin   /* multiply the local matrix */
979e176bc59SStefano Zampini   ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr);
9802e74eeadSLisandro Dalcin 
9812e74eeadSLisandro Dalcin   /* scatter product back into global vector */
982e176bc59SStefano Zampini   ierr = VecSet(x,0);CHKERRQ(ierr);
983e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
984e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9852e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9862e74eeadSLisandro Dalcin }
9872e74eeadSLisandro Dalcin 
9882e74eeadSLisandro Dalcin #undef __FUNCT__
9892e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
990a8116848SStefano Zampini static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
9912e74eeadSLisandro Dalcin {
992650997f4SStefano Zampini   Vec            temp_vec;
9932e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9942e74eeadSLisandro Dalcin 
9952e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
996650997f4SStefano Zampini   if (v3 != v2) {
997650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
998650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
999650997f4SStefano Zampini   } else {
1000650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
1001650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
1002650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
1003650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
1004650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
1005650997f4SStefano Zampini   }
10062e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
10072e74eeadSLisandro Dalcin }
10082e74eeadSLisandro Dalcin 
10092e74eeadSLisandro Dalcin #undef __FUNCT__
1010b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
1011a8116848SStefano Zampini static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
1012b4319ba4SBarry Smith {
1013b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
1014dfbe8321SBarry Smith   PetscErrorCode ierr;
1015b4319ba4SBarry Smith   PetscViewer    sviewer;
1016ee2491ecSStefano Zampini   PetscBool      isascii,view = PETSC_TRUE;
1017b4319ba4SBarry Smith 
1018b4319ba4SBarry Smith   PetscFunctionBegin;
1019ee2491ecSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
1020ee2491ecSStefano Zampini   if (isascii)  {
1021ee2491ecSStefano Zampini     PetscViewerFormat format;
1022ee2491ecSStefano Zampini 
1023ee2491ecSStefano Zampini     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
1024ee2491ecSStefano Zampini     if (format == PETSC_VIEWER_ASCII_INFO) view = PETSC_FALSE;
1025ee2491ecSStefano Zampini   }
1026ee2491ecSStefano Zampini   if (!view) PetscFunctionReturn(0);
10273f08860eSBarry Smith   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
1028b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
10293f08860eSBarry Smith   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
10306e520ac8SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1031b4319ba4SBarry Smith   PetscFunctionReturn(0);
1032b4319ba4SBarry Smith }
1033b4319ba4SBarry Smith 
1034b4319ba4SBarry Smith #undef __FUNCT__
1035b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
1036a8116848SStefano Zampini static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
1037b4319ba4SBarry Smith {
1038dfbe8321SBarry Smith   PetscErrorCode ierr;
1039e176bc59SStefano Zampini   PetscInt       nr,rbs,nc,cbs;
1040b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1041b4319ba4SBarry Smith   IS             from,to;
1042e176bc59SStefano Zampini   Vec            cglobal,rglobal;
1043b4319ba4SBarry Smith 
1044b4319ba4SBarry Smith   PetscFunctionBegin;
1045784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
1046e176bc59SStefano Zampini   PetscCheckSameComm(A,1,cmapping,3);
10473bbff08aSStefano Zampini   /* Destroy any previous data */
104870cf5478SStefano Zampini   ierr = VecDestroy(&is->x);CHKERRQ(ierr);
104970cf5478SStefano Zampini   ierr = VecDestroy(&is->y);CHKERRQ(ierr);
10503fd1c9e7SStefano Zampini   ierr = VecDestroy(&is->counter);CHKERRQ(ierr);
1051e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr);
1052e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr);
10531c47cb0fSStefano Zampini   ierr = MatDestroy(&is->A);CHKERRQ(ierr);
105428f4e0baSStefano Zampini   ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
105528f4e0baSStefano Zampini   ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
10563bbff08aSStefano Zampini 
10573bbff08aSStefano Zampini   /* Setup Layout and set local to global maps */
1058fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1059fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1060fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
1061fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
1062b4319ba4SBarry Smith 
1063b4319ba4SBarry Smith   /* Create the local matrix A */
1064e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr);
1065e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr);
1066e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr);
1067e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr);
1068f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
1069e176bc59SStefano Zampini   ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr);
1070e176bc59SStefano Zampini   ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr);
1071e176bc59SStefano Zampini   ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr);
1072ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
1073ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
1074b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
1075b4319ba4SBarry Smith 
1076f26d0771SStefano Zampini   if (!is->islocalref) { /* setup scatters and local vectors for MatMult */
1077b4319ba4SBarry Smith     /* Create the local work vectors */
10783bbff08aSStefano Zampini     ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr);
1079b4319ba4SBarry Smith 
1080e176bc59SStefano Zampini     /* setup the global to local scatters */
1081e176bc59SStefano Zampini     ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr);
1082e176bc59SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr);
1083784ac674SJed Brown     ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
1084e176bc59SStefano Zampini     ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr);
1085e176bc59SStefano Zampini     if (rmapping != cmapping) {
1086e176bc59SStefano Zampini       ierr = ISDestroy(&to);CHKERRQ(ierr);
1087e176bc59SStefano Zampini       ierr = ISDestroy(&from);CHKERRQ(ierr);
1088e176bc59SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr);
1089e176bc59SStefano Zampini       ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr);
1090e176bc59SStefano Zampini       ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr);
1091e176bc59SStefano Zampini     } else {
1092e176bc59SStefano Zampini       ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr);
1093e176bc59SStefano Zampini       is->cctx = is->rctx;
1094e176bc59SStefano Zampini     }
10953fd1c9e7SStefano Zampini 
10963fd1c9e7SStefano Zampini     /* interface counter vector (local) */
10973fd1c9e7SStefano Zampini     ierr = VecDuplicate(is->y,&is->counter);CHKERRQ(ierr);
10983fd1c9e7SStefano Zampini     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
10993fd1c9e7SStefano Zampini     ierr = VecScatterBegin(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11003fd1c9e7SStefano Zampini     ierr = VecScatterEnd(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11013fd1c9e7SStefano Zampini     ierr = VecScatterBegin(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11023fd1c9e7SStefano Zampini     ierr = VecScatterEnd(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11033fd1c9e7SStefano Zampini 
11043fd1c9e7SStefano Zampini     /* free workspace */
1105e176bc59SStefano Zampini     ierr = VecDestroy(&rglobal);CHKERRQ(ierr);
1106e176bc59SStefano Zampini     ierr = VecDestroy(&cglobal);CHKERRQ(ierr);
11076bf464f9SBarry Smith     ierr = ISDestroy(&to);CHKERRQ(ierr);
11086bf464f9SBarry Smith     ierr = ISDestroy(&from);CHKERRQ(ierr);
1109f26d0771SStefano Zampini   }
11109c0b00dbSStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
1111b4319ba4SBarry Smith   PetscFunctionReturn(0);
1112b4319ba4SBarry Smith }
1113b4319ba4SBarry Smith 
11142e74eeadSLisandro Dalcin #undef __FUNCT__
11152e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
1116a8116848SStefano Zampini static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
11172e74eeadSLisandro Dalcin {
11182e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
11192e74eeadSLisandro Dalcin   PetscErrorCode ierr;
112097563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
112197563a80SStefano Zampini   PetscInt       i,zm,zn;
112297563a80SStefano Zampini #endif
1123f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
11242e74eeadSLisandro Dalcin 
11252e74eeadSLisandro Dalcin   PetscFunctionBegin;
11262e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
1127f26d0771SStefano Zampini   if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n);
112897563a80SStefano Zampini   /* count negative indices */
112997563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
113097563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
11312e74eeadSLisandro Dalcin #endif
113297563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
113397563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
113497563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
113597563a80SStefano Zampini   /* count negative indices (should be the same as before) */
113697563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
113797563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
113897563a80SStefano Zampini   if (zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the row indices can not be mapped! Maybe you should not use MATIS");
113997563a80SStefano Zampini   if (zn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the column indices can not be mapped! Maybe you should not use MATIS");
114097563a80SStefano Zampini #endif
11412e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
11422e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
11432e74eeadSLisandro Dalcin }
11442e74eeadSLisandro Dalcin 
1145b4319ba4SBarry Smith #undef __FUNCT__
114697563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS"
1147a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
114897563a80SStefano Zampini {
114997563a80SStefano Zampini   Mat_IS         *is = (Mat_IS*)mat->data;
115097563a80SStefano Zampini   PetscErrorCode ierr;
115197563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
115297563a80SStefano Zampini   PetscInt       i,zm,zn;
115397563a80SStefano Zampini #endif
1154f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
115597563a80SStefano Zampini 
115697563a80SStefano Zampini   PetscFunctionBegin;
115797563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
1158f26d0771SStefano Zampini   if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column block indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n);
115997563a80SStefano Zampini   /* count negative indices */
116097563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
116197563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
116297563a80SStefano Zampini #endif
116397563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
116497563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
116597563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
116697563a80SStefano Zampini   /* count negative indices (should be the same as before) */
116797563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
116897563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
116997563a80SStefano Zampini   if (zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the row indices can not be mapped! Maybe you should not use MATIS");
117097563a80SStefano Zampini   if (zn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the column indices can not be mapped! Maybe you should not use MATIS");
117197563a80SStefano Zampini #endif
117297563a80SStefano Zampini   ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
117397563a80SStefano Zampini   PetscFunctionReturn(0);
117497563a80SStefano Zampini }
117597563a80SStefano Zampini 
117697563a80SStefano Zampini #undef __FUNCT__
1177b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
1178a8116848SStefano Zampini static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
1179b4319ba4SBarry Smith {
1180dfbe8321SBarry Smith   PetscErrorCode ierr;
1181b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1182b4319ba4SBarry Smith 
1183b4319ba4SBarry Smith   PetscFunctionBegin;
1184b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
1185b4319ba4SBarry Smith   PetscFunctionReturn(0);
1186b4319ba4SBarry Smith }
1187b4319ba4SBarry Smith 
1188b4319ba4SBarry Smith #undef __FUNCT__
1189f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
1190a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
1191f0006bf2SLisandro Dalcin {
1192f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
1193f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
1194f0006bf2SLisandro Dalcin 
1195f0006bf2SLisandro Dalcin   PetscFunctionBegin;
1196f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
1197f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
1198f0006bf2SLisandro Dalcin }
1199f0006bf2SLisandro Dalcin 
1200f0006bf2SLisandro Dalcin #undef __FUNCT__
1201f0ae7da4SStefano Zampini #define __FUNCT__ "MatISZeroRowsColumnsLocal_Private"
1202f0ae7da4SStefano Zampini static PetscErrorCode MatISZeroRowsColumnsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,PetscBool columns)
12032e74eeadSLisandro Dalcin {
1204f0ae7da4SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
12052e74eeadSLisandro Dalcin   PetscErrorCode ierr;
12062e74eeadSLisandro Dalcin 
12072e74eeadSLisandro Dalcin   PetscFunctionBegin;
1208f0ae7da4SStefano Zampini   if (!n) {
1209f0ae7da4SStefano Zampini     is->pure_neumann = PETSC_TRUE;
1210f0ae7da4SStefano Zampini   } else {
1211f0ae7da4SStefano Zampini     PetscInt i;
1212f0ae7da4SStefano Zampini     is->pure_neumann = PETSC_FALSE;
1213f0ae7da4SStefano Zampini 
1214f0ae7da4SStefano Zampini     if (columns) {
1215f0ae7da4SStefano Zampini       ierr = MatZeroRowsColumns(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
1216f0ae7da4SStefano Zampini     } else {
1217f0ae7da4SStefano Zampini       ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
12182e74eeadSLisandro Dalcin     }
1219f0ae7da4SStefano Zampini     if (diag != 0.) {
1220f0ae7da4SStefano Zampini       const PetscScalar *array;
1221f0ae7da4SStefano Zampini       ierr = VecGetArrayRead(is->counter,&array);CHKERRQ(ierr);
1222f0ae7da4SStefano Zampini       for (i=0; i<n; i++) {
1223f0ae7da4SStefano Zampini         ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
1224f0ae7da4SStefano Zampini       }
1225f0ae7da4SStefano Zampini       ierr = VecRestoreArrayRead(is->counter,&array);CHKERRQ(ierr);
1226f0ae7da4SStefano Zampini     }
1227f0ae7da4SStefano Zampini     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1228f0ae7da4SStefano Zampini     ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1229f0ae7da4SStefano Zampini   }
12302e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
12312e74eeadSLisandro Dalcin }
12322e74eeadSLisandro Dalcin 
12332e74eeadSLisandro Dalcin #undef __FUNCT__
1234f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_Private_IS"
1235f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Private_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b,PetscBool columns)
1236b4319ba4SBarry Smith {
12376e520ac8SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
12386e520ac8SStefano Zampini   PetscInt       nr,nl,len,i;
12396e520ac8SStefano Zampini   PetscInt       *lrows;
1240dfbe8321SBarry Smith   PetscErrorCode ierr;
1241b4319ba4SBarry Smith 
1242b4319ba4SBarry Smith   PetscFunctionBegin;
1243f0ae7da4SStefano Zampini #if defined(PETSC_USE_DEBUG)
1244f0ae7da4SStefano Zampini   if (columns || diag != 0. || (x && b)) {
1245f0ae7da4SStefano Zampini     PetscBool cong;
1246f0ae7da4SStefano Zampini     ierr = PetscLayoutCompare(A->rmap,A->cmap,&cong);CHKERRQ(ierr);
1247f0ae7da4SStefano Zampini     if (!cong && columns) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Columns can be zeroed if and only if A->rmap and A->cmap are congruent for MATIS");
1248f0ae7da4SStefano Zampini     if (!cong && diag != 0.) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Nonzero diagonal value supported if and only if A->rmap and A->cmap are congruent for MATIS");
1249f0ae7da4SStefano Zampini     if (!cong && x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"A->rmap and A->cmap need to be congruent");
1250b4319ba4SBarry Smith   }
1251f0ae7da4SStefano Zampini #endif
12526e520ac8SStefano Zampini   /* get locally owned rows */
1253f0ae7da4SStefano Zampini   ierr = PetscLayoutMapLocal_Private(A->rmap,n,rows,&len,&lrows,NULL);CHKERRQ(ierr);
12546e520ac8SStefano Zampini   /* fix right hand side if needed */
12556e520ac8SStefano Zampini   if (x && b) {
12566e520ac8SStefano Zampini     const PetscScalar *xx;
12576e520ac8SStefano Zampini     PetscScalar       *bb;
12582205254eSKarl Rupp 
12596e520ac8SStefano Zampini     ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr);
12606e520ac8SStefano Zampini     ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
12616e520ac8SStefano Zampini     for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]];
12626e520ac8SStefano Zampini     ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr);
12636e520ac8SStefano Zampini     ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr);
1264b4319ba4SBarry Smith   }
12656e520ac8SStefano Zampini   /* get rows associated to the local matrices */
12663d996552SStefano Zampini   ierr = MatISSetUpSF(A);CHKERRQ(ierr);
12676e520ac8SStefano Zampini   ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr);
12686e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr);
12696e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
12706e520ac8SStefano Zampini   for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1;
12716e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
12726e520ac8SStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
12736e520ac8SStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
12746e520ac8SStefano Zampini   ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr);
12756e520ac8SStefano Zampini   for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i;
1276f0ae7da4SStefano Zampini   ierr = MatISZeroRowsColumnsLocal_Private(A,nr,lrows,diag,columns);CHKERRQ(ierr);
12776e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
1278b4319ba4SBarry Smith   PetscFunctionReturn(0);
1279b4319ba4SBarry Smith }
1280b4319ba4SBarry Smith 
1281b4319ba4SBarry Smith #undef __FUNCT__
1282f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRows_IS"
1283f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1284b4319ba4SBarry Smith {
1285b4319ba4SBarry Smith   PetscErrorCode ierr;
1286b4319ba4SBarry Smith 
1287b4319ba4SBarry Smith   PetscFunctionBegin;
1288f0ae7da4SStefano Zampini   ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_FALSE);CHKERRQ(ierr);
1289f0ae7da4SStefano Zampini   PetscFunctionReturn(0);
1290f0ae7da4SStefano Zampini }
1291b4319ba4SBarry Smith 
1292f0ae7da4SStefano Zampini #undef __FUNCT__
1293f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_IS"
1294f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1295f0ae7da4SStefano Zampini {
1296f0ae7da4SStefano Zampini   PetscErrorCode ierr;
1297f0ae7da4SStefano Zampini 
1298f0ae7da4SStefano Zampini   PetscFunctionBegin;
1299f0ae7da4SStefano Zampini   ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_TRUE);CHKERRQ(ierr);
1300b4319ba4SBarry Smith   PetscFunctionReturn(0);
1301b4319ba4SBarry Smith }
1302b4319ba4SBarry Smith 
1303b4319ba4SBarry Smith #undef __FUNCT__
1304b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
1305a8116848SStefano Zampini static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
1306b4319ba4SBarry Smith {
1307b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1308dfbe8321SBarry Smith   PetscErrorCode ierr;
1309b4319ba4SBarry Smith 
1310b4319ba4SBarry Smith   PetscFunctionBegin;
1311b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
1312b4319ba4SBarry Smith   PetscFunctionReturn(0);
1313b4319ba4SBarry Smith }
1314b4319ba4SBarry Smith 
1315b4319ba4SBarry Smith #undef __FUNCT__
1316b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
1317a8116848SStefano Zampini static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
1318b4319ba4SBarry Smith {
1319b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1320dfbe8321SBarry Smith   PetscErrorCode ierr;
1321b4319ba4SBarry Smith 
1322b4319ba4SBarry Smith   PetscFunctionBegin;
1323b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
1324b4319ba4SBarry Smith   PetscFunctionReturn(0);
1325b4319ba4SBarry Smith }
1326b4319ba4SBarry Smith 
1327b4319ba4SBarry Smith #undef __FUNCT__
1328b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
1329a8116848SStefano Zampini static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
1330b4319ba4SBarry Smith {
1331b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
1332b4319ba4SBarry Smith 
1333b4319ba4SBarry Smith   PetscFunctionBegin;
1334b4319ba4SBarry Smith   *local = is->A;
1335b4319ba4SBarry Smith   PetscFunctionReturn(0);
1336b4319ba4SBarry Smith }
1337b4319ba4SBarry Smith 
1338b4319ba4SBarry Smith #undef __FUNCT__
1339b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
1340b4319ba4SBarry Smith /*@
1341b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
1342b4319ba4SBarry Smith 
1343b4319ba4SBarry Smith   Input Parameter:
1344b4319ba4SBarry Smith .  mat - the matrix
1345b4319ba4SBarry Smith 
1346b4319ba4SBarry Smith   Output Parameter:
1347eb82efa4SStefano Zampini .  local - the local matrix
1348b4319ba4SBarry Smith 
1349b4319ba4SBarry Smith   Level: advanced
1350b4319ba4SBarry Smith 
1351b4319ba4SBarry Smith   Notes:
1352b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
1353b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
1354b4319ba4SBarry Smith   of the MatSetValues() operation.
1355b4319ba4SBarry Smith 
1356b4319ba4SBarry Smith .seealso: MATIS
1357b4319ba4SBarry Smith @*/
13587087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
1359b4319ba4SBarry Smith {
13604ac538c5SBarry Smith   PetscErrorCode ierr;
1361b4319ba4SBarry Smith 
1362b4319ba4SBarry Smith   PetscFunctionBegin;
13630700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1364b4319ba4SBarry Smith   PetscValidPointer(local,2);
13654ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
1366b4319ba4SBarry Smith   PetscFunctionReturn(0);
1367b4319ba4SBarry Smith }
1368b4319ba4SBarry Smith 
13693b03a366Sstefano_zampini #undef __FUNCT__
13703b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
1371a8116848SStefano Zampini static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
13723b03a366Sstefano_zampini {
13733b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
13743b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
13753b03a366Sstefano_zampini   PetscErrorCode ierr;
13763b03a366Sstefano_zampini 
13773b03a366Sstefano_zampini   PetscFunctionBegin;
13784e4c7dbeSStefano Zampini   if (is->A) {
13793b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
13803b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
1381f0ae7da4SStefano Zampini     if (orows != nrows || ocols != ncols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local MATIS matrix should be of size %Dx%D (you passed a %Dx%D matrix)",orows,ocols,nrows,ncols);
13824e4c7dbeSStefano Zampini   }
13833b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
13843b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
13853b03a366Sstefano_zampini   is->A = local;
13863b03a366Sstefano_zampini   PetscFunctionReturn(0);
13873b03a366Sstefano_zampini }
13883b03a366Sstefano_zampini 
13893b03a366Sstefano_zampini #undef __FUNCT__
13903b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
13913b03a366Sstefano_zampini /*@
1392eb82efa4SStefano Zampini     MatISSetLocalMat - Replace the local matrix stored inside a MATIS object.
13933b03a366Sstefano_zampini 
13943b03a366Sstefano_zampini   Input Parameter:
13953b03a366Sstefano_zampini .  mat - the matrix
1396eb82efa4SStefano Zampini .  local - the local matrix
13973b03a366Sstefano_zampini 
13983b03a366Sstefano_zampini   Output Parameter:
13993b03a366Sstefano_zampini 
14003b03a366Sstefano_zampini   Level: advanced
14013b03a366Sstefano_zampini 
14023b03a366Sstefano_zampini   Notes:
14033b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
14043b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
14053b03a366Sstefano_zampini 
14063b03a366Sstefano_zampini .seealso: MATIS
14073b03a366Sstefano_zampini @*/
14083b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
14093b03a366Sstefano_zampini {
14103b03a366Sstefano_zampini   PetscErrorCode ierr;
14113b03a366Sstefano_zampini 
14123b03a366Sstefano_zampini   PetscFunctionBegin;
14133b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1414b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
14153b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
14163b03a366Sstefano_zampini   PetscFunctionReturn(0);
14173b03a366Sstefano_zampini }
14183b03a366Sstefano_zampini 
14196726f965SBarry Smith #undef __FUNCT__
14206726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
1421a8116848SStefano Zampini static PetscErrorCode MatZeroEntries_IS(Mat A)
14226726f965SBarry Smith {
14236726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
14246726f965SBarry Smith   PetscErrorCode ierr;
14256726f965SBarry Smith 
14266726f965SBarry Smith   PetscFunctionBegin;
14276726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
14286726f965SBarry Smith   PetscFunctionReturn(0);
14296726f965SBarry Smith }
14306726f965SBarry Smith 
14316726f965SBarry Smith #undef __FUNCT__
14322e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
1433a8116848SStefano Zampini static PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
14342e74eeadSLisandro Dalcin {
14352e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
14362e74eeadSLisandro Dalcin   PetscErrorCode ierr;
14372e74eeadSLisandro Dalcin 
14382e74eeadSLisandro Dalcin   PetscFunctionBegin;
14392e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
14402e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
14412e74eeadSLisandro Dalcin }
14422e74eeadSLisandro Dalcin 
14432e74eeadSLisandro Dalcin #undef __FUNCT__
14442e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
1445a8116848SStefano Zampini static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
14462e74eeadSLisandro Dalcin {
14472e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
14482e74eeadSLisandro Dalcin   PetscErrorCode ierr;
14492e74eeadSLisandro Dalcin 
14502e74eeadSLisandro Dalcin   PetscFunctionBegin;
14512e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
1452e176bc59SStefano Zampini   ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr);
14532e74eeadSLisandro Dalcin 
14542e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
14552e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
1456e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1457e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14582e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
14592e74eeadSLisandro Dalcin }
14602e74eeadSLisandro Dalcin 
14612e74eeadSLisandro Dalcin #undef __FUNCT__
14626726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
1463a8116848SStefano Zampini static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
14646726f965SBarry Smith {
14656726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
14666726f965SBarry Smith   PetscErrorCode ierr;
14676726f965SBarry Smith 
14686726f965SBarry Smith   PetscFunctionBegin;
14694e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
14706726f965SBarry Smith   PetscFunctionReturn(0);
14716726f965SBarry Smith }
14726726f965SBarry Smith 
1473284134d9SBarry Smith #undef __FUNCT__
1474f26d0771SStefano Zampini #define __FUNCT__ "MatAXPY_IS"
1475f26d0771SStefano Zampini static PetscErrorCode MatAXPY_IS(Mat Y,PetscScalar a,Mat X,MatStructure str)
1476f26d0771SStefano Zampini {
1477f26d0771SStefano Zampini   Mat_IS         *y = (Mat_IS*)Y->data;
1478f26d0771SStefano Zampini   Mat_IS         *x;
1479f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1480f26d0771SStefano Zampini   PetscBool      ismatis;
1481f26d0771SStefano Zampini #endif
1482f26d0771SStefano Zampini   PetscErrorCode ierr;
1483f26d0771SStefano Zampini 
1484f26d0771SStefano Zampini   PetscFunctionBegin;
1485f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1486f26d0771SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)X,MATIS,&ismatis);CHKERRQ(ierr);
1487f26d0771SStefano Zampini   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS");
1488f26d0771SStefano Zampini #endif
1489f26d0771SStefano Zampini   x = (Mat_IS*)X->data;
1490f26d0771SStefano Zampini   ierr = MatAXPY(y->A,a,x->A,str);CHKERRQ(ierr);
1491f26d0771SStefano Zampini   PetscFunctionReturn(0);
1492f26d0771SStefano Zampini }
1493f26d0771SStefano Zampini 
1494f26d0771SStefano Zampini #undef __FUNCT__
1495f26d0771SStefano Zampini #define __FUNCT__ "MatGetLocalSubMatrix_IS"
1496f26d0771SStefano Zampini static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat)
1497f26d0771SStefano Zampini {
1498f26d0771SStefano Zampini   Mat                    lA;
1499f26d0771SStefano Zampini   Mat_IS                 *matis;
1500f26d0771SStefano Zampini   ISLocalToGlobalMapping rl2g,cl2g;
1501f26d0771SStefano Zampini   IS                     is;
1502f26d0771SStefano Zampini   const PetscInt         *rg,*rl;
1503f26d0771SStefano Zampini   PetscInt               nrg;
1504f26d0771SStefano Zampini   PetscInt               N,M,nrl,i,*idxs;
1505f26d0771SStefano Zampini   PetscErrorCode         ierr;
1506f26d0771SStefano Zampini 
1507f26d0771SStefano Zampini   PetscFunctionBegin;
1508f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr);
1509f26d0771SStefano Zampini   ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr);
1510f26d0771SStefano Zampini   ierr = ISGetIndices(row,&rl);CHKERRQ(ierr);
1511f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr);
1512f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1513f0ae7da4SStefano Zampini   for (i=0;i<nrl;i++) if (rl[i]>=nrg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local row index %D -> %D greater than maximum possible %D",i,rl[i],nrg);
1514f26d0771SStefano Zampini #endif
1515f26d0771SStefano Zampini   ierr = PetscMalloc1(nrg,&idxs);CHKERRQ(ierr);
1516f26d0771SStefano Zampini   /* map from [0,nrl) to row */
1517f26d0771SStefano Zampini   for (i=0;i<nrl;i++) idxs[i] = rl[i];
1518f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1519f26d0771SStefano Zampini   for (i=nrl;i<nrg;i++) idxs[i] = nrg;
1520f26d0771SStefano Zampini #else
1521f26d0771SStefano Zampini   for (i=nrl;i<nrg;i++) idxs[i] = -1;
1522f26d0771SStefano Zampini #endif
1523f26d0771SStefano Zampini   ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr);
1524f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr);
1525f26d0771SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
1526f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
1527f26d0771SStefano Zampini   ierr = ISDestroy(&is);CHKERRQ(ierr);
1528f26d0771SStefano Zampini   /* compute new l2g map for columns */
1529f26d0771SStefano Zampini   if (col != row || A->rmap->mapping != A->cmap->mapping) {
1530f26d0771SStefano Zampini     const PetscInt *cg,*cl;
1531f26d0771SStefano Zampini     PetscInt       ncg;
1532f26d0771SStefano Zampini     PetscInt       ncl;
1533f26d0771SStefano Zampini 
1534f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr);
1535f26d0771SStefano Zampini     ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr);
1536f26d0771SStefano Zampini     ierr = ISGetIndices(col,&cl);CHKERRQ(ierr);
1537f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr);
1538f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1539f0ae7da4SStefano Zampini     for (i=0;i<ncl;i++) if (cl[i]>=ncg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local column index %D -> %D greater than maximum possible %D",i,cl[i],ncg);
1540f26d0771SStefano Zampini #endif
1541f26d0771SStefano Zampini     ierr = PetscMalloc1(ncg,&idxs);CHKERRQ(ierr);
1542f26d0771SStefano Zampini     /* map from [0,ncl) to col */
1543f26d0771SStefano Zampini     for (i=0;i<ncl;i++) idxs[i] = cl[i];
1544f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1545f26d0771SStefano Zampini     for (i=ncl;i<ncg;i++) idxs[i] = ncg;
1546f26d0771SStefano Zampini #else
1547f26d0771SStefano Zampini     for (i=ncl;i<ncg;i++) idxs[i] = -1;
1548f26d0771SStefano Zampini #endif
1549f26d0771SStefano Zampini     ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr);
1550f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr);
1551f26d0771SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
1552f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
1553f26d0771SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
1554f26d0771SStefano Zampini   } else {
1555f26d0771SStefano Zampini     ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr);
1556f26d0771SStefano Zampini     cl2g = rl2g;
1557f26d0771SStefano Zampini   }
1558f26d0771SStefano Zampini   /* create the MATIS submatrix */
1559f26d0771SStefano Zampini   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1560f26d0771SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr);
1561f26d0771SStefano Zampini   ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1562f26d0771SStefano Zampini   ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr);
1563b0aa3428SStefano Zampini   matis = (Mat_IS*)((*submat)->data);
1564f26d0771SStefano Zampini   matis->islocalref = PETSC_TRUE;
1565f26d0771SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr);
1566f26d0771SStefano Zampini   ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr);
1567f26d0771SStefano Zampini   ierr = MatISSetLocalMat(*submat,lA);CHKERRQ(ierr);
1568f26d0771SStefano Zampini   ierr = MatSetUp(*submat);CHKERRQ(ierr);
1569f26d0771SStefano Zampini   ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1570f26d0771SStefano Zampini   ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1571f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
1572f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
1573f26d0771SStefano Zampini   /* remove unsupported ops */
1574f26d0771SStefano Zampini   ierr = PetscMemzero((*submat)->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1575f26d0771SStefano Zampini   (*submat)->ops->destroy               = MatDestroy_IS;
1576f26d0771SStefano Zampini   (*submat)->ops->setvalueslocal        = MatSetValuesLocal_SubMat_IS;
1577f26d0771SStefano Zampini   (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS;
1578f26d0771SStefano Zampini   (*submat)->ops->assemblybegin         = MatAssemblyBegin_IS;
1579f26d0771SStefano Zampini   (*submat)->ops->assemblyend           = MatAssemblyEnd_IS;
1580f26d0771SStefano Zampini   PetscFunctionReturn(0);
1581f26d0771SStefano Zampini }
1582f26d0771SStefano Zampini 
1583f26d0771SStefano Zampini #undef __FUNCT__
1584284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
1585284134d9SBarry Smith /*@
15863c212e90SHong Zhang     MatCreateIS - Creates a "process" unassembled matrix, assembled on each
1587284134d9SBarry Smith        process but not across processes.
1588284134d9SBarry Smith 
1589284134d9SBarry Smith    Input Parameters:
1590284134d9SBarry Smith +     comm    - MPI communicator that will share the matrix
1591e176bc59SStefano Zampini .     bs      - block size of the matrix
1592df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
1593e176bc59SStefano Zampini .     rmap    - local to global map for rows
1594e176bc59SStefano Zampini -     cmap    - local to global map for cols
1595284134d9SBarry Smith 
1596284134d9SBarry Smith    Output Parameter:
1597284134d9SBarry Smith .    A - the resulting matrix
1598284134d9SBarry Smith 
15998e6c10adSSatish Balay    Level: advanced
16008e6c10adSSatish Balay 
16013c212e90SHong Zhang    Notes: See MATIS for more details.
16023c212e90SHong Zhang           m and n are NOT related to the size of the map; they are the size of the part of the vector owned
1603e176bc59SStefano Zampini           by that process. The sizes of rmap and cmap define the size of the local matrices.
16043c212e90SHong Zhang           If either rmap or cmap are NULL, then the matrix is assumed to be square.
1605284134d9SBarry Smith 
1606284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
1607284134d9SBarry Smith @*/
1608e176bc59SStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A)
1609284134d9SBarry Smith {
1610284134d9SBarry Smith   PetscErrorCode ierr;
1611284134d9SBarry Smith 
1612284134d9SBarry Smith   PetscFunctionBegin;
1613e176bc59SStefano Zampini   if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping");
1614284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1615d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
1616284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1617284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
1618e176bc59SStefano Zampini   if (rmap && cmap) {
1619e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr);
1620e176bc59SStefano Zampini   } else if (!rmap) {
1621e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr);
1622e176bc59SStefano Zampini   } else {
1623e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr);
1624e176bc59SStefano Zampini   }
1625284134d9SBarry Smith   PetscFunctionReturn(0);
1626284134d9SBarry Smith }
1627284134d9SBarry Smith 
1628b4319ba4SBarry Smith /*MC
1629f26d0771SStefano Zampini    MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition methods (e.g. PCBDDC or KSPFETIDP).
1630b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
1631b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
1632b4319ba4SBarry Smith    product is handled "implicitly".
1633b4319ba4SBarry Smith 
1634b4319ba4SBarry Smith    Operations Provided:
16356726f965SBarry Smith +  MatMult()
16362e74eeadSLisandro Dalcin .  MatMultAdd()
16372e74eeadSLisandro Dalcin .  MatMultTranspose()
16382e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
16396726f965SBarry Smith .  MatZeroEntries()
16406726f965SBarry Smith .  MatSetOption()
16412e74eeadSLisandro Dalcin .  MatZeroRows()
16422e74eeadSLisandro Dalcin .  MatSetValues()
164397563a80SStefano Zampini .  MatSetValuesBlocked()
16446726f965SBarry Smith .  MatSetValuesLocal()
164597563a80SStefano Zampini .  MatSetValuesBlockedLocal()
16462e74eeadSLisandro Dalcin .  MatScale()
16472e74eeadSLisandro Dalcin .  MatGetDiagonal()
16482b404112SStefano Zampini .  MatMissingDiagonal()
16492b404112SStefano Zampini .  MatDuplicate()
16502b404112SStefano Zampini .  MatCopy()
1651f26d0771SStefano Zampini .  MatAXPY()
1652f26d0771SStefano Zampini .  MatGetSubMatrix()
1653f26d0771SStefano Zampini .  MatGetLocalSubMatrix()
16546726f965SBarry Smith -  MatSetLocalToGlobalMapping()
1655b4319ba4SBarry Smith 
1656b4319ba4SBarry Smith    Options Database Keys:
1657b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
1658b4319ba4SBarry Smith 
1659b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
1660b4319ba4SBarry Smith 
1661b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
1662b4319ba4SBarry Smith 
1663b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
1664eb82efa4SStefano Zampini           MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation()
1665b4319ba4SBarry Smith 
1666b4319ba4SBarry Smith   Level: advanced
1667b4319ba4SBarry Smith 
1668f26d0771SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), MatCreateIS(), PCBDDC, KSPFETIDP
1669b4319ba4SBarry Smith 
1670b4319ba4SBarry Smith M*/
1671b4319ba4SBarry Smith 
1672b4319ba4SBarry Smith #undef __FUNCT__
1673b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
16748cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
1675b4319ba4SBarry Smith {
1676dfbe8321SBarry Smith   PetscErrorCode ierr;
1677b4319ba4SBarry Smith   Mat_IS         *b;
1678b4319ba4SBarry Smith 
1679b4319ba4SBarry Smith   PetscFunctionBegin;
1680b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1681b4319ba4SBarry Smith   A->data = (void*)b;
1682b4319ba4SBarry Smith 
1683e176bc59SStefano Zampini   /* matrix ops */
1684e176bc59SStefano Zampini   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1685b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
16862e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
16872e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
16882e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1689b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1690b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
16912e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
169298921651SStefano Zampini   A->ops->setvaluesblocked        = MatSetValuesBlocked_IS;
1693b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1694f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
16952e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1696f0ae7da4SStefano Zampini   A->ops->zerorowscolumns         = MatZeroRowsColumns_IS;
1697b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1698b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1699b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
17006726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
17012e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
17022e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
17036726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
170469796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
170569796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1706*45471136SStefano Zampini   A->ops->isstructurallysymmetric = MatIsStructurallySymmetric_IS;
1707ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
17086bd84002SStefano Zampini   A->ops->missingdiagonal         = MatMissingDiagonal_IS;
17092b404112SStefano Zampini   A->ops->copy                    = MatCopy_IS;
1710659959c5SStefano Zampini   A->ops->getlocalsubmatrix       = MatGetLocalSubMatrix_IS;
1711a8116848SStefano Zampini   A->ops->getsubmatrix            = MatGetSubMatrix_IS;
1712f26d0771SStefano Zampini   A->ops->axpy                    = MatAXPY_IS;
17133fd1c9e7SStefano Zampini   A->ops->diagonalset             = MatDiagonalSet_IS;
17143fd1c9e7SStefano Zampini   A->ops->shift                   = MatShift_IS;
1715b4319ba4SBarry Smith 
1716b7ce53b6SStefano Zampini   /* special MATIS functions */
1717bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1718bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1719b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
17202e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
1721cf0a3239SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",MatISSetUpSF_IS);CHKERRQ(ierr);
172217667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1723b4319ba4SBarry Smith   PetscFunctionReturn(0);
1724b4319ba4SBarry Smith }
1725