xref: /petsc/src/mat/impls/is/matis.c (revision f0ae7da42e5dc6744970261bd164770f012c0037)
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*/
1128f4e0baSStefano Zampini 
12f26d0771SStefano Zampini #define MATIS_MAX_ENTRIES_INSERTION 2048
13f26d0771SStefano Zampini 
14a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat);
15a72627d2SStefano Zampini 
1628f4e0baSStefano Zampini #undef __FUNCT__
173fd1c9e7SStefano Zampini #define __FUNCT__ "MatDiagonalSet_IS"
183fd1c9e7SStefano Zampini PetscErrorCode  MatDiagonalSet_IS(Mat A,Vec D,InsertMode insmode)
193fd1c9e7SStefano Zampini {
203fd1c9e7SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
213fd1c9e7SStefano Zampini   PetscErrorCode ierr;
223fd1c9e7SStefano Zampini 
233fd1c9e7SStefano Zampini   PetscFunctionBegin;
243fd1c9e7SStefano Zampini   if (!D) { /* this code branch is used by MatShift_IS */
253fd1c9e7SStefano Zampini     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
263fd1c9e7SStefano Zampini   } else {
273fd1c9e7SStefano Zampini     ierr = VecScatterBegin(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
283fd1c9e7SStefano Zampini     ierr = VecScatterEnd(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
293fd1c9e7SStefano Zampini   }
303fd1c9e7SStefano Zampini   ierr = VecPointwiseDivide(is->y,is->y,is->counter);CHKERRQ(ierr);
313fd1c9e7SStefano Zampini   ierr = MatDiagonalSet(is->A,is->y,insmode);CHKERRQ(ierr);
323fd1c9e7SStefano Zampini   PetscFunctionReturn(0);
333fd1c9e7SStefano Zampini }
343fd1c9e7SStefano Zampini 
353fd1c9e7SStefano Zampini #undef __FUNCT__
363fd1c9e7SStefano Zampini #define __FUNCT__ "MatShift_IS"
373fd1c9e7SStefano Zampini PetscErrorCode  MatShift_IS(Mat A,PetscScalar a)
383fd1c9e7SStefano Zampini {
393fd1c9e7SStefano Zampini   PetscErrorCode ierr;
403fd1c9e7SStefano Zampini 
413fd1c9e7SStefano Zampini   PetscFunctionBegin;
423fd1c9e7SStefano Zampini   ierr = MatDiagonalSet_IS(A,NULL,ADD_VALUES);CHKERRQ(ierr);
433fd1c9e7SStefano Zampini   PetscFunctionReturn(0);
443fd1c9e7SStefano Zampini }
453fd1c9e7SStefano Zampini 
463fd1c9e7SStefano Zampini #undef __FUNCT__
47f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesLocal_SubMat_IS"
48f26d0771SStefano Zampini static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
49f26d0771SStefano Zampini {
50f26d0771SStefano Zampini   PetscErrorCode ierr;
51f26d0771SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
52f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
53f26d0771SStefano Zampini 
54f26d0771SStefano Zampini   PetscFunctionBegin;
55f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
56f26d0771SStefano 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);
57f26d0771SStefano Zampini #endif
58f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr);
59f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr);
60f26d0771SStefano Zampini   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
61f26d0771SStefano Zampini   PetscFunctionReturn(0);
62f26d0771SStefano Zampini }
63f26d0771SStefano Zampini 
64f26d0771SStefano Zampini #undef __FUNCT__
65f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesBlockedLocal_SubMat_IS"
66f26d0771SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
67f26d0771SStefano Zampini {
68f26d0771SStefano Zampini   PetscErrorCode ierr;
69f26d0771SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
70f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
71f26d0771SStefano Zampini 
72f26d0771SStefano Zampini   PetscFunctionBegin;
73f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
74f26d0771SStefano 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);
75f26d0771SStefano Zampini #endif
76f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr);
77f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr);
78f26d0771SStefano Zampini   ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
79f26d0771SStefano Zampini   PetscFunctionReturn(0);
80f26d0771SStefano Zampini }
81f26d0771SStefano Zampini 
82f26d0771SStefano Zampini #undef __FUNCT__
83a8116848SStefano Zampini #define __FUNCT__ "PetscLayoutMapLocal_Private"
84*f0ae7da4SStefano Zampini static PetscErrorCode PetscLayoutMapLocal_Private(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs)
85a8116848SStefano Zampini {
86a8116848SStefano Zampini   PetscInt      *owners = map->range;
87a8116848SStefano Zampini   PetscInt       n      = map->n;
88a8116848SStefano Zampini   PetscSF        sf;
89a8116848SStefano Zampini   PetscInt      *lidxs,*work = NULL;
90a8116848SStefano Zampini   PetscSFNode   *ridxs;
91a8116848SStefano Zampini   PetscMPIInt    rank;
92a8116848SStefano Zampini   PetscInt       r, p = 0, len = 0;
93a8116848SStefano Zampini   PetscErrorCode ierr;
94a8116848SStefano Zampini 
95a8116848SStefano Zampini   PetscFunctionBegin;
96a8116848SStefano Zampini   /* Create SF where leaves are input idxs and roots are owned idxs (code adapted from MatZeroRowsMapLocal_Private) */
97a8116848SStefano Zampini   ierr = MPI_Comm_rank(map->comm,&rank);CHKERRQ(ierr);
98a8116848SStefano Zampini   ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr);
99a8116848SStefano Zampini   for (r = 0; r < n; ++r) lidxs[r] = -1;
100a8116848SStefano Zampini   ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr);
101a8116848SStefano Zampini   for (r = 0; r < N; ++r) {
102a8116848SStefano Zampini     const PetscInt idx = idxs[r];
103a8116848SStefano 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);
104a8116848SStefano Zampini     if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */
105a8116848SStefano Zampini       ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr);
106a8116848SStefano Zampini     }
107a8116848SStefano Zampini     ridxs[r].rank = p;
108a8116848SStefano Zampini     ridxs[r].index = idxs[r] - owners[p];
109a8116848SStefano Zampini   }
110a8116848SStefano Zampini   ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr);
111a8116848SStefano Zampini   ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr);
112a8116848SStefano Zampini   ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr);
113a8116848SStefano Zampini   ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr);
114*f0ae7da4SStefano Zampini   if (ogidxs) { /* communicate global idxs */
115a8116848SStefano Zampini     PetscInt cum = 0,start,*work2;
116*f0ae7da4SStefano Zampini 
117*f0ae7da4SStefano Zampini     ierr = PetscMalloc1(n,&work);CHKERRQ(ierr);
118a8116848SStefano Zampini     ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr);
119a8116848SStefano Zampini     for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++;
120a8116848SStefano Zampini     ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRQ(ierr);
121a8116848SStefano Zampini     start -= cum;
122a8116848SStefano Zampini     cum = 0;
123a8116848SStefano Zampini     for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++;
124a8116848SStefano Zampini     ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr);
125a8116848SStefano Zampini     ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr);
126a8116848SStefano Zampini     ierr = PetscFree(work2);CHKERRQ(ierr);
127a8116848SStefano Zampini   }
128a8116848SStefano Zampini   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
129a8116848SStefano Zampini   /* Compress and put in indices */
130a8116848SStefano Zampini   for (r = 0; r < n; ++r)
131a8116848SStefano Zampini     if (lidxs[r] >= 0) {
132a8116848SStefano Zampini       if (work) work[len] = work[r];
133a8116848SStefano Zampini       lidxs[len++] = r;
134a8116848SStefano Zampini     }
135a8116848SStefano Zampini   if (on) *on = len;
136a8116848SStefano Zampini   if (oidxs) *oidxs = lidxs;
137a8116848SStefano Zampini   if (ogidxs) *ogidxs = work;
138a8116848SStefano Zampini   PetscFunctionReturn(0);
139a8116848SStefano Zampini }
140a8116848SStefano Zampini 
141a8116848SStefano Zampini #undef __FUNCT__
142a8116848SStefano Zampini #define __FUNCT__ "MatGetSubMatrix_IS"
143a8116848SStefano Zampini static PetscErrorCode MatGetSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat)
144a8116848SStefano Zampini {
145a8116848SStefano Zampini   Mat               locmat,newlocmat;
146a8116848SStefano Zampini   Mat_IS            *newmatis;
147a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG)
148a8116848SStefano Zampini   Vec               rtest,ltest;
149a8116848SStefano Zampini   const PetscScalar *array;
150a8116848SStefano Zampini #endif
151a8116848SStefano Zampini   const PetscInt    *idxs;
152a8116848SStefano Zampini   PetscInt          i,m,n;
153a8116848SStefano Zampini   PetscErrorCode    ierr;
154a8116848SStefano Zampini 
155a8116848SStefano Zampini   PetscFunctionBegin;
156a8116848SStefano Zampini   if (scall == MAT_REUSE_MATRIX) {
157a8116848SStefano Zampini     PetscBool ismatis;
158a8116848SStefano Zampini 
159a8116848SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr);
160a8116848SStefano Zampini     if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type");
161a8116848SStefano Zampini     newmatis = (Mat_IS*)(*newmat)->data;
162a8116848SStefano Zampini     if (!newmatis->getsub_ris) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS");
163a8116848SStefano Zampini     if (!newmatis->getsub_cis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS");
164a8116848SStefano Zampini   }
165a8116848SStefano Zampini   /* irow and icol may not have duplicate entries */
166a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG)
167a8116848SStefano Zampini   ierr = MatCreateVecs(mat,&ltest,&rtest);CHKERRQ(ierr);
168a8116848SStefano Zampini   ierr = ISGetLocalSize(irow,&n);CHKERRQ(ierr);
169a8116848SStefano Zampini   ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr);
170a8116848SStefano Zampini   for (i=0;i<n;i++) {
171a8116848SStefano Zampini     ierr = VecSetValue(rtest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr);
172a8116848SStefano Zampini   }
173a8116848SStefano Zampini   ierr = VecAssemblyBegin(rtest);CHKERRQ(ierr);
174a8116848SStefano Zampini   ierr = VecAssemblyEnd(rtest);CHKERRQ(ierr);
175a8116848SStefano Zampini   ierr = VecGetLocalSize(rtest,&n);CHKERRQ(ierr);
176a8116848SStefano Zampini   ierr = VecGetOwnershipRange(rtest,&m,NULL);CHKERRQ(ierr);
177a8116848SStefano Zampini   ierr = VecGetArrayRead(rtest,&array);CHKERRQ(ierr);
178fd479f66SMatthew 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]));
179a8116848SStefano Zampini   ierr = VecRestoreArrayRead(rtest,&array);CHKERRQ(ierr);
180a8116848SStefano Zampini   ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr);
181a8116848SStefano Zampini   ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr);
182a8116848SStefano Zampini   ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr);
183a8116848SStefano Zampini   for (i=0;i<n;i++) {
184a8116848SStefano Zampini     ierr = VecSetValue(ltest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr);
185a8116848SStefano Zampini   }
186a8116848SStefano Zampini   ierr = VecAssemblyBegin(ltest);CHKERRQ(ierr);
187a8116848SStefano Zampini   ierr = VecAssemblyEnd(ltest);CHKERRQ(ierr);
188a8116848SStefano Zampini   ierr = VecGetLocalSize(ltest,&n);CHKERRQ(ierr);
189a8116848SStefano Zampini   ierr = VecGetOwnershipRange(ltest,&m,NULL);CHKERRQ(ierr);
190a8116848SStefano Zampini   ierr = VecGetArrayRead(ltest,&array);CHKERRQ(ierr);
191fd479f66SMatthew 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]));
192a8116848SStefano Zampini   ierr = VecRestoreArrayRead(ltest,&array);CHKERRQ(ierr);
193a8116848SStefano Zampini   ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr);
194a8116848SStefano Zampini   ierr = VecDestroy(&rtest);CHKERRQ(ierr);
195a8116848SStefano Zampini   ierr = VecDestroy(&ltest);CHKERRQ(ierr);
196a8116848SStefano Zampini #endif
197a8116848SStefano Zampini   if (scall == MAT_INITIAL_MATRIX) {
198a8116848SStefano Zampini     Mat_IS                 *matis = (Mat_IS*)mat->data;
199a8116848SStefano Zampini     ISLocalToGlobalMapping rl2g;
200a8116848SStefano Zampini     IS                     is;
201a8116848SStefano Zampini     PetscInt               *lidxs,*lgidxs,*newgidxs;
2026dd40735SStefano Zampini     PetscInt               ll,newloc;
203a8116848SStefano Zampini     MPI_Comm               comm;
204a8116848SStefano Zampini 
205a8116848SStefano Zampini     ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
206a8116848SStefano Zampini     ierr = ISGetLocalSize(irow,&m);CHKERRQ(ierr);
207a8116848SStefano Zampini     ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr);
208a8116848SStefano Zampini     ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
209a8116848SStefano Zampini     ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr);
210a8116848SStefano Zampini     ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
211a8116848SStefano Zampini     /* communicate irow to their owners in the layout */
212a8116848SStefano Zampini     ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr);
213*f0ae7da4SStefano Zampini     ierr = PetscLayoutMapLocal_Private(mat->rmap,m,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr);
214a8116848SStefano Zampini     ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr);
215a8116848SStefano Zampini     if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
216a8116848SStefano Zampini       ierr = MatISComputeSF_Private(mat);CHKERRQ(ierr);
217a8116848SStefano Zampini     }
218a8116848SStefano Zampini     ierr = PetscMemzero(matis->sf_rootdata,matis->sf_nroots*sizeof(PetscInt));CHKERRQ(ierr);
219a8116848SStefano Zampini     for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1;
220a8116848SStefano Zampini     ierr = PetscFree(lidxs);CHKERRQ(ierr);
221a8116848SStefano Zampini     ierr = PetscFree(lgidxs);CHKERRQ(ierr);
222a8116848SStefano Zampini     ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
223a8116848SStefano Zampini     ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
224a8116848SStefano Zampini     for (i=0,newloc=0;i<matis->sf_nleaves;i++) if (matis->sf_leafdata[i]) newloc++;
225a8116848SStefano Zampini     ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr);
226a8116848SStefano Zampini     ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr);
227a8116848SStefano Zampini     for (i=0,newloc=0;i<matis->sf_nleaves;i++)
228a8116848SStefano Zampini       if (matis->sf_leafdata[i]) {
229a8116848SStefano Zampini         lidxs[newloc] = i;
230a8116848SStefano Zampini         newgidxs[newloc++] = matis->sf_leafdata[i]-1;
231a8116848SStefano Zampini       }
232a8116848SStefano Zampini     ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
233a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
234a8116848SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
235a8116848SStefano Zampini     /* local is to extract local submatrix */
236a8116848SStefano Zampini     newmatis = (Mat_IS*)(*newmat)->data;
237a8116848SStefano Zampini     ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris);CHKERRQ(ierr);
238a8116848SStefano Zampini     if (mat->congruentlayouts == -1) { /* first time we compare rows and cols layouts */
239a8116848SStefano Zampini       PetscBool cong;
240a8116848SStefano Zampini       ierr = PetscLayoutCompare(mat->rmap,mat->cmap,&cong);CHKERRQ(ierr);
241a8116848SStefano Zampini       if (cong) mat->congruentlayouts = 1;
242a8116848SStefano Zampini       else      mat->congruentlayouts = 0;
243a8116848SStefano Zampini     }
244a8116848SStefano Zampini     if (mat->congruentlayouts && irow == icol) {
245a8116848SStefano Zampini       ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g);CHKERRQ(ierr);
246a8116848SStefano Zampini       ierr = PetscObjectReference((PetscObject)newmatis->getsub_ris);CHKERRQ(ierr);
247a8116848SStefano Zampini       newmatis->getsub_cis = newmatis->getsub_ris;
248a8116848SStefano Zampini     } else {
249a8116848SStefano Zampini       ISLocalToGlobalMapping cl2g;
250a8116848SStefano Zampini 
251a8116848SStefano Zampini       /* communicate icol to their owners in the layout */
252a8116848SStefano Zampini       ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr);
253*f0ae7da4SStefano Zampini       ierr = PetscLayoutMapLocal_Private(mat->cmap,n,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr);
254a8116848SStefano Zampini       ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr);
255a8116848SStefano Zampini       ierr = PetscMemzero(matis->csf_rootdata,matis->csf_nroots*sizeof(PetscInt));CHKERRQ(ierr);
256a8116848SStefano Zampini       for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1;
257a8116848SStefano Zampini       ierr = PetscFree(lidxs);CHKERRQ(ierr);
258a8116848SStefano Zampini       ierr = PetscFree(lgidxs);CHKERRQ(ierr);
259a8116848SStefano Zampini       ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr);
260a8116848SStefano Zampini       ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr);
261a8116848SStefano Zampini       for (i=0,newloc=0;i<matis->csf_nleaves;i++) if (matis->csf_leafdata[i]) newloc++;
262a8116848SStefano Zampini       ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr);
263a8116848SStefano Zampini       ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr);
264a8116848SStefano Zampini       for (i=0,newloc=0;i<matis->csf_nleaves;i++)
265a8116848SStefano Zampini         if (matis->csf_leafdata[i]) {
266a8116848SStefano Zampini           lidxs[newloc] = i;
267a8116848SStefano Zampini           newgidxs[newloc++] = matis->csf_leafdata[i]-1;
268a8116848SStefano Zampini         }
269a8116848SStefano Zampini       ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
270a8116848SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
271a8116848SStefano Zampini       ierr = ISDestroy(&is);CHKERRQ(ierr);
272a8116848SStefano Zampini       /* local is to extract local submatrix */
273a8116848SStefano Zampini       ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis);CHKERRQ(ierr);
274a8116848SStefano Zampini       ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr);
275a8116848SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
276a8116848SStefano Zampini     }
277a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
278a8116848SStefano Zampini   } else {
279a8116848SStefano Zampini     ierr = MatISGetLocalMat(*newmat,&newlocmat);CHKERRQ(ierr);
280a8116848SStefano Zampini   }
281a8116848SStefano Zampini   ierr = MatISGetLocalMat(mat,&locmat);CHKERRQ(ierr);
282a8116848SStefano Zampini   newmatis = (Mat_IS*)(*newmat)->data;
283a8116848SStefano Zampini   ierr = MatGetSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat);CHKERRQ(ierr);
284a8116848SStefano Zampini   if (scall == MAT_INITIAL_MATRIX) {
285a8116848SStefano Zampini     ierr = MatISSetLocalMat(*newmat,newlocmat);CHKERRQ(ierr);
286a8116848SStefano Zampini     ierr = MatDestroy(&newlocmat);CHKERRQ(ierr);
287a8116848SStefano Zampini   }
288a8116848SStefano Zampini   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
289a8116848SStefano Zampini   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
290a8116848SStefano Zampini   PetscFunctionReturn(0);
291a8116848SStefano Zampini }
292a8116848SStefano Zampini 
293a8116848SStefano Zampini #undef __FUNCT__
2942b404112SStefano Zampini #define __FUNCT__ "MatCopy_IS"
295a8116848SStefano Zampini static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str)
2962b404112SStefano Zampini {
2972b404112SStefano Zampini   Mat_IS         *a = (Mat_IS*)A->data,*b;
2982b404112SStefano Zampini   PetscBool      ismatis;
2992b404112SStefano Zampini   PetscErrorCode ierr;
3002b404112SStefano Zampini 
3012b404112SStefano Zampini   PetscFunctionBegin;
3022b404112SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr);
3032b404112SStefano Zampini   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented");
3042b404112SStefano Zampini   b = (Mat_IS*)B->data;
3052b404112SStefano Zampini   ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr);
3062b404112SStefano Zampini   PetscFunctionReturn(0);
3072b404112SStefano Zampini }
3082b404112SStefano Zampini 
3092b404112SStefano Zampini #undef __FUNCT__
3106bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS"
311a8116848SStefano Zampini static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool  *missing,PetscInt *d)
3126bd84002SStefano Zampini {
313527b2640SStefano Zampini   Vec               v;
314527b2640SStefano Zampini   const PetscScalar *array;
315527b2640SStefano Zampini   PetscInt          i,n;
3166bd84002SStefano Zampini   PetscErrorCode    ierr;
3176bd84002SStefano Zampini 
3186bd84002SStefano Zampini   PetscFunctionBegin;
319527b2640SStefano Zampini   *missing = PETSC_FALSE;
320527b2640SStefano Zampini   ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr);
321527b2640SStefano Zampini   ierr = MatGetDiagonal(A,v);CHKERRQ(ierr);
322527b2640SStefano Zampini   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
323527b2640SStefano Zampini   ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr);
324527b2640SStefano Zampini   for (i=0;i<n;i++) if (array[i] == 0.) break;
325527b2640SStefano Zampini   ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr);
326527b2640SStefano Zampini   ierr = VecDestroy(&v);CHKERRQ(ierr);
327527b2640SStefano Zampini   if (i != n) *missing = PETSC_TRUE;
328527b2640SStefano Zampini   if (d) {
329527b2640SStefano Zampini     *d = -1;
330527b2640SStefano Zampini     if (*missing) {
331527b2640SStefano Zampini       PetscInt rstart;
332527b2640SStefano Zampini       ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
333527b2640SStefano Zampini       *d = i+rstart;
334527b2640SStefano Zampini     }
335527b2640SStefano Zampini   }
3366bd84002SStefano Zampini   PetscFunctionReturn(0);
3376bd84002SStefano Zampini }
3386bd84002SStefano Zampini 
3396bd84002SStefano Zampini #undef __FUNCT__
34028f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private"
341a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B)
34228f4e0baSStefano Zampini {
34328f4e0baSStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
34428f4e0baSStefano Zampini   const PetscInt *gidxs;
34528f4e0baSStefano Zampini   PetscErrorCode ierr;
34628f4e0baSStefano Zampini 
34728f4e0baSStefano Zampini   PetscFunctionBegin;
34828f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr);
34928f4e0baSStefano Zampini   ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr);
35028f4e0baSStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr);
3513bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
35228f4e0baSStefano Zampini   ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
3533bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
35428f4e0baSStefano Zampini   ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr);
355a8116848SStefano Zampini   if (B->rmap->mapping != B->cmap->mapping) { /* setup SF for columns */
356a8116848SStefano Zampini     ierr = MatGetSize(matis->A,NULL,&matis->csf_nleaves);CHKERRQ(ierr);
357a8116848SStefano Zampini     ierr = MatGetLocalSize(B,NULL,&matis->csf_nroots);CHKERRQ(ierr);
358a8116848SStefano Zampini     ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf);CHKERRQ(ierr);
359a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr);
360a8116848SStefano Zampini     ierr = PetscSFSetGraphLayout(matis->csf,B->cmap,matis->csf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
361a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr);
362a8116848SStefano Zampini     ierr = PetscMalloc2(matis->csf_nroots,&matis->csf_rootdata,matis->csf_nleaves,&matis->csf_leafdata);CHKERRQ(ierr);
363a8116848SStefano Zampini   } else {
364a8116848SStefano Zampini     matis->csf = matis->sf;
365a8116848SStefano Zampini     matis->csf_nleaves = matis->sf_nleaves;
366a8116848SStefano Zampini     matis->csf_nroots = matis->sf_nroots;
367a8116848SStefano Zampini     matis->csf_leafdata = matis->sf_leafdata;
368a8116848SStefano Zampini     matis->csf_rootdata = matis->sf_rootdata;
369a8116848SStefano Zampini   }
37028f4e0baSStefano Zampini   PetscFunctionReturn(0);
37128f4e0baSStefano Zampini }
3722e1947a5SStefano Zampini 
3732e1947a5SStefano Zampini #undef __FUNCT__
3742e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation"
375eb82efa4SStefano Zampini /*@
376a88811baSStefano Zampini    MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix.
377a88811baSStefano Zampini 
378a88811baSStefano Zampini    Collective on MPI_Comm
379a88811baSStefano Zampini 
380a88811baSStefano Zampini    Input Parameters:
381a88811baSStefano Zampini +  B - the matrix
382a88811baSStefano Zampini .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
383a88811baSStefano Zampini            (same value is used for all local rows)
384a88811baSStefano Zampini .  d_nnz - array containing the number of nonzeros in the various rows of the
385a88811baSStefano Zampini            DIAGONAL portion of the local submatrix (possibly different for each row)
386a88811baSStefano Zampini            or NULL, if d_nz is used to specify the nonzero structure.
387a88811baSStefano Zampini            The size of this array is equal to the number of local rows, i.e 'm'.
388a88811baSStefano Zampini            For matrices that will be factored, you must leave room for (and set)
389a88811baSStefano Zampini            the diagonal entry even if it is zero.
390a88811baSStefano Zampini .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
391a88811baSStefano Zampini            submatrix (same value is used for all local rows).
392a88811baSStefano Zampini -  o_nnz - array containing the number of nonzeros in the various rows of the
393a88811baSStefano Zampini            OFF-DIAGONAL portion of the local submatrix (possibly different for
394a88811baSStefano Zampini            each row) or NULL, if o_nz is used to specify the nonzero
395a88811baSStefano Zampini            structure. The size of this array is equal to the number
396a88811baSStefano Zampini            of local rows, i.e 'm'.
397a88811baSStefano Zampini 
398a88811baSStefano Zampini    If the *_nnz parameter is given then the *_nz parameter is ignored
399a88811baSStefano Zampini 
400a88811baSStefano Zampini    Level: intermediate
401a88811baSStefano Zampini 
402a88811baSStefano Zampini    Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition
403a88811baSStefano Zampini           from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local
404a88811baSStefano Zampini           matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects.
405a88811baSStefano Zampini 
406a88811baSStefano Zampini .keywords: matrix
407a88811baSStefano Zampini 
4083c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS
409a88811baSStefano Zampini @*/
4102e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
4112e1947a5SStefano Zampini {
4122e1947a5SStefano Zampini   PetscErrorCode ierr;
4132e1947a5SStefano Zampini 
4142e1947a5SStefano Zampini   PetscFunctionBegin;
4152e1947a5SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
4162e1947a5SStefano Zampini   PetscValidType(B,1);
4172e1947a5SStefano Zampini   ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
4182e1947a5SStefano Zampini   PetscFunctionReturn(0);
4192e1947a5SStefano Zampini }
4202e1947a5SStefano Zampini 
4212e1947a5SStefano Zampini #undef __FUNCT__
4222e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS"
4237230de76SStefano Zampini static PetscErrorCode  MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
4242e1947a5SStefano Zampini {
4252e1947a5SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
42628f4e0baSStefano Zampini   PetscInt       bs,i,nlocalcols;
4272e1947a5SStefano Zampini   PetscErrorCode ierr;
4282e1947a5SStefano Zampini 
4292e1947a5SStefano Zampini   PetscFunctionBegin;
4306c4ed002SBarry Smith   if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
43128f4e0baSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
43228f4e0baSStefano Zampini     ierr = MatISComputeSF_Private(B);CHKERRQ(ierr);
43328f4e0baSStefano Zampini   }
4342e1947a5SStefano Zampini   if (!d_nnz) {
43528f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz;
4362e1947a5SStefano Zampini   } else {
43728f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i];
4382e1947a5SStefano Zampini   }
4392e1947a5SStefano Zampini   if (!o_nnz) {
44028f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz;
4412e1947a5SStefano Zampini   } else {
44228f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i];
4432e1947a5SStefano Zampini   }
44428f4e0baSStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
44528f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr);
44628f4e0baSStefano Zampini   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
44728f4e0baSStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
44828f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves;i++) {
44928f4e0baSStefano Zampini     matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols);
4502e1947a5SStefano Zampini   }
45128f4e0baSStefano Zampini   ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr);
45228f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
45328f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs;
4542e1947a5SStefano Zampini   }
45528f4e0baSStefano Zampini   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
45628f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
45728f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i]-i;
4582e1947a5SStefano Zampini   }
45928f4e0baSStefano Zampini   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
4602e1947a5SStefano Zampini   PetscFunctionReturn(0);
4612e1947a5SStefano Zampini }
462b4319ba4SBarry Smith 
463b4319ba4SBarry Smith #undef __FUNCT__
4643927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private"
4653927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
4663927de2eSStefano Zampini {
4673927de2eSStefano Zampini   Mat_IS          *matis = (Mat_IS*)(A->data);
4683927de2eSStefano Zampini   PetscInt        *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership;
469ecf5a873SStefano Zampini   const PetscInt  *global_indices_r,*global_indices_c;
4703927de2eSStefano Zampini   PetscInt        i,j,bs,rows,cols;
4713927de2eSStefano Zampini   PetscInt        lrows,lcols;
4723927de2eSStefano Zampini   PetscInt        local_rows,local_cols;
4733927de2eSStefano Zampini   PetscMPIInt     nsubdomains;
4743927de2eSStefano Zampini   PetscBool       isdense,issbaij;
4753927de2eSStefano Zampini   PetscErrorCode  ierr;
4763927de2eSStefano Zampini 
4773927de2eSStefano Zampini   PetscFunctionBegin;
4783927de2eSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr);
4793927de2eSStefano Zampini   ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr);
4803927de2eSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
4813927de2eSStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
4823927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
4833927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
484ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
485ecf5a873SStefano Zampini   if (A->rmap->mapping != A->cmap->mapping) {
4867230de76SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
487ecf5a873SStefano Zampini   } else {
488ecf5a873SStefano Zampini     global_indices_c = global_indices_r;
489ecf5a873SStefano Zampini   }
490ecf5a873SStefano Zampini 
4913927de2eSStefano Zampini   if (issbaij) {
4923927de2eSStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
4933927de2eSStefano Zampini   }
4943927de2eSStefano Zampini   /*
495ecf5a873SStefano Zampini      An SF reduce is needed to sum up properly on shared rows.
4963927de2eSStefano Zampini      Note that generally preallocation is not exact, since it overestimates nonzeros
4973927de2eSStefano Zampini   */
4983927de2eSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
4993927de2eSStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
5003927de2eSStefano Zampini   }
5013927de2eSStefano Zampini   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
5023927de2eSStefano Zampini   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
5033927de2eSStefano Zampini   /* All processes need to compute entire row ownership */
5043927de2eSStefano Zampini   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
5053927de2eSStefano Zampini   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
5063927de2eSStefano Zampini   for (i=0;i<nsubdomains;i++) {
5073927de2eSStefano Zampini     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
5083927de2eSStefano Zampini       row_ownership[j] = i;
5093927de2eSStefano Zampini     }
5103927de2eSStefano Zampini   }
5117230de76SStefano Zampini   ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
5123927de2eSStefano Zampini 
5133927de2eSStefano Zampini   /*
5143927de2eSStefano Zampini      my_dnz and my_onz contains exact contribution to preallocation from each local mat
5153927de2eSStefano Zampini      then, they will be summed up properly. This way, preallocation is always sufficient
5163927de2eSStefano Zampini   */
5173927de2eSStefano Zampini   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
5183927de2eSStefano Zampini   /* preallocation as a MATAIJ */
5193927de2eSStefano Zampini   if (isdense) { /* special case for dense local matrices */
5203927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
521ecf5a873SStefano Zampini       PetscInt index_row = global_indices_r[i];
5223927de2eSStefano Zampini       for (j=i;j<local_rows;j++) {
5233927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
524ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[j];
5253927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
5263927de2eSStefano Zampini           my_dnz[i] += 1;
5273927de2eSStefano Zampini         } else { /* offdiag block */
5283927de2eSStefano Zampini           my_onz[i] += 1;
5293927de2eSStefano Zampini         }
5303927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
5313927de2eSStefano Zampini         if (i != j) {
5323927de2eSStefano Zampini           owner = row_ownership[index_col];
5333927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
5343927de2eSStefano Zampini             my_dnz[j] += 1;
5353927de2eSStefano Zampini           } else {
5363927de2eSStefano Zampini             my_onz[j] += 1;
5373927de2eSStefano Zampini           }
5383927de2eSStefano Zampini         }
5393927de2eSStefano Zampini       }
5403927de2eSStefano Zampini     }
541bb1015c3SStefano Zampini   } else if (matis->A->ops->getrowij) {
542bb1015c3SStefano Zampini     const PetscInt *ii,*jj,*jptr;
543bb1015c3SStefano Zampini     PetscBool      done;
544bb1015c3SStefano Zampini     ierr = MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr);
545bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__);
546bb1015c3SStefano Zampini     jptr = jj;
547bb1015c3SStefano Zampini     for (i=0;i<local_rows;i++) {
548bb1015c3SStefano Zampini       PetscInt index_row = global_indices_r[i];
549bb1015c3SStefano Zampini       for (j=0;j<ii[i+1]-ii[i];j++,jptr++) {
550bb1015c3SStefano Zampini         PetscInt owner = row_ownership[index_row];
551bb1015c3SStefano Zampini         PetscInt index_col = global_indices_c[*jptr];
552bb1015c3SStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
553bb1015c3SStefano Zampini           my_dnz[i] += 1;
554bb1015c3SStefano Zampini         } else { /* offdiag block */
555bb1015c3SStefano Zampini           my_onz[i] += 1;
556bb1015c3SStefano Zampini         }
557bb1015c3SStefano Zampini         /* same as before, interchanging rows and cols */
558bb1015c3SStefano Zampini         if (issbaij && index_col != index_row) {
559bb1015c3SStefano Zampini           owner = row_ownership[index_col];
560bb1015c3SStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
561bb1015c3SStefano Zampini             my_dnz[*jptr] += 1;
562bb1015c3SStefano Zampini           } else {
563bb1015c3SStefano Zampini             my_onz[*jptr] += 1;
564bb1015c3SStefano Zampini           }
565bb1015c3SStefano Zampini         }
566bb1015c3SStefano Zampini       }
567bb1015c3SStefano Zampini     }
568bb1015c3SStefano Zampini     ierr = MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr);
569bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__);
570bb1015c3SStefano Zampini   } else { /* loop over rows and use MatGetRow */
5713927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
5723927de2eSStefano Zampini       const PetscInt *cols;
573ecf5a873SStefano Zampini       PetscInt       ncols,index_row = global_indices_r[i];
5743927de2eSStefano Zampini       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
5753927de2eSStefano Zampini       for (j=0;j<ncols;j++) {
5763927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
577ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[cols[j]];
5783927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
5793927de2eSStefano Zampini           my_dnz[i] += 1;
5803927de2eSStefano Zampini         } else { /* offdiag block */
5813927de2eSStefano Zampini           my_onz[i] += 1;
5823927de2eSStefano Zampini         }
5833927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
584d9a9e74cSStefano Zampini         if (issbaij && index_col != index_row) {
5853927de2eSStefano Zampini           owner = row_ownership[index_col];
5863927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
587d9a9e74cSStefano Zampini             my_dnz[cols[j]] += 1;
5883927de2eSStefano Zampini           } else {
589d9a9e74cSStefano Zampini             my_onz[cols[j]] += 1;
5903927de2eSStefano Zampini           }
5913927de2eSStefano Zampini         }
5923927de2eSStefano Zampini       }
5933927de2eSStefano Zampini       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
5943927de2eSStefano Zampini     }
5953927de2eSStefano Zampini   }
596ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
597ecf5a873SStefano Zampini   if (global_indices_c != global_indices_r) {
5987230de76SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
599ecf5a873SStefano Zampini   }
6003927de2eSStefano Zampini   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
601ecf5a873SStefano Zampini 
602ecf5a873SStefano Zampini   /* Reduce my_dnz and my_onz */
6033927de2eSStefano Zampini   if (maxreduce) {
6043927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
6053927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
606bb1015c3SStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
6073927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
6083927de2eSStefano Zampini   } else {
6093927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
6103927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
611bb1015c3SStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
6123927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
6133927de2eSStefano Zampini   }
6143927de2eSStefano Zampini   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
6153927de2eSStefano Zampini 
6163927de2eSStefano Zampini   /* Resize preallocation if overestimated */
6173927de2eSStefano Zampini   for (i=0;i<lrows;i++) {
6183927de2eSStefano Zampini     dnz[i] = PetscMin(dnz[i],lcols);
6193927de2eSStefano Zampini     onz[i] = PetscMin(onz[i],cols-lcols);
6203927de2eSStefano Zampini   }
6213927de2eSStefano Zampini   /* set preallocation */
6223927de2eSStefano Zampini   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
6233927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
6243927de2eSStefano Zampini     dnz[i] = dnz[i*bs]/bs;
6253927de2eSStefano Zampini     onz[i] = onz[i*bs]/bs;
6263927de2eSStefano Zampini   }
6273927de2eSStefano Zampini   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
6283927de2eSStefano Zampini   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
6293927de2eSStefano Zampini   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
6303927de2eSStefano Zampini   if (issbaij) {
6313927de2eSStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
6323927de2eSStefano Zampini   }
6333927de2eSStefano Zampini   PetscFunctionReturn(0);
6343927de2eSStefano Zampini }
6353927de2eSStefano Zampini 
6363927de2eSStefano Zampini #undef __FUNCT__
637b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
6387230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
639b7ce53b6SStefano Zampini {
640b7ce53b6SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
641d9a9e74cSStefano Zampini   Mat            local_mat;
642b7ce53b6SStefano Zampini   /* info on mat */
6433cfa4ea4SStefano Zampini   PetscInt       bs,rows,cols,lrows,lcols;
644b7ce53b6SStefano Zampini   PetscInt       local_rows,local_cols;
645686e3a49SStefano Zampini   PetscBool      isdense,issbaij,isseqaij;
6467c03b4e8SStefano Zampini   PetscMPIInt    nsubdomains;
647b7ce53b6SStefano Zampini   /* values insertion */
648b7ce53b6SStefano Zampini   PetscScalar    *array;
649b7ce53b6SStefano Zampini   /* work */
650b7ce53b6SStefano Zampini   PetscErrorCode ierr;
651b7ce53b6SStefano Zampini 
652b7ce53b6SStefano Zampini   PetscFunctionBegin;
653b7ce53b6SStefano Zampini   /* get info from mat */
6547c03b4e8SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
6557c03b4e8SStefano Zampini   if (nsubdomains == 1) {
6567c03b4e8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
6577c03b4e8SStefano Zampini       ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr);
6587c03b4e8SStefano Zampini     } else {
6597c03b4e8SStefano Zampini       ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
6607c03b4e8SStefano Zampini     }
6617c03b4e8SStefano Zampini     PetscFunctionReturn(0);
6627c03b4e8SStefano Zampini   }
663b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
664b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
6653cfa4ea4SStefano Zampini   ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr);
666b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
667b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
668686e3a49SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
669b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
670b7ce53b6SStefano Zampini 
671b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
672b7ce53b6SStefano Zampini     MatType     new_mat_type;
6733927de2eSStefano Zampini     PetscBool   issbaij_red;
674b7ce53b6SStefano Zampini 
675b7ce53b6SStefano Zampini     /* determining new matrix type */
676b2566f29SBarry Smith     ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
677b7ce53b6SStefano Zampini     if (issbaij_red) {
678b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
679b7ce53b6SStefano Zampini     } else {
680b7ce53b6SStefano Zampini       if (bs>1) {
681b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
682b7ce53b6SStefano Zampini       } else {
683b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
684b7ce53b6SStefano Zampini       }
685b7ce53b6SStefano Zampini     }
686b7ce53b6SStefano Zampini 
6873927de2eSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
6883cfa4ea4SStefano Zampini     ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr);
6893927de2eSStefano Zampini     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
6903927de2eSStefano Zampini     ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr);
6913927de2eSStefano Zampini     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
692b7ce53b6SStefano Zampini   } else {
6933cfa4ea4SStefano Zampini     PetscInt mbs,mrows,mcols,mlrows,mlcols;
694b7ce53b6SStefano Zampini     /* some checks */
695b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
696b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
6973cfa4ea4SStefano Zampini     ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr);
6986c4ed002SBarry Smith     if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
6996c4ed002SBarry Smith     if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
7006c4ed002SBarry Smith     if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows);
7016c4ed002SBarry Smith     if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols);
7026c4ed002SBarry Smith     if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
703b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
704b7ce53b6SStefano Zampini   }
705d9a9e74cSStefano Zampini 
706d9a9e74cSStefano Zampini   if (issbaij) {
707d9a9e74cSStefano Zampini     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr);
708d9a9e74cSStefano Zampini   } else {
709d9a9e74cSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
710d9a9e74cSStefano Zampini     local_mat = matis->A;
711d9a9e74cSStefano Zampini   }
712686e3a49SStefano Zampini 
713b7ce53b6SStefano Zampini   /* Set values */
714ecf5a873SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr);
715b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
716ecf5a873SStefano Zampini     PetscInt i,*dummy_rows,*dummy_cols;
717ecf5a873SStefano Zampini 
718ecf5a873SStefano Zampini     if (local_rows != local_cols) {
719ecf5a873SStefano Zampini       ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr);
720ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
721ecf5a873SStefano Zampini       for (i=0;i<local_cols;i++) dummy_cols[i] = i;
722ecf5a873SStefano Zampini     } else {
723ecf5a873SStefano Zampini       ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr);
724ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
725ecf5a873SStefano Zampini       dummy_cols = dummy_rows;
726ecf5a873SStefano Zampini     }
727b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
728d9a9e74cSStefano Zampini     ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr);
729ecf5a873SStefano Zampini     ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr);
730d9a9e74cSStefano Zampini     ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr);
731ecf5a873SStefano Zampini     if (dummy_rows != dummy_cols) {
732ecf5a873SStefano Zampini       ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr);
733ecf5a873SStefano Zampini     } else {
734ecf5a873SStefano Zampini       ierr = PetscFree(dummy_rows);CHKERRQ(ierr);
735ecf5a873SStefano Zampini     }
736686e3a49SStefano Zampini   } else if (isseqaij) {
737ecf5a873SStefano Zampini     PetscInt  i,nvtxs,*xadj,*adjncy;
738686e3a49SStefano Zampini     PetscBool done;
739686e3a49SStefano Zampini 
740d9a9e74cSStefano Zampini     ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
741bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__);
742d9a9e74cSStefano Zampini     ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr);
743686e3a49SStefano Zampini     for (i=0;i<nvtxs;i++) {
744ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr);
745686e3a49SStefano Zampini     }
746d9a9e74cSStefano Zampini     ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
747bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__);
748d9a9e74cSStefano Zampini     ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr);
749686e3a49SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
750ecf5a873SStefano Zampini     PetscInt i;
751c0962df8SStefano Zampini 
752686e3a49SStefano Zampini     for (i=0;i<local_rows;i++) {
753686e3a49SStefano Zampini       PetscInt       j;
754ecf5a873SStefano Zampini       const PetscInt *local_indices_cols;
755686e3a49SStefano Zampini 
756ecf5a873SStefano Zampini       ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
757ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr);
758ecf5a873SStefano Zampini       ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
759686e3a49SStefano Zampini     }
760b7ce53b6SStefano Zampini   }
761b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
762d9a9e74cSStefano Zampini   ierr = MatDestroy(&local_mat);CHKERRQ(ierr);
763b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
764b7ce53b6SStefano Zampini   if (isdense) {
765b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
766b7ce53b6SStefano Zampini   }
767b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
768b7ce53b6SStefano Zampini }
769b7ce53b6SStefano Zampini 
770b7ce53b6SStefano Zampini #undef __FUNCT__
771b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
772b7ce53b6SStefano Zampini /*@
773b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
774b7ce53b6SStefano Zampini 
775b7ce53b6SStefano Zampini   Input Parameter:
776b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
777b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
778b7ce53b6SStefano Zampini 
779b7ce53b6SStefano Zampini   Output Parameter:
780b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
781b7ce53b6SStefano Zampini 
782b7ce53b6SStefano Zampini   Level: developer
783b7ce53b6SStefano Zampini 
784eb82efa4SStefano Zampini   Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested.
785b7ce53b6SStefano Zampini 
786b7ce53b6SStefano Zampini .seealso: MATIS
787b7ce53b6SStefano Zampini @*/
788b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
789b7ce53b6SStefano Zampini {
790b7ce53b6SStefano Zampini   PetscErrorCode ierr;
791b7ce53b6SStefano Zampini 
792b7ce53b6SStefano Zampini   PetscFunctionBegin;
793b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
794b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
795b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
796b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
797b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
798b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
7996c4ed002SBarry Smith     if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix");
800b7ce53b6SStefano Zampini   }
801b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
802b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
803b7ce53b6SStefano Zampini }
804b7ce53b6SStefano Zampini 
805b7ce53b6SStefano Zampini #undef __FUNCT__
806ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
807ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
808ad6194a2SStefano Zampini {
809ad6194a2SStefano Zampini   PetscErrorCode ierr;
810ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
811ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
812ad6194a2SStefano Zampini   Mat            B,localmat;
813ad6194a2SStefano Zampini 
814ad6194a2SStefano Zampini   PetscFunctionBegin;
815ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
816ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
817ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
818e176bc59SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr);
819ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
820ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
821b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
822ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
823ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
824ad6194a2SStefano Zampini   *newmat = B;
825ad6194a2SStefano Zampini   PetscFunctionReturn(0);
826ad6194a2SStefano Zampini }
827ad6194a2SStefano Zampini 
828ad6194a2SStefano Zampini #undef __FUNCT__
82969796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
830a8116848SStefano Zampini static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
83169796d55SStefano Zampini {
83269796d55SStefano Zampini   PetscErrorCode ierr;
83369796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
83469796d55SStefano Zampini   PetscBool      local_sym;
83569796d55SStefano Zampini 
83669796d55SStefano Zampini   PetscFunctionBegin;
83769796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
838b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
83969796d55SStefano Zampini   PetscFunctionReturn(0);
84069796d55SStefano Zampini }
84169796d55SStefano Zampini 
84269796d55SStefano Zampini #undef __FUNCT__
84369796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
844a8116848SStefano Zampini static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
84569796d55SStefano Zampini {
84669796d55SStefano Zampini   PetscErrorCode ierr;
84769796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
84869796d55SStefano Zampini   PetscBool      local_sym;
84969796d55SStefano Zampini 
85069796d55SStefano Zampini   PetscFunctionBegin;
85169796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
852b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
85369796d55SStefano Zampini   PetscFunctionReturn(0);
85469796d55SStefano Zampini }
85569796d55SStefano Zampini 
85669796d55SStefano Zampini #undef __FUNCT__
857b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
858a8116848SStefano Zampini static PetscErrorCode MatDestroy_IS(Mat A)
859b4319ba4SBarry Smith {
860dfbe8321SBarry Smith   PetscErrorCode ierr;
861b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
862b4319ba4SBarry Smith 
863b4319ba4SBarry Smith   PetscFunctionBegin;
8646bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
865e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr);
866e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr);
8676bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
8686bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
8693fd1c9e7SStefano Zampini   ierr = VecDestroy(&b->counter);CHKERRQ(ierr);
870a8116848SStefano Zampini   ierr = ISDestroy(&b->getsub_ris);CHKERRQ(ierr);
871a8116848SStefano Zampini   ierr = ISDestroy(&b->getsub_cis);CHKERRQ(ierr);
872a8116848SStefano Zampini   if (b->sf != b->csf) {
873a8116848SStefano Zampini     ierr = PetscSFDestroy(&b->csf);CHKERRQ(ierr);
874a8116848SStefano Zampini     ierr = PetscFree2(b->csf_rootdata,b->csf_leafdata);CHKERRQ(ierr);
875a8116848SStefano Zampini   }
87628f4e0baSStefano Zampini   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
87728f4e0baSStefano Zampini   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
878bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
879dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
880bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
881b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
882b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
8832e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
884b4319ba4SBarry Smith   PetscFunctionReturn(0);
885b4319ba4SBarry Smith }
886b4319ba4SBarry Smith 
887b4319ba4SBarry Smith #undef __FUNCT__
888b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
889a8116848SStefano Zampini static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
890b4319ba4SBarry Smith {
891dfbe8321SBarry Smith   PetscErrorCode ierr;
892b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
893b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
894b4319ba4SBarry Smith 
895b4319ba4SBarry Smith   PetscFunctionBegin;
896b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
897e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
898e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
899b4319ba4SBarry Smith 
900b4319ba4SBarry Smith   /* multiply the local matrix */
901b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
902b4319ba4SBarry Smith 
903b4319ba4SBarry Smith   /* scatter product back into global memory */
9042dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
905e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
906e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
907b4319ba4SBarry Smith   PetscFunctionReturn(0);
908b4319ba4SBarry Smith }
909b4319ba4SBarry Smith 
910b4319ba4SBarry Smith #undef __FUNCT__
9112e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
912a8116848SStefano Zampini static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
9132e74eeadSLisandro Dalcin {
914650997f4SStefano Zampini   Vec            temp_vec;
9152e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9162e74eeadSLisandro Dalcin 
9172e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
918650997f4SStefano Zampini   if (v3 != v2) {
919650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
920650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
921650997f4SStefano Zampini   } else {
922650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
923650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
924650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
925650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
926650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
927650997f4SStefano Zampini   }
9282e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9292e74eeadSLisandro Dalcin }
9302e74eeadSLisandro Dalcin 
9312e74eeadSLisandro Dalcin #undef __FUNCT__
9322e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
933a8116848SStefano Zampini static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x)
9342e74eeadSLisandro Dalcin {
9352e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9362e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9372e74eeadSLisandro Dalcin 
938e176bc59SStefano Zampini   PetscFunctionBegin;
9392e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
940e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
941e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9422e74eeadSLisandro Dalcin 
9432e74eeadSLisandro Dalcin   /* multiply the local matrix */
944e176bc59SStefano Zampini   ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr);
9452e74eeadSLisandro Dalcin 
9462e74eeadSLisandro Dalcin   /* scatter product back into global vector */
947e176bc59SStefano Zampini   ierr = VecSet(x,0);CHKERRQ(ierr);
948e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
949e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9502e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9512e74eeadSLisandro Dalcin }
9522e74eeadSLisandro Dalcin 
9532e74eeadSLisandro Dalcin #undef __FUNCT__
9542e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
955a8116848SStefano Zampini static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
9562e74eeadSLisandro Dalcin {
957650997f4SStefano Zampini   Vec            temp_vec;
9582e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9592e74eeadSLisandro Dalcin 
9602e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
961650997f4SStefano Zampini   if (v3 != v2) {
962650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
963650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
964650997f4SStefano Zampini   } else {
965650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
966650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
967650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
968650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
969650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
970650997f4SStefano Zampini   }
9712e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9722e74eeadSLisandro Dalcin }
9732e74eeadSLisandro Dalcin 
9742e74eeadSLisandro Dalcin #undef __FUNCT__
975b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
976a8116848SStefano Zampini static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
977b4319ba4SBarry Smith {
978b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
979dfbe8321SBarry Smith   PetscErrorCode ierr;
980b4319ba4SBarry Smith   PetscViewer    sviewer;
981b4319ba4SBarry Smith 
982b4319ba4SBarry Smith   PetscFunctionBegin;
9833f08860eSBarry Smith   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
984b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
9853f08860eSBarry Smith   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
9866e520ac8SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
987b4319ba4SBarry Smith   PetscFunctionReturn(0);
988b4319ba4SBarry Smith }
989b4319ba4SBarry Smith 
990b4319ba4SBarry Smith #undef __FUNCT__
991b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
992a8116848SStefano Zampini static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
993b4319ba4SBarry Smith {
994dfbe8321SBarry Smith   PetscErrorCode ierr;
995e176bc59SStefano Zampini   PetscInt       nr,rbs,nc,cbs;
996b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
997b4319ba4SBarry Smith   IS             from,to;
998e176bc59SStefano Zampini   Vec            cglobal,rglobal;
999b4319ba4SBarry Smith 
1000b4319ba4SBarry Smith   PetscFunctionBegin;
1001784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
1002e176bc59SStefano Zampini   PetscCheckSameComm(A,1,cmapping,3);
10033bbff08aSStefano Zampini   /* Destroy any previous data */
100470cf5478SStefano Zampini   ierr = VecDestroy(&is->x);CHKERRQ(ierr);
100570cf5478SStefano Zampini   ierr = VecDestroy(&is->y);CHKERRQ(ierr);
10063fd1c9e7SStefano Zampini   ierr = VecDestroy(&is->counter);CHKERRQ(ierr);
1007e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr);
1008e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr);
10091c47cb0fSStefano Zampini   ierr = MatDestroy(&is->A);CHKERRQ(ierr);
101028f4e0baSStefano Zampini   ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
101128f4e0baSStefano Zampini   ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
10123bbff08aSStefano Zampini 
10133bbff08aSStefano Zampini   /* Setup Layout and set local to global maps */
1014fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1015fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1016fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
1017fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
1018b4319ba4SBarry Smith 
1019b4319ba4SBarry Smith   /* Create the local matrix A */
1020e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr);
1021e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr);
1022e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr);
1023e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr);
1024f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
1025e176bc59SStefano Zampini   ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr);
1026e176bc59SStefano Zampini   ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr);
1027e176bc59SStefano Zampini   ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr);
1028ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
1029ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
1030b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
1031b4319ba4SBarry Smith 
1032f26d0771SStefano Zampini   if (!is->islocalref) { /* setup scatters and local vectors for MatMult */
1033b4319ba4SBarry Smith     /* Create the local work vectors */
10343bbff08aSStefano Zampini     ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr);
1035b4319ba4SBarry Smith 
1036e176bc59SStefano Zampini     /* setup the global to local scatters */
1037e176bc59SStefano Zampini     ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr);
1038e176bc59SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr);
1039784ac674SJed Brown     ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
1040e176bc59SStefano Zampini     ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr);
1041e176bc59SStefano Zampini     if (rmapping != cmapping) {
1042e176bc59SStefano Zampini       ierr = ISDestroy(&to);CHKERRQ(ierr);
1043e176bc59SStefano Zampini       ierr = ISDestroy(&from);CHKERRQ(ierr);
1044e176bc59SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr);
1045e176bc59SStefano Zampini       ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr);
1046e176bc59SStefano Zampini       ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr);
1047e176bc59SStefano Zampini     } else {
1048e176bc59SStefano Zampini       ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr);
1049e176bc59SStefano Zampini       is->cctx = is->rctx;
1050e176bc59SStefano Zampini     }
10513fd1c9e7SStefano Zampini 
10523fd1c9e7SStefano Zampini     /* interface counter vector (local) */
10533fd1c9e7SStefano Zampini     ierr = VecDuplicate(is->y,&is->counter);CHKERRQ(ierr);
10543fd1c9e7SStefano Zampini     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
10553fd1c9e7SStefano Zampini     ierr = VecScatterBegin(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10563fd1c9e7SStefano Zampini     ierr = VecScatterEnd(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10573fd1c9e7SStefano Zampini     ierr = VecScatterBegin(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10583fd1c9e7SStefano Zampini     ierr = VecScatterEnd(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10593fd1c9e7SStefano Zampini 
10603fd1c9e7SStefano Zampini     /* free workspace */
1061e176bc59SStefano Zampini     ierr = VecDestroy(&rglobal);CHKERRQ(ierr);
1062e176bc59SStefano Zampini     ierr = VecDestroy(&cglobal);CHKERRQ(ierr);
10636bf464f9SBarry Smith     ierr = ISDestroy(&to);CHKERRQ(ierr);
10646bf464f9SBarry Smith     ierr = ISDestroy(&from);CHKERRQ(ierr);
1065f26d0771SStefano Zampini   }
106648ff6bf3SStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
1067b4319ba4SBarry Smith   PetscFunctionReturn(0);
1068b4319ba4SBarry Smith }
1069b4319ba4SBarry Smith 
10702e74eeadSLisandro Dalcin #undef __FUNCT__
10712e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
1072a8116848SStefano Zampini static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
10732e74eeadSLisandro Dalcin {
10742e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
10752e74eeadSLisandro Dalcin   PetscErrorCode ierr;
107697563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
107797563a80SStefano Zampini   PetscInt       i,zm,zn;
107897563a80SStefano Zampini #endif
1079f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
10802e74eeadSLisandro Dalcin 
10812e74eeadSLisandro Dalcin   PetscFunctionBegin;
10822e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
1083f26d0771SStefano 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);
108497563a80SStefano Zampini   /* count negative indices */
108597563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
108697563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
10872e74eeadSLisandro Dalcin #endif
108897563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
108997563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
109097563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
109197563a80SStefano Zampini   /* count negative indices (should be the same as before) */
109297563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
109397563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
109497563a80SStefano 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");
109597563a80SStefano 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");
109697563a80SStefano Zampini #endif
10972e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
10982e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
10992e74eeadSLisandro Dalcin }
11002e74eeadSLisandro Dalcin 
1101b4319ba4SBarry Smith #undef __FUNCT__
110297563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS"
1103a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
110497563a80SStefano Zampini {
110597563a80SStefano Zampini   Mat_IS         *is = (Mat_IS*)mat->data;
110697563a80SStefano Zampini   PetscErrorCode ierr;
110797563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
110897563a80SStefano Zampini   PetscInt       i,zm,zn;
110997563a80SStefano Zampini #endif
1110f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
111197563a80SStefano Zampini 
111297563a80SStefano Zampini   PetscFunctionBegin;
111397563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
1114f26d0771SStefano 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);
111597563a80SStefano Zampini   /* count negative indices */
111697563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
111797563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
111897563a80SStefano Zampini #endif
111997563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
112097563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
112197563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
112297563a80SStefano Zampini   /* count negative indices (should be the same as before) */
112397563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
112497563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
112597563a80SStefano 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");
112697563a80SStefano 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");
112797563a80SStefano Zampini #endif
112897563a80SStefano Zampini   ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
112997563a80SStefano Zampini   PetscFunctionReturn(0);
113097563a80SStefano Zampini }
113197563a80SStefano Zampini 
113297563a80SStefano Zampini #undef __FUNCT__
1133b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
1134a8116848SStefano Zampini static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
1135b4319ba4SBarry Smith {
1136dfbe8321SBarry Smith   PetscErrorCode ierr;
1137b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1138b4319ba4SBarry Smith 
1139b4319ba4SBarry Smith   PetscFunctionBegin;
1140b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
1141b4319ba4SBarry Smith   PetscFunctionReturn(0);
1142b4319ba4SBarry Smith }
1143b4319ba4SBarry Smith 
1144b4319ba4SBarry Smith #undef __FUNCT__
1145f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
1146a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
1147f0006bf2SLisandro Dalcin {
1148f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
1149f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
1150f0006bf2SLisandro Dalcin 
1151f0006bf2SLisandro Dalcin   PetscFunctionBegin;
1152f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
1153f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
1154f0006bf2SLisandro Dalcin }
1155f0006bf2SLisandro Dalcin 
1156f0006bf2SLisandro Dalcin #undef __FUNCT__
1157*f0ae7da4SStefano Zampini #define __FUNCT__ "MatISZeroRowsColumnsLocal_Private"
1158*f0ae7da4SStefano Zampini static PetscErrorCode MatISZeroRowsColumnsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,PetscBool columns)
1159*f0ae7da4SStefano Zampini {
1160*f0ae7da4SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
1161*f0ae7da4SStefano Zampini   PetscErrorCode ierr;
1162*f0ae7da4SStefano Zampini 
1163*f0ae7da4SStefano Zampini   PetscFunctionBegin;
1164*f0ae7da4SStefano Zampini   if (!n) {
1165*f0ae7da4SStefano Zampini     is->pure_neumann = PETSC_TRUE;
1166*f0ae7da4SStefano Zampini   } else {
1167*f0ae7da4SStefano Zampini     PetscInt i;
1168*f0ae7da4SStefano Zampini     is->pure_neumann = PETSC_FALSE;
1169*f0ae7da4SStefano Zampini 
1170*f0ae7da4SStefano Zampini     if (columns) {
1171*f0ae7da4SStefano Zampini       ierr = MatZeroRowsColumns(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
1172*f0ae7da4SStefano Zampini     } else {
1173*f0ae7da4SStefano Zampini       ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
1174*f0ae7da4SStefano Zampini     }
1175*f0ae7da4SStefano Zampini     if (diag != 0.) {
1176*f0ae7da4SStefano Zampini       const PetscScalar *array;
1177*f0ae7da4SStefano Zampini       ierr = VecGetArrayRead(is->counter,&array);CHKERRQ(ierr);
1178*f0ae7da4SStefano Zampini       for (i=0; i<n; i++) {
1179*f0ae7da4SStefano Zampini         ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
1180*f0ae7da4SStefano Zampini       }
1181*f0ae7da4SStefano Zampini       ierr = VecRestoreArrayRead(is->counter,&array);CHKERRQ(ierr);
1182*f0ae7da4SStefano Zampini     }
1183*f0ae7da4SStefano Zampini     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1184*f0ae7da4SStefano Zampini     ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1185*f0ae7da4SStefano Zampini   }
1186*f0ae7da4SStefano Zampini   PetscFunctionReturn(0);
1187*f0ae7da4SStefano Zampini }
1188*f0ae7da4SStefano Zampini 
1189*f0ae7da4SStefano Zampini #undef __FUNCT__
1190*f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_Private_IS"
1191*f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Private_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b,PetscBool columns)
11922e74eeadSLisandro Dalcin {
11936e520ac8SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
11946e520ac8SStefano Zampini   PetscInt       nr,nl,len,i;
11956e520ac8SStefano Zampini   PetscInt       *lrows;
11962e74eeadSLisandro Dalcin   PetscErrorCode ierr;
11972e74eeadSLisandro Dalcin 
11982e74eeadSLisandro Dalcin   PetscFunctionBegin;
1199*f0ae7da4SStefano Zampini #if defined(PETSC_USE_DEBUG)
1200*f0ae7da4SStefano Zampini   if (columns || diag != 0. || (x && b)) {
1201*f0ae7da4SStefano Zampini     PetscBool cong;
1202*f0ae7da4SStefano Zampini     ierr = PetscLayoutCompare(A->rmap,A->cmap,&cong);CHKERRQ(ierr);
1203*f0ae7da4SStefano 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");
1204*f0ae7da4SStefano 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");
1205*f0ae7da4SStefano Zampini     if (!cong && x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"A->rmap and A->cmap need to be congruent");
1206*f0ae7da4SStefano Zampini   }
1207*f0ae7da4SStefano Zampini #endif
12086e520ac8SStefano Zampini   /* get locally owned rows */
1209*f0ae7da4SStefano Zampini   ierr = PetscLayoutMapLocal_Private(A->rmap,n,rows,&len,&lrows,NULL);CHKERRQ(ierr);
12106e520ac8SStefano Zampini   /* fix right hand side if needed */
12116e520ac8SStefano Zampini   if (x && b) {
12126e520ac8SStefano Zampini     const PetscScalar *xx;
12136e520ac8SStefano Zampini     PetscScalar       *bb;
12146e520ac8SStefano Zampini 
12156e520ac8SStefano Zampini     ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr);
12166e520ac8SStefano Zampini     ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
12176e520ac8SStefano Zampini     for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]];
12186e520ac8SStefano Zampini     ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr);
12196e520ac8SStefano Zampini     ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr);
12202e74eeadSLisandro Dalcin   }
12216e520ac8SStefano Zampini   /* get rows associated to the local matrices */
12226e520ac8SStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
12236e520ac8SStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
12246e520ac8SStefano Zampini   }
12256e520ac8SStefano Zampini   ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr);
12266e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr);
12276e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
12286e520ac8SStefano Zampini   for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1;
12296e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
12306e520ac8SStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
12316e520ac8SStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
12326e520ac8SStefano Zampini   ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr);
12336e520ac8SStefano Zampini   for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i;
1234*f0ae7da4SStefano Zampini   ierr = MatISZeroRowsColumnsLocal_Private(A,nr,lrows,diag,columns);CHKERRQ(ierr);
12356e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
12362e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
12372e74eeadSLisandro Dalcin }
12382e74eeadSLisandro Dalcin 
12392e74eeadSLisandro Dalcin #undef __FUNCT__
1240*f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRows_IS"
1241*f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1242b4319ba4SBarry Smith {
1243dfbe8321SBarry Smith   PetscErrorCode ierr;
1244b4319ba4SBarry Smith 
1245b4319ba4SBarry Smith   PetscFunctionBegin;
1246*f0ae7da4SStefano Zampini   ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_FALSE);CHKERRQ(ierr);
1247*f0ae7da4SStefano Zampini   PetscFunctionReturn(0);
1248*f0ae7da4SStefano Zampini }
12492205254eSKarl Rupp 
1250*f0ae7da4SStefano Zampini #undef __FUNCT__
1251*f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_IS"
1252*f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1253*f0ae7da4SStefano Zampini {
1254*f0ae7da4SStefano Zampini   PetscErrorCode ierr;
1255*f0ae7da4SStefano Zampini 
1256*f0ae7da4SStefano Zampini   PetscFunctionBegin;
1257*f0ae7da4SStefano Zampini   ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_TRUE);CHKERRQ(ierr);
1258b4319ba4SBarry Smith   PetscFunctionReturn(0);
1259b4319ba4SBarry Smith }
1260b4319ba4SBarry Smith 
1261b4319ba4SBarry Smith #undef __FUNCT__
1262b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
1263a8116848SStefano Zampini static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
1264b4319ba4SBarry Smith {
1265b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1266dfbe8321SBarry Smith   PetscErrorCode ierr;
1267b4319ba4SBarry Smith 
1268b4319ba4SBarry Smith   PetscFunctionBegin;
1269b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
1270b4319ba4SBarry Smith   PetscFunctionReturn(0);
1271b4319ba4SBarry Smith }
1272b4319ba4SBarry Smith 
1273b4319ba4SBarry Smith #undef __FUNCT__
1274b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
1275a8116848SStefano Zampini static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
1276b4319ba4SBarry Smith {
1277b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1278dfbe8321SBarry Smith   PetscErrorCode ierr;
1279b4319ba4SBarry Smith 
1280b4319ba4SBarry Smith   PetscFunctionBegin;
1281b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
1282b4319ba4SBarry Smith   PetscFunctionReturn(0);
1283b4319ba4SBarry Smith }
1284b4319ba4SBarry Smith 
1285b4319ba4SBarry Smith #undef __FUNCT__
1286b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
1287a8116848SStefano Zampini static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
1288b4319ba4SBarry Smith {
1289b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
1290b4319ba4SBarry Smith 
1291b4319ba4SBarry Smith   PetscFunctionBegin;
1292b4319ba4SBarry Smith   *local = is->A;
1293b4319ba4SBarry Smith   PetscFunctionReturn(0);
1294b4319ba4SBarry Smith }
1295b4319ba4SBarry Smith 
1296b4319ba4SBarry Smith #undef __FUNCT__
1297b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
1298b4319ba4SBarry Smith /*@
1299b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
1300b4319ba4SBarry Smith 
1301b4319ba4SBarry Smith   Input Parameter:
1302b4319ba4SBarry Smith .  mat - the matrix
1303b4319ba4SBarry Smith 
1304b4319ba4SBarry Smith   Output Parameter:
1305eb82efa4SStefano Zampini .  local - the local matrix
1306b4319ba4SBarry Smith 
1307b4319ba4SBarry Smith   Level: advanced
1308b4319ba4SBarry Smith 
1309b4319ba4SBarry Smith   Notes:
1310b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
1311b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
1312b4319ba4SBarry Smith   of the MatSetValues() operation.
1313b4319ba4SBarry Smith 
1314b4319ba4SBarry Smith .seealso: MATIS
1315b4319ba4SBarry Smith @*/
13167087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
1317b4319ba4SBarry Smith {
13184ac538c5SBarry Smith   PetscErrorCode ierr;
1319b4319ba4SBarry Smith 
1320b4319ba4SBarry Smith   PetscFunctionBegin;
13210700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1322b4319ba4SBarry Smith   PetscValidPointer(local,2);
13234ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
1324b4319ba4SBarry Smith   PetscFunctionReturn(0);
1325b4319ba4SBarry Smith }
1326b4319ba4SBarry Smith 
13273b03a366Sstefano_zampini #undef __FUNCT__
13283b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
1329a8116848SStefano Zampini static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
13303b03a366Sstefano_zampini {
13313b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
13323b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
13333b03a366Sstefano_zampini   PetscErrorCode ierr;
13343b03a366Sstefano_zampini 
13353b03a366Sstefano_zampini   PetscFunctionBegin;
13364e4c7dbeSStefano Zampini   if (is->A) {
13373b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
13383b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
1339*f0ae7da4SStefano 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);
13404e4c7dbeSStefano Zampini   }
13413b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
13423b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
13433b03a366Sstefano_zampini   is->A = local;
13443b03a366Sstefano_zampini   PetscFunctionReturn(0);
13453b03a366Sstefano_zampini }
13463b03a366Sstefano_zampini 
13473b03a366Sstefano_zampini #undef __FUNCT__
13483b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
13493b03a366Sstefano_zampini /*@
1350eb82efa4SStefano Zampini     MatISSetLocalMat - Replace the local matrix stored inside a MATIS object.
13513b03a366Sstefano_zampini 
13523b03a366Sstefano_zampini   Input Parameter:
13533b03a366Sstefano_zampini .  mat - the matrix
1354eb82efa4SStefano Zampini .  local - the local matrix
13553b03a366Sstefano_zampini 
13563b03a366Sstefano_zampini   Output Parameter:
13573b03a366Sstefano_zampini 
13583b03a366Sstefano_zampini   Level: advanced
13593b03a366Sstefano_zampini 
13603b03a366Sstefano_zampini   Notes:
13613b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
13623b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
13633b03a366Sstefano_zampini 
13643b03a366Sstefano_zampini .seealso: MATIS
13653b03a366Sstefano_zampini @*/
13663b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
13673b03a366Sstefano_zampini {
13683b03a366Sstefano_zampini   PetscErrorCode ierr;
13693b03a366Sstefano_zampini 
13703b03a366Sstefano_zampini   PetscFunctionBegin;
13713b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1372b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
13733b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
13743b03a366Sstefano_zampini   PetscFunctionReturn(0);
13753b03a366Sstefano_zampini }
13763b03a366Sstefano_zampini 
13776726f965SBarry Smith #undef __FUNCT__
13786726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
1379a8116848SStefano Zampini static PetscErrorCode MatZeroEntries_IS(Mat A)
13806726f965SBarry Smith {
13816726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
13826726f965SBarry Smith   PetscErrorCode ierr;
13836726f965SBarry Smith 
13846726f965SBarry Smith   PetscFunctionBegin;
13856726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
13866726f965SBarry Smith   PetscFunctionReturn(0);
13876726f965SBarry Smith }
13886726f965SBarry Smith 
13896726f965SBarry Smith #undef __FUNCT__
13902e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
1391a8116848SStefano Zampini static PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
13922e74eeadSLisandro Dalcin {
13932e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
13942e74eeadSLisandro Dalcin   PetscErrorCode ierr;
13952e74eeadSLisandro Dalcin 
13962e74eeadSLisandro Dalcin   PetscFunctionBegin;
13972e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
13982e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
13992e74eeadSLisandro Dalcin }
14002e74eeadSLisandro Dalcin 
14012e74eeadSLisandro Dalcin #undef __FUNCT__
14022e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
1403a8116848SStefano Zampini static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
14042e74eeadSLisandro Dalcin {
14052e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
14062e74eeadSLisandro Dalcin   PetscErrorCode ierr;
14072e74eeadSLisandro Dalcin 
14082e74eeadSLisandro Dalcin   PetscFunctionBegin;
14092e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
1410e176bc59SStefano Zampini   ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr);
14112e74eeadSLisandro Dalcin 
14122e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
14132e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
1414e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1415e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
14162e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
14172e74eeadSLisandro Dalcin }
14182e74eeadSLisandro Dalcin 
14192e74eeadSLisandro Dalcin #undef __FUNCT__
14206726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
1421a8116848SStefano Zampini static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
14226726f965SBarry Smith {
14236726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
14246726f965SBarry Smith   PetscErrorCode ierr;
14256726f965SBarry Smith 
14266726f965SBarry Smith   PetscFunctionBegin;
14274e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
14286726f965SBarry Smith   PetscFunctionReturn(0);
14296726f965SBarry Smith }
14306726f965SBarry Smith 
1431284134d9SBarry Smith #undef __FUNCT__
1432f26d0771SStefano Zampini #define __FUNCT__ "MatAXPY_IS"
1433f26d0771SStefano Zampini static PetscErrorCode MatAXPY_IS(Mat Y,PetscScalar a,Mat X,MatStructure str)
1434f26d0771SStefano Zampini {
1435f26d0771SStefano Zampini   Mat_IS         *y = (Mat_IS*)Y->data;
1436f26d0771SStefano Zampini   Mat_IS         *x;
1437f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1438f26d0771SStefano Zampini   PetscBool      ismatis;
1439f26d0771SStefano Zampini #endif
1440f26d0771SStefano Zampini   PetscErrorCode ierr;
1441f26d0771SStefano Zampini 
1442f26d0771SStefano Zampini   PetscFunctionBegin;
1443f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1444f26d0771SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)X,MATIS,&ismatis);CHKERRQ(ierr);
1445f26d0771SStefano Zampini   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS");
1446f26d0771SStefano Zampini #endif
1447f26d0771SStefano Zampini   x = (Mat_IS*)X->data;
1448f26d0771SStefano Zampini   ierr = MatAXPY(y->A,a,x->A,str);CHKERRQ(ierr);
1449f26d0771SStefano Zampini   PetscFunctionReturn(0);
1450f26d0771SStefano Zampini }
1451f26d0771SStefano Zampini 
1452f26d0771SStefano Zampini #undef __FUNCT__
1453f26d0771SStefano Zampini #define __FUNCT__ "MatGetLocalSubMatrix_IS"
1454f26d0771SStefano Zampini static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat)
1455f26d0771SStefano Zampini {
1456f26d0771SStefano Zampini   Mat                    lA;
1457f26d0771SStefano Zampini   Mat_IS                 *matis;
1458f26d0771SStefano Zampini   ISLocalToGlobalMapping rl2g,cl2g;
1459f26d0771SStefano Zampini   IS                     is;
1460f26d0771SStefano Zampini   const PetscInt         *rg,*rl;
1461f26d0771SStefano Zampini   PetscInt               nrg;
1462f26d0771SStefano Zampini   PetscInt               N,M,nrl,i,*idxs;
1463f26d0771SStefano Zampini   PetscErrorCode         ierr;
1464f26d0771SStefano Zampini 
1465f26d0771SStefano Zampini   PetscFunctionBegin;
1466f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr);
1467f26d0771SStefano Zampini   ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr);
1468f26d0771SStefano Zampini   ierr = ISGetIndices(row,&rl);CHKERRQ(ierr);
1469f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr);
1470f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1471*f0ae7da4SStefano 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);
1472f26d0771SStefano Zampini #endif
1473f26d0771SStefano Zampini   ierr = PetscMalloc1(nrg,&idxs);CHKERRQ(ierr);
1474f26d0771SStefano Zampini   /* map from [0,nrl) to row */
1475f26d0771SStefano Zampini   for (i=0;i<nrl;i++) idxs[i] = rl[i];
1476f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1477f26d0771SStefano Zampini   for (i=nrl;i<nrg;i++) idxs[i] = nrg;
1478f26d0771SStefano Zampini #else
1479f26d0771SStefano Zampini   for (i=nrl;i<nrg;i++) idxs[i] = -1;
1480f26d0771SStefano Zampini #endif
1481f26d0771SStefano Zampini   ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr);
1482f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr);
1483f26d0771SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
1484f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
1485f26d0771SStefano Zampini   ierr = ISDestroy(&is);CHKERRQ(ierr);
1486f26d0771SStefano Zampini   /* compute new l2g map for columns */
1487f26d0771SStefano Zampini   if (col != row || A->rmap->mapping != A->cmap->mapping) {
1488f26d0771SStefano Zampini     const PetscInt *cg,*cl;
1489f26d0771SStefano Zampini     PetscInt       ncg;
1490f26d0771SStefano Zampini     PetscInt       ncl;
1491f26d0771SStefano Zampini 
1492f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr);
1493f26d0771SStefano Zampini     ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr);
1494f26d0771SStefano Zampini     ierr = ISGetIndices(col,&cl);CHKERRQ(ierr);
1495f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr);
1496f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1497*f0ae7da4SStefano 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);
1498f26d0771SStefano Zampini #endif
1499f26d0771SStefano Zampini     ierr = PetscMalloc1(ncg,&idxs);CHKERRQ(ierr);
1500f26d0771SStefano Zampini     /* map from [0,ncl) to col */
1501f26d0771SStefano Zampini     for (i=0;i<ncl;i++) idxs[i] = cl[i];
1502f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1503f26d0771SStefano Zampini     for (i=ncl;i<ncg;i++) idxs[i] = ncg;
1504f26d0771SStefano Zampini #else
1505f26d0771SStefano Zampini     for (i=ncl;i<ncg;i++) idxs[i] = -1;
1506f26d0771SStefano Zampini #endif
1507f26d0771SStefano Zampini     ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr);
1508f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr);
1509f26d0771SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
1510f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
1511f26d0771SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
1512f26d0771SStefano Zampini   } else {
1513f26d0771SStefano Zampini     ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr);
1514f26d0771SStefano Zampini     cl2g = rl2g;
1515f26d0771SStefano Zampini   }
1516f26d0771SStefano Zampini   /* create the MATIS submatrix */
1517f26d0771SStefano Zampini   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1518f26d0771SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr);
1519f26d0771SStefano Zampini   ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1520f26d0771SStefano Zampini   ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr);
1521b0aa3428SStefano Zampini   matis = (Mat_IS*)((*submat)->data);
1522f26d0771SStefano Zampini   matis->islocalref = PETSC_TRUE;
1523f26d0771SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr);
1524f26d0771SStefano Zampini   ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr);
1525f26d0771SStefano Zampini   ierr = MatISSetLocalMat(*submat,lA);CHKERRQ(ierr);
1526f26d0771SStefano Zampini   ierr = MatSetUp(*submat);CHKERRQ(ierr);
1527f26d0771SStefano Zampini   ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1528f26d0771SStefano Zampini   ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1529f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
1530f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
1531f26d0771SStefano Zampini   /* remove unsupported ops */
1532f26d0771SStefano Zampini   ierr = PetscMemzero((*submat)->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1533f26d0771SStefano Zampini   (*submat)->ops->destroy               = MatDestroy_IS;
1534f26d0771SStefano Zampini   (*submat)->ops->setvalueslocal        = MatSetValuesLocal_SubMat_IS;
1535f26d0771SStefano Zampini   (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS;
1536f26d0771SStefano Zampini   (*submat)->ops->assemblybegin         = MatAssemblyBegin_IS;
1537f26d0771SStefano Zampini   (*submat)->ops->assemblyend           = MatAssemblyEnd_IS;
1538f26d0771SStefano Zampini   PetscFunctionReturn(0);
1539f26d0771SStefano Zampini }
1540f26d0771SStefano Zampini 
1541f26d0771SStefano Zampini #undef __FUNCT__
1542284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
1543284134d9SBarry Smith /*@
15443c212e90SHong Zhang     MatCreateIS - Creates a "process" unassembled matrix, assembled on each
1545284134d9SBarry Smith        process but not across processes.
1546284134d9SBarry Smith 
1547284134d9SBarry Smith    Input Parameters:
1548284134d9SBarry Smith +     comm    - MPI communicator that will share the matrix
1549e176bc59SStefano Zampini .     bs      - block size of the matrix
1550df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
1551e176bc59SStefano Zampini .     rmap    - local to global map for rows
1552e176bc59SStefano Zampini -     cmap    - local to global map for cols
1553284134d9SBarry Smith 
1554284134d9SBarry Smith    Output Parameter:
1555284134d9SBarry Smith .    A - the resulting matrix
1556284134d9SBarry Smith 
15578e6c10adSSatish Balay    Level: advanced
15588e6c10adSSatish Balay 
15593c212e90SHong Zhang    Notes: See MATIS for more details.
15603c212e90SHong Zhang           m and n are NOT related to the size of the map; they are the size of the part of the vector owned
1561e176bc59SStefano Zampini           by that process. The sizes of rmap and cmap define the size of the local matrices.
15623c212e90SHong Zhang           If either rmap or cmap are NULL, then the matrix is assumed to be square.
1563284134d9SBarry Smith 
1564284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
1565284134d9SBarry Smith @*/
1566e176bc59SStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A)
1567284134d9SBarry Smith {
1568284134d9SBarry Smith   PetscErrorCode ierr;
1569284134d9SBarry Smith 
1570284134d9SBarry Smith   PetscFunctionBegin;
1571e176bc59SStefano Zampini   if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping");
1572284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1573d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
1574284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1575284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
1576e176bc59SStefano Zampini   if (rmap && cmap) {
1577e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr);
1578e176bc59SStefano Zampini   } else if (!rmap) {
1579e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr);
1580e176bc59SStefano Zampini   } else {
1581e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr);
1582e176bc59SStefano Zampini   }
1583284134d9SBarry Smith   PetscFunctionReturn(0);
1584284134d9SBarry Smith }
1585284134d9SBarry Smith 
1586b4319ba4SBarry Smith /*MC
1587f26d0771SStefano Zampini    MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition methods (e.g. PCBDDC or KSPFETIDP).
1588b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
1589b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
1590b4319ba4SBarry Smith    product is handled "implicitly".
1591b4319ba4SBarry Smith 
1592b4319ba4SBarry Smith    Operations Provided:
15936726f965SBarry Smith +  MatMult()
15942e74eeadSLisandro Dalcin .  MatMultAdd()
15952e74eeadSLisandro Dalcin .  MatMultTranspose()
15962e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
15976726f965SBarry Smith .  MatZeroEntries()
15986726f965SBarry Smith .  MatSetOption()
15992e74eeadSLisandro Dalcin .  MatZeroRows()
16002e74eeadSLisandro Dalcin .  MatSetValues()
160197563a80SStefano Zampini .  MatSetValuesBlocked()
16026726f965SBarry Smith .  MatSetValuesLocal()
160397563a80SStefano Zampini .  MatSetValuesBlockedLocal()
16042e74eeadSLisandro Dalcin .  MatScale()
16052e74eeadSLisandro Dalcin .  MatGetDiagonal()
16062b404112SStefano Zampini .  MatMissingDiagonal()
16072b404112SStefano Zampini .  MatDuplicate()
16082b404112SStefano Zampini .  MatCopy()
1609f26d0771SStefano Zampini .  MatAXPY()
1610f26d0771SStefano Zampini .  MatGetSubMatrix()
1611f26d0771SStefano Zampini .  MatGetLocalSubMatrix()
16126726f965SBarry Smith -  MatSetLocalToGlobalMapping()
1613b4319ba4SBarry Smith 
1614b4319ba4SBarry Smith    Options Database Keys:
1615b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
1616b4319ba4SBarry Smith 
1617b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
1618b4319ba4SBarry Smith 
1619b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
1620b4319ba4SBarry Smith 
1621b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
1622eb82efa4SStefano Zampini           MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation()
1623b4319ba4SBarry Smith 
1624b4319ba4SBarry Smith   Level: advanced
1625b4319ba4SBarry Smith 
1626f26d0771SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), MatCreateIS(), PCBDDC, KSPFETIDP
1627b4319ba4SBarry Smith 
1628b4319ba4SBarry Smith M*/
1629b4319ba4SBarry Smith 
1630b4319ba4SBarry Smith #undef __FUNCT__
1631b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
16328cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
1633b4319ba4SBarry Smith {
1634dfbe8321SBarry Smith   PetscErrorCode ierr;
1635b4319ba4SBarry Smith   Mat_IS         *b;
1636b4319ba4SBarry Smith 
1637b4319ba4SBarry Smith   PetscFunctionBegin;
1638b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1639b4319ba4SBarry Smith   A->data = (void*)b;
1640b4319ba4SBarry Smith 
1641e176bc59SStefano Zampini   /* matrix ops */
1642e176bc59SStefano Zampini   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1643b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
16442e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
16452e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
16462e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1647b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1648b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
16492e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
165098921651SStefano Zampini   A->ops->setvaluesblocked        = MatSetValuesBlocked_IS;
1651b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1652f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
16532e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1654*f0ae7da4SStefano Zampini   A->ops->zerorowscolumns         = MatZeroRowsColumns_IS;
1655b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1656b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1657b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
16586726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
16592e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
16602e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
16616726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
166269796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
166369796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1664ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
16656bd84002SStefano Zampini   A->ops->missingdiagonal         = MatMissingDiagonal_IS;
16662b404112SStefano Zampini   A->ops->copy                    = MatCopy_IS;
1667659959c5SStefano Zampini   A->ops->getlocalsubmatrix       = MatGetLocalSubMatrix_IS;
1668a8116848SStefano Zampini   A->ops->getsubmatrix            = MatGetSubMatrix_IS;
1669f26d0771SStefano Zampini   A->ops->axpy                    = MatAXPY_IS;
16703fd1c9e7SStefano Zampini   A->ops->diagonalset             = MatDiagonalSet_IS;
16713fd1c9e7SStefano Zampini   A->ops->shift                   = MatShift_IS;
1672b4319ba4SBarry Smith 
1673b7ce53b6SStefano Zampini   /* special MATIS functions */
1674bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1675bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1676b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
16772e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
167817667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1679b4319ba4SBarry Smith   PetscFunctionReturn(0);
1680b4319ba4SBarry Smith }
1681