xref: /petsc/src/mat/impls/is/matis.c (revision a2ccb5f9ee22d16282353a7d1b8730ab22a81dc7)
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__
177fa8f2d3SStefano Zampini #define __FUNCT__ "MatGetInfo_IS"
187fa8f2d3SStefano Zampini static PetscErrorCode MatGetInfo_IS(Mat A,MatInfoType flag,MatInfo *ginfo)
197fa8f2d3SStefano Zampini {
207fa8f2d3SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
217fa8f2d3SStefano Zampini   MatInfo        info;
22314ce898Sstefano_zampini   PetscReal      isend[6],irecv[6];
237fa8f2d3SStefano Zampini   PetscInt       bs;
247fa8f2d3SStefano Zampini   PetscErrorCode ierr;
257fa8f2d3SStefano Zampini 
267fa8f2d3SStefano Zampini   PetscFunctionBegin;
277fa8f2d3SStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
28*a2ccb5f9Sstefano_zampini   if (matis->A->ops->getinfo) {
297fa8f2d3SStefano Zampini     ierr     = MatGetInfo(matis->A,MAT_LOCAL,&info);CHKERRQ(ierr);
307fa8f2d3SStefano Zampini     isend[0] = info.nz_used;
317fa8f2d3SStefano Zampini     isend[1] = info.nz_allocated;
327fa8f2d3SStefano Zampini     isend[2] = info.nz_unneeded;
337fa8f2d3SStefano Zampini     isend[3] = info.memory;
347fa8f2d3SStefano Zampini     isend[4] = info.mallocs;
35*a2ccb5f9Sstefano_zampini   } else {
36*a2ccb5f9Sstefano_zampini     isend[0] = 0.;
37*a2ccb5f9Sstefano_zampini     isend[1] = 0.;
38*a2ccb5f9Sstefano_zampini     isend[2] = 0.;
39*a2ccb5f9Sstefano_zampini     isend[3] = 0.;
40*a2ccb5f9Sstefano_zampini     isend[4] = 0.;
41*a2ccb5f9Sstefano_zampini   }
42314ce898Sstefano_zampini   isend[5] = matis->A->num_ass;
437fa8f2d3SStefano Zampini   if (flag == MAT_LOCAL) {
447fa8f2d3SStefano Zampini     ginfo->nz_used      = isend[0];
457fa8f2d3SStefano Zampini     ginfo->nz_allocated = isend[1];
467fa8f2d3SStefano Zampini     ginfo->nz_unneeded  = isend[2];
477fa8f2d3SStefano Zampini     ginfo->memory       = isend[3];
487fa8f2d3SStefano Zampini     ginfo->mallocs      = isend[4];
49314ce898Sstefano_zampini     ginfo->assemblies   = isend[5];
507fa8f2d3SStefano Zampini   } else if (flag == MAT_GLOBAL_MAX) {
51314ce898Sstefano_zampini     ierr = MPIU_Allreduce(isend,irecv,6,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
527fa8f2d3SStefano Zampini 
537fa8f2d3SStefano Zampini     ginfo->nz_used      = irecv[0];
547fa8f2d3SStefano Zampini     ginfo->nz_allocated = irecv[1];
557fa8f2d3SStefano Zampini     ginfo->nz_unneeded  = irecv[2];
567fa8f2d3SStefano Zampini     ginfo->memory       = irecv[3];
577fa8f2d3SStefano Zampini     ginfo->mallocs      = irecv[4];
58314ce898Sstefano_zampini     ginfo->assemblies   = irecv[5];
597fa8f2d3SStefano Zampini   } else if (flag == MAT_GLOBAL_SUM) {
607fa8f2d3SStefano Zampini     ierr = MPIU_Allreduce(isend,irecv,5,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
617fa8f2d3SStefano Zampini 
627fa8f2d3SStefano Zampini     ginfo->nz_used      = irecv[0];
637fa8f2d3SStefano Zampini     ginfo->nz_allocated = irecv[1];
647fa8f2d3SStefano Zampini     ginfo->nz_unneeded  = irecv[2];
657fa8f2d3SStefano Zampini     ginfo->memory       = irecv[3];
667fa8f2d3SStefano Zampini     ginfo->mallocs      = irecv[4];
677fa8f2d3SStefano Zampini     ginfo->assemblies   = A->num_ass;
687fa8f2d3SStefano Zampini   }
697fa8f2d3SStefano Zampini   ginfo->block_size        = bs;
707fa8f2d3SStefano Zampini   ginfo->fill_ratio_given  = 0;
717fa8f2d3SStefano Zampini   ginfo->fill_ratio_needed = 0;
727fa8f2d3SStefano Zampini   ginfo->factor_mallocs    = 0;
737fa8f2d3SStefano Zampini   PetscFunctionReturn(0);
747fa8f2d3SStefano Zampini }
757fa8f2d3SStefano Zampini 
767fa8f2d3SStefano Zampini #undef __FUNCT__
77d7f69cd0SStefano Zampini #define __FUNCT__ "MatTranspose_IS"
78d7f69cd0SStefano Zampini PetscErrorCode MatTranspose_IS(Mat A,MatReuse reuse,Mat *B)
79d7f69cd0SStefano Zampini {
80d7f69cd0SStefano Zampini   Mat                    C,lC,lA;
81d7f69cd0SStefano Zampini   ISLocalToGlobalMapping rl2g,cl2g;
82d7f69cd0SStefano Zampini   PetscBool              notset = PETSC_FALSE,cong = PETSC_TRUE;
83d7f69cd0SStefano Zampini   PetscErrorCode         ierr;
84d7f69cd0SStefano Zampini 
85d7f69cd0SStefano Zampini   PetscFunctionBegin;
86d7f69cd0SStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
87d7f69cd0SStefano Zampini     PetscBool rcong,ccong;
88d7f69cd0SStefano Zampini 
89d7f69cd0SStefano Zampini     ierr = PetscLayoutCompare((*B)->cmap,A->rmap,&rcong);CHKERRQ(ierr);
90d7f69cd0SStefano Zampini     ierr = PetscLayoutCompare((*B)->rmap,A->cmap,&ccong);CHKERRQ(ierr);
91fa520230SStefano Zampini     cong = (PetscBool)(rcong || ccong);
92d7f69cd0SStefano Zampini   }
93d7f69cd0SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX || *B == A || !cong) {
94d7f69cd0SStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
95d7f69cd0SStefano Zampini   } else {
96d7f69cd0SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATIS,&notset);CHKERRQ(ierr);
97d7f69cd0SStefano Zampini     C = *B;
98d7f69cd0SStefano Zampini   }
99d7f69cd0SStefano Zampini 
100d7f69cd0SStefano Zampini   if (!notset) {
101d7f69cd0SStefano Zampini     ierr = MatSetSizes(C,A->cmap->n,A->rmap->n,A->cmap->N,A->rmap->N);CHKERRQ(ierr);
102d7f69cd0SStefano Zampini     ierr = MatSetBlockSizes(C,PetscAbs(A->cmap->bs),PetscAbs(A->rmap->bs));CHKERRQ(ierr);
103d7f69cd0SStefano Zampini     ierr = MatSetType(C,MATIS);CHKERRQ(ierr);
104d7f69cd0SStefano Zampini     ierr = MatGetLocalToGlobalMapping(A,&rl2g,&cl2g);CHKERRQ(ierr);
105d7f69cd0SStefano Zampini     ierr = MatSetLocalToGlobalMapping(C,cl2g,rl2g);CHKERRQ(ierr);
106d7f69cd0SStefano Zampini   }
107d7f69cd0SStefano Zampini 
108d7f69cd0SStefano Zampini   /* perform local transposition */
109d7f69cd0SStefano Zampini   ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr);
110d7f69cd0SStefano Zampini   ierr = MatTranspose(lA,MAT_INITIAL_MATRIX,&lC);CHKERRQ(ierr);
111d7f69cd0SStefano Zampini   ierr = MatISSetLocalMat(C,lC);CHKERRQ(ierr);
112d7f69cd0SStefano Zampini   ierr = MatDestroy(&lC);CHKERRQ(ierr);
113d7f69cd0SStefano Zampini 
114d7f69cd0SStefano Zampini   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
115d7f69cd0SStefano Zampini   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
116d7f69cd0SStefano Zampini 
117d7f69cd0SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX || *B != A) {
118d7f69cd0SStefano Zampini     if (!cong) {
119d7f69cd0SStefano Zampini       ierr = MatDestroy(B);CHKERRQ(ierr);
120d7f69cd0SStefano Zampini     }
121d7f69cd0SStefano Zampini     *B = C;
122d7f69cd0SStefano Zampini   } else {
123d7f69cd0SStefano Zampini     ierr = MatHeaderMerge(A,&C);CHKERRQ(ierr);
124d7f69cd0SStefano Zampini   }
125d7f69cd0SStefano Zampini   PetscFunctionReturn(0);
126d7f69cd0SStefano Zampini }
127d7f69cd0SStefano Zampini 
128d7f69cd0SStefano Zampini #undef __FUNCT__
1293fd1c9e7SStefano Zampini #define __FUNCT__ "MatDiagonalSet_IS"
1303fd1c9e7SStefano Zampini PetscErrorCode  MatDiagonalSet_IS(Mat A,Vec D,InsertMode insmode)
1313fd1c9e7SStefano Zampini {
1323fd1c9e7SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
1333fd1c9e7SStefano Zampini   PetscErrorCode ierr;
1343fd1c9e7SStefano Zampini 
1353fd1c9e7SStefano Zampini   PetscFunctionBegin;
1363fd1c9e7SStefano Zampini   if (!D) { /* this code branch is used by MatShift_IS */
1373fd1c9e7SStefano Zampini     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
1383fd1c9e7SStefano Zampini   } else {
1393fd1c9e7SStefano Zampini     ierr = VecScatterBegin(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1403fd1c9e7SStefano Zampini     ierr = VecScatterEnd(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1413fd1c9e7SStefano Zampini   }
1423fd1c9e7SStefano Zampini   ierr = VecPointwiseDivide(is->y,is->y,is->counter);CHKERRQ(ierr);
1433fd1c9e7SStefano Zampini   ierr = MatDiagonalSet(is->A,is->y,insmode);CHKERRQ(ierr);
1443fd1c9e7SStefano Zampini   PetscFunctionReturn(0);
1453fd1c9e7SStefano Zampini }
1463fd1c9e7SStefano Zampini 
1473fd1c9e7SStefano Zampini #undef __FUNCT__
1483fd1c9e7SStefano Zampini #define __FUNCT__ "MatShift_IS"
1493fd1c9e7SStefano Zampini PetscErrorCode  MatShift_IS(Mat A,PetscScalar a)
1503fd1c9e7SStefano Zampini {
1513fd1c9e7SStefano Zampini   PetscErrorCode ierr;
1523fd1c9e7SStefano Zampini 
1533fd1c9e7SStefano Zampini   PetscFunctionBegin;
1543fd1c9e7SStefano Zampini   ierr = MatDiagonalSet_IS(A,NULL,ADD_VALUES);CHKERRQ(ierr);
1553fd1c9e7SStefano Zampini   PetscFunctionReturn(0);
1563fd1c9e7SStefano Zampini }
1573fd1c9e7SStefano Zampini 
1583fd1c9e7SStefano Zampini #undef __FUNCT__
159f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesLocal_SubMat_IS"
160f26d0771SStefano Zampini static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
161f26d0771SStefano Zampini {
162f26d0771SStefano Zampini   PetscErrorCode ierr;
163f26d0771SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
164f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
165f26d0771SStefano Zampini 
166f26d0771SStefano Zampini   PetscFunctionBegin;
167f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
168f26d0771SStefano 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);
169f26d0771SStefano Zampini #endif
170f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr);
171f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr);
172f26d0771SStefano Zampini   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
173f26d0771SStefano Zampini   PetscFunctionReturn(0);
174f26d0771SStefano Zampini }
175f26d0771SStefano Zampini 
176f26d0771SStefano Zampini #undef __FUNCT__
177f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesBlockedLocal_SubMat_IS"
178f26d0771SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
179f26d0771SStefano Zampini {
180f26d0771SStefano Zampini   PetscErrorCode ierr;
181f26d0771SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
182f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
183f26d0771SStefano Zampini 
184f26d0771SStefano Zampini   PetscFunctionBegin;
185f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
186f26d0771SStefano 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);
187f26d0771SStefano Zampini #endif
188f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr);
189f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr);
190f26d0771SStefano Zampini   ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
191f26d0771SStefano Zampini   PetscFunctionReturn(0);
192f26d0771SStefano Zampini }
193f26d0771SStefano Zampini 
194f26d0771SStefano Zampini #undef __FUNCT__
195a8116848SStefano Zampini #define __FUNCT__ "PetscLayoutMapLocal_Private"
196f0ae7da4SStefano Zampini static PetscErrorCode PetscLayoutMapLocal_Private(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs)
197a8116848SStefano Zampini {
198a8116848SStefano Zampini   PetscInt      *owners = map->range;
199a8116848SStefano Zampini   PetscInt       n      = map->n;
200a8116848SStefano Zampini   PetscSF        sf;
201a8116848SStefano Zampini   PetscInt      *lidxs,*work = NULL;
202a8116848SStefano Zampini   PetscSFNode   *ridxs;
203a8116848SStefano Zampini   PetscMPIInt    rank;
204a8116848SStefano Zampini   PetscInt       r, p = 0, len = 0;
205a8116848SStefano Zampini   PetscErrorCode ierr;
206a8116848SStefano Zampini 
207a8116848SStefano Zampini   PetscFunctionBegin;
208a8116848SStefano Zampini   /* Create SF where leaves are input idxs and roots are owned idxs (code adapted from MatZeroRowsMapLocal_Private) */
209a8116848SStefano Zampini   ierr = MPI_Comm_rank(map->comm,&rank);CHKERRQ(ierr);
210a8116848SStefano Zampini   ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr);
211a8116848SStefano Zampini   for (r = 0; r < n; ++r) lidxs[r] = -1;
212a8116848SStefano Zampini   ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr);
213a8116848SStefano Zampini   for (r = 0; r < N; ++r) {
214a8116848SStefano Zampini     const PetscInt idx = idxs[r];
215a8116848SStefano 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);
216a8116848SStefano Zampini     if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */
217a8116848SStefano Zampini       ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr);
218a8116848SStefano Zampini     }
219a8116848SStefano Zampini     ridxs[r].rank = p;
220a8116848SStefano Zampini     ridxs[r].index = idxs[r] - owners[p];
221a8116848SStefano Zampini   }
222a8116848SStefano Zampini   ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr);
223a8116848SStefano Zampini   ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr);
224a8116848SStefano Zampini   ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr);
225a8116848SStefano Zampini   ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr);
226f0ae7da4SStefano Zampini   if (ogidxs) { /* communicate global idxs */
227a8116848SStefano Zampini     PetscInt cum = 0,start,*work2;
228f0ae7da4SStefano Zampini 
229f0ae7da4SStefano Zampini     ierr = PetscMalloc1(n,&work);CHKERRQ(ierr);
230a8116848SStefano Zampini     ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr);
231a8116848SStefano Zampini     for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++;
232a8116848SStefano Zampini     ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRQ(ierr);
233a8116848SStefano Zampini     start -= cum;
234a8116848SStefano Zampini     cum = 0;
235a8116848SStefano Zampini     for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++;
236a8116848SStefano Zampini     ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr);
237a8116848SStefano Zampini     ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr);
238a8116848SStefano Zampini     ierr = PetscFree(work2);CHKERRQ(ierr);
239a8116848SStefano Zampini   }
240a8116848SStefano Zampini   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
241a8116848SStefano Zampini   /* Compress and put in indices */
242a8116848SStefano Zampini   for (r = 0; r < n; ++r)
243a8116848SStefano Zampini     if (lidxs[r] >= 0) {
244a8116848SStefano Zampini       if (work) work[len] = work[r];
245a8116848SStefano Zampini       lidxs[len++] = r;
246a8116848SStefano Zampini     }
247a8116848SStefano Zampini   if (on) *on = len;
248a8116848SStefano Zampini   if (oidxs) *oidxs = lidxs;
249a8116848SStefano Zampini   if (ogidxs) *ogidxs = work;
250a8116848SStefano Zampini   PetscFunctionReturn(0);
251a8116848SStefano Zampini }
252a8116848SStefano Zampini 
253a8116848SStefano Zampini #undef __FUNCT__
254a8116848SStefano Zampini #define __FUNCT__ "MatGetSubMatrix_IS"
255a8116848SStefano Zampini static PetscErrorCode MatGetSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat)
256a8116848SStefano Zampini {
257a8116848SStefano Zampini   Mat               locmat,newlocmat;
258a8116848SStefano Zampini   Mat_IS            *newmatis;
259a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG)
260a8116848SStefano Zampini   Vec               rtest,ltest;
261a8116848SStefano Zampini   const PetscScalar *array;
262a8116848SStefano Zampini #endif
263a8116848SStefano Zampini   const PetscInt    *idxs;
264a8116848SStefano Zampini   PetscInt          i,m,n;
265a8116848SStefano Zampini   PetscErrorCode    ierr;
266a8116848SStefano Zampini 
267a8116848SStefano Zampini   PetscFunctionBegin;
268a8116848SStefano Zampini   if (scall == MAT_REUSE_MATRIX) {
269a8116848SStefano Zampini     PetscBool ismatis;
270a8116848SStefano Zampini 
271a8116848SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr);
272a8116848SStefano Zampini     if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type");
273a8116848SStefano Zampini     newmatis = (Mat_IS*)(*newmat)->data;
274a8116848SStefano Zampini     if (!newmatis->getsub_ris) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS");
275a8116848SStefano Zampini     if (!newmatis->getsub_cis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS");
276a8116848SStefano Zampini   }
277a8116848SStefano Zampini   /* irow and icol may not have duplicate entries */
278a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG)
279a8116848SStefano Zampini   ierr = MatCreateVecs(mat,&ltest,&rtest);CHKERRQ(ierr);
280a8116848SStefano Zampini   ierr = ISGetLocalSize(irow,&n);CHKERRQ(ierr);
281a8116848SStefano Zampini   ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr);
282a8116848SStefano Zampini   for (i=0;i<n;i++) {
283a8116848SStefano Zampini     ierr = VecSetValue(rtest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr);
284a8116848SStefano Zampini   }
285a8116848SStefano Zampini   ierr = VecAssemblyBegin(rtest);CHKERRQ(ierr);
286a8116848SStefano Zampini   ierr = VecAssemblyEnd(rtest);CHKERRQ(ierr);
287a8116848SStefano Zampini   ierr = VecGetLocalSize(rtest,&n);CHKERRQ(ierr);
288a8116848SStefano Zampini   ierr = VecGetOwnershipRange(rtest,&m,NULL);CHKERRQ(ierr);
289a8116848SStefano Zampini   ierr = VecGetArrayRead(rtest,&array);CHKERRQ(ierr);
290fd479f66SMatthew 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]));
291a8116848SStefano Zampini   ierr = VecRestoreArrayRead(rtest,&array);CHKERRQ(ierr);
292a8116848SStefano Zampini   ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr);
293a8116848SStefano Zampini   ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr);
294a8116848SStefano Zampini   ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr);
295a8116848SStefano Zampini   for (i=0;i<n;i++) {
296a8116848SStefano Zampini     ierr = VecSetValue(ltest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr);
297a8116848SStefano Zampini   }
298a8116848SStefano Zampini   ierr = VecAssemblyBegin(ltest);CHKERRQ(ierr);
299a8116848SStefano Zampini   ierr = VecAssemblyEnd(ltest);CHKERRQ(ierr);
300a8116848SStefano Zampini   ierr = VecGetLocalSize(ltest,&n);CHKERRQ(ierr);
301a8116848SStefano Zampini   ierr = VecGetOwnershipRange(ltest,&m,NULL);CHKERRQ(ierr);
302a8116848SStefano Zampini   ierr = VecGetArrayRead(ltest,&array);CHKERRQ(ierr);
303fd479f66SMatthew 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]));
304a8116848SStefano Zampini   ierr = VecRestoreArrayRead(ltest,&array);CHKERRQ(ierr);
305a8116848SStefano Zampini   ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr);
306a8116848SStefano Zampini   ierr = VecDestroy(&rtest);CHKERRQ(ierr);
307a8116848SStefano Zampini   ierr = VecDestroy(&ltest);CHKERRQ(ierr);
308a8116848SStefano Zampini #endif
309a8116848SStefano Zampini   if (scall == MAT_INITIAL_MATRIX) {
310a8116848SStefano Zampini     Mat_IS                 *matis = (Mat_IS*)mat->data;
311a8116848SStefano Zampini     ISLocalToGlobalMapping rl2g;
312a8116848SStefano Zampini     IS                     is;
313a8116848SStefano Zampini     PetscInt               *lidxs,*lgidxs,*newgidxs;
3146dd40735SStefano Zampini     PetscInt               ll,newloc;
315a8116848SStefano Zampini     MPI_Comm               comm;
316a8116848SStefano Zampini 
317a8116848SStefano Zampini     ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
318a8116848SStefano Zampini     ierr = ISGetLocalSize(irow,&m);CHKERRQ(ierr);
319a8116848SStefano Zampini     ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr);
320a8116848SStefano Zampini     ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
321a8116848SStefano Zampini     ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr);
322a8116848SStefano Zampini     ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
323a8116848SStefano Zampini     /* communicate irow to their owners in the layout */
324a8116848SStefano Zampini     ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr);
325f0ae7da4SStefano Zampini     ierr = PetscLayoutMapLocal_Private(mat->rmap,m,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr);
326a8116848SStefano Zampini     ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr);
327a8116848SStefano Zampini     if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
328a8116848SStefano Zampini       ierr = MatISComputeSF_Private(mat);CHKERRQ(ierr);
329a8116848SStefano Zampini     }
330a8116848SStefano Zampini     ierr = PetscMemzero(matis->sf_rootdata,matis->sf_nroots*sizeof(PetscInt));CHKERRQ(ierr);
331a8116848SStefano Zampini     for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1;
332a8116848SStefano Zampini     ierr = PetscFree(lidxs);CHKERRQ(ierr);
333a8116848SStefano Zampini     ierr = PetscFree(lgidxs);CHKERRQ(ierr);
334a8116848SStefano Zampini     ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
335a8116848SStefano Zampini     ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
336a8116848SStefano Zampini     for (i=0,newloc=0;i<matis->sf_nleaves;i++) if (matis->sf_leafdata[i]) newloc++;
337a8116848SStefano Zampini     ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr);
338a8116848SStefano Zampini     ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr);
339a8116848SStefano Zampini     for (i=0,newloc=0;i<matis->sf_nleaves;i++)
340a8116848SStefano Zampini       if (matis->sf_leafdata[i]) {
341a8116848SStefano Zampini         lidxs[newloc] = i;
342a8116848SStefano Zampini         newgidxs[newloc++] = matis->sf_leafdata[i]-1;
343a8116848SStefano Zampini       }
344a8116848SStefano Zampini     ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
345a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
346a8116848SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
347a8116848SStefano Zampini     /* local is to extract local submatrix */
348a8116848SStefano Zampini     newmatis = (Mat_IS*)(*newmat)->data;
349a8116848SStefano Zampini     ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris);CHKERRQ(ierr);
350a8116848SStefano Zampini     if (mat->congruentlayouts == -1) { /* first time we compare rows and cols layouts */
351a8116848SStefano Zampini       PetscBool cong;
352a8116848SStefano Zampini       ierr = PetscLayoutCompare(mat->rmap,mat->cmap,&cong);CHKERRQ(ierr);
353a8116848SStefano Zampini       if (cong) mat->congruentlayouts = 1;
354a8116848SStefano Zampini       else      mat->congruentlayouts = 0;
355a8116848SStefano Zampini     }
356a8116848SStefano Zampini     if (mat->congruentlayouts && irow == icol) {
357a8116848SStefano Zampini       ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g);CHKERRQ(ierr);
358a8116848SStefano Zampini       ierr = PetscObjectReference((PetscObject)newmatis->getsub_ris);CHKERRQ(ierr);
359a8116848SStefano Zampini       newmatis->getsub_cis = newmatis->getsub_ris;
360a8116848SStefano Zampini     } else {
361a8116848SStefano Zampini       ISLocalToGlobalMapping cl2g;
362a8116848SStefano Zampini 
363a8116848SStefano Zampini       /* communicate icol to their owners in the layout */
364a8116848SStefano Zampini       ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr);
365f0ae7da4SStefano Zampini       ierr = PetscLayoutMapLocal_Private(mat->cmap,n,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr);
366a8116848SStefano Zampini       ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr);
367a8116848SStefano Zampini       ierr = PetscMemzero(matis->csf_rootdata,matis->csf_nroots*sizeof(PetscInt));CHKERRQ(ierr);
368a8116848SStefano Zampini       for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1;
369a8116848SStefano Zampini       ierr = PetscFree(lidxs);CHKERRQ(ierr);
370a8116848SStefano Zampini       ierr = PetscFree(lgidxs);CHKERRQ(ierr);
371a8116848SStefano Zampini       ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr);
372a8116848SStefano Zampini       ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr);
373a8116848SStefano Zampini       for (i=0,newloc=0;i<matis->csf_nleaves;i++) if (matis->csf_leafdata[i]) newloc++;
374a8116848SStefano Zampini       ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr);
375a8116848SStefano Zampini       ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr);
376a8116848SStefano Zampini       for (i=0,newloc=0;i<matis->csf_nleaves;i++)
377a8116848SStefano Zampini         if (matis->csf_leafdata[i]) {
378a8116848SStefano Zampini           lidxs[newloc] = i;
379a8116848SStefano Zampini           newgidxs[newloc++] = matis->csf_leafdata[i]-1;
380a8116848SStefano Zampini         }
381a8116848SStefano Zampini       ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
382a8116848SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
383a8116848SStefano Zampini       ierr = ISDestroy(&is);CHKERRQ(ierr);
384a8116848SStefano Zampini       /* local is to extract local submatrix */
385a8116848SStefano Zampini       ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis);CHKERRQ(ierr);
386a8116848SStefano Zampini       ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr);
387a8116848SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
388a8116848SStefano Zampini     }
389a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
390a8116848SStefano Zampini   } else {
391a8116848SStefano Zampini     ierr = MatISGetLocalMat(*newmat,&newlocmat);CHKERRQ(ierr);
392a8116848SStefano Zampini   }
393a8116848SStefano Zampini   ierr = MatISGetLocalMat(mat,&locmat);CHKERRQ(ierr);
394a8116848SStefano Zampini   newmatis = (Mat_IS*)(*newmat)->data;
395a8116848SStefano Zampini   ierr = MatGetSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat);CHKERRQ(ierr);
396a8116848SStefano Zampini   if (scall == MAT_INITIAL_MATRIX) {
397a8116848SStefano Zampini     ierr = MatISSetLocalMat(*newmat,newlocmat);CHKERRQ(ierr);
398a8116848SStefano Zampini     ierr = MatDestroy(&newlocmat);CHKERRQ(ierr);
399a8116848SStefano Zampini   }
400a8116848SStefano Zampini   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
401a8116848SStefano Zampini   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
402a8116848SStefano Zampini   PetscFunctionReturn(0);
403a8116848SStefano Zampini }
404a8116848SStefano Zampini 
405a8116848SStefano Zampini #undef __FUNCT__
4062b404112SStefano Zampini #define __FUNCT__ "MatCopy_IS"
407a8116848SStefano Zampini static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str)
4082b404112SStefano Zampini {
4092b404112SStefano Zampini   Mat_IS         *a = (Mat_IS*)A->data,*b;
4102b404112SStefano Zampini   PetscBool      ismatis;
4112b404112SStefano Zampini   PetscErrorCode ierr;
4122b404112SStefano Zampini 
4132b404112SStefano Zampini   PetscFunctionBegin;
4142b404112SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr);
4152b404112SStefano Zampini   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented");
4162b404112SStefano Zampini   b = (Mat_IS*)B->data;
4172b404112SStefano Zampini   ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr);
4182b404112SStefano Zampini   PetscFunctionReturn(0);
4192b404112SStefano Zampini }
4202b404112SStefano Zampini 
4212b404112SStefano Zampini #undef __FUNCT__
4226bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS"
423a8116848SStefano Zampini static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool  *missing,PetscInt *d)
4246bd84002SStefano Zampini {
425527b2640SStefano Zampini   Vec               v;
426527b2640SStefano Zampini   const PetscScalar *array;
427527b2640SStefano Zampini   PetscInt          i,n;
4286bd84002SStefano Zampini   PetscErrorCode    ierr;
4296bd84002SStefano Zampini 
4306bd84002SStefano Zampini   PetscFunctionBegin;
431527b2640SStefano Zampini   *missing = PETSC_FALSE;
432527b2640SStefano Zampini   ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr);
433527b2640SStefano Zampini   ierr = MatGetDiagonal(A,v);CHKERRQ(ierr);
434527b2640SStefano Zampini   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
435527b2640SStefano Zampini   ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr);
436527b2640SStefano Zampini   for (i=0;i<n;i++) if (array[i] == 0.) break;
437527b2640SStefano Zampini   ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr);
438527b2640SStefano Zampini   ierr = VecDestroy(&v);CHKERRQ(ierr);
439527b2640SStefano Zampini   if (i != n) *missing = PETSC_TRUE;
440527b2640SStefano Zampini   if (d) {
441527b2640SStefano Zampini     *d = -1;
442527b2640SStefano Zampini     if (*missing) {
443527b2640SStefano Zampini       PetscInt rstart;
444527b2640SStefano Zampini       ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
445527b2640SStefano Zampini       *d = i+rstart;
446527b2640SStefano Zampini     }
447527b2640SStefano Zampini   }
4486bd84002SStefano Zampini   PetscFunctionReturn(0);
4496bd84002SStefano Zampini }
4506bd84002SStefano Zampini 
4516bd84002SStefano Zampini #undef __FUNCT__
45228f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private"
453a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B)
45428f4e0baSStefano Zampini {
45528f4e0baSStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
45628f4e0baSStefano Zampini   const PetscInt *gidxs;
45728f4e0baSStefano Zampini   PetscErrorCode ierr;
45828f4e0baSStefano Zampini 
45928f4e0baSStefano Zampini   PetscFunctionBegin;
46028f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr);
46128f4e0baSStefano Zampini   ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr);
46228f4e0baSStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr);
4633bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
46428f4e0baSStefano Zampini   ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
4653bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
46628f4e0baSStefano Zampini   ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr);
467a8116848SStefano Zampini   if (B->rmap->mapping != B->cmap->mapping) { /* setup SF for columns */
468a8116848SStefano Zampini     ierr = MatGetSize(matis->A,NULL,&matis->csf_nleaves);CHKERRQ(ierr);
469a8116848SStefano Zampini     ierr = MatGetLocalSize(B,NULL,&matis->csf_nroots);CHKERRQ(ierr);
470a8116848SStefano Zampini     ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf);CHKERRQ(ierr);
471a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr);
472a8116848SStefano Zampini     ierr = PetscSFSetGraphLayout(matis->csf,B->cmap,matis->csf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
473a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr);
474a8116848SStefano Zampini     ierr = PetscMalloc2(matis->csf_nroots,&matis->csf_rootdata,matis->csf_nleaves,&matis->csf_leafdata);CHKERRQ(ierr);
475a8116848SStefano Zampini   } else {
476a8116848SStefano Zampini     matis->csf = matis->sf;
477a8116848SStefano Zampini     matis->csf_nleaves = matis->sf_nleaves;
478a8116848SStefano Zampini     matis->csf_nroots = matis->sf_nroots;
479a8116848SStefano Zampini     matis->csf_leafdata = matis->sf_leafdata;
480a8116848SStefano Zampini     matis->csf_rootdata = matis->sf_rootdata;
481a8116848SStefano Zampini   }
48228f4e0baSStefano Zampini   PetscFunctionReturn(0);
48328f4e0baSStefano Zampini }
4842e1947a5SStefano Zampini 
4852e1947a5SStefano Zampini #undef __FUNCT__
4862e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation"
487eb82efa4SStefano Zampini /*@
488a88811baSStefano Zampini    MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix.
489a88811baSStefano Zampini 
490a88811baSStefano Zampini    Collective on MPI_Comm
491a88811baSStefano Zampini 
492a88811baSStefano Zampini    Input Parameters:
493a88811baSStefano Zampini +  B - the matrix
494a88811baSStefano Zampini .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
495a88811baSStefano Zampini            (same value is used for all local rows)
496a88811baSStefano Zampini .  d_nnz - array containing the number of nonzeros in the various rows of the
497a88811baSStefano Zampini            DIAGONAL portion of the local submatrix (possibly different for each row)
498a88811baSStefano Zampini            or NULL, if d_nz is used to specify the nonzero structure.
499a88811baSStefano Zampini            The size of this array is equal to the number of local rows, i.e 'm'.
500a88811baSStefano Zampini            For matrices that will be factored, you must leave room for (and set)
501a88811baSStefano Zampini            the diagonal entry even if it is zero.
502a88811baSStefano Zampini .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
503a88811baSStefano Zampini            submatrix (same value is used for all local rows).
504a88811baSStefano Zampini -  o_nnz - array containing the number of nonzeros in the various rows of the
505a88811baSStefano Zampini            OFF-DIAGONAL portion of the local submatrix (possibly different for
506a88811baSStefano Zampini            each row) or NULL, if o_nz is used to specify the nonzero
507a88811baSStefano Zampini            structure. The size of this array is equal to the number
508a88811baSStefano Zampini            of local rows, i.e 'm'.
509a88811baSStefano Zampini 
510a88811baSStefano Zampini    If the *_nnz parameter is given then the *_nz parameter is ignored
511a88811baSStefano Zampini 
512a88811baSStefano Zampini    Level: intermediate
513a88811baSStefano Zampini 
514a88811baSStefano Zampini    Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition
515a88811baSStefano Zampini           from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local
516a88811baSStefano Zampini           matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects.
517a88811baSStefano Zampini 
518a88811baSStefano Zampini .keywords: matrix
519a88811baSStefano Zampini 
5203c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS
521a88811baSStefano Zampini @*/
5222e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
5232e1947a5SStefano Zampini {
5242e1947a5SStefano Zampini   PetscErrorCode ierr;
5252e1947a5SStefano Zampini 
5262e1947a5SStefano Zampini   PetscFunctionBegin;
5272e1947a5SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
5282e1947a5SStefano Zampini   PetscValidType(B,1);
5292e1947a5SStefano Zampini   ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
5302e1947a5SStefano Zampini   PetscFunctionReturn(0);
5312e1947a5SStefano Zampini }
5322e1947a5SStefano Zampini 
5332e1947a5SStefano Zampini #undef __FUNCT__
5342e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS"
5357230de76SStefano Zampini static PetscErrorCode  MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
5362e1947a5SStefano Zampini {
5372e1947a5SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
53828f4e0baSStefano Zampini   PetscInt       bs,i,nlocalcols;
5392e1947a5SStefano Zampini   PetscErrorCode ierr;
5402e1947a5SStefano Zampini 
5412e1947a5SStefano Zampini   PetscFunctionBegin;
5426c4ed002SBarry Smith   if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
54328f4e0baSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
54428f4e0baSStefano Zampini     ierr = MatISComputeSF_Private(B);CHKERRQ(ierr);
54528f4e0baSStefano Zampini   }
5462e1947a5SStefano Zampini   if (!d_nnz) {
54728f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz;
5482e1947a5SStefano Zampini   } else {
54928f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i];
5502e1947a5SStefano Zampini   }
5512e1947a5SStefano Zampini   if (!o_nnz) {
55228f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz;
5532e1947a5SStefano Zampini   } else {
55428f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i];
5552e1947a5SStefano Zampini   }
55628f4e0baSStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
55728f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr);
55828f4e0baSStefano Zampini   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
55928f4e0baSStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
56028f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves;i++) {
56128f4e0baSStefano Zampini     matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols);
5622e1947a5SStefano Zampini   }
56328f4e0baSStefano Zampini   ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr);
56428f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
56528f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs;
5662e1947a5SStefano Zampini   }
56728f4e0baSStefano Zampini   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
56828f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
56928f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i]-i;
5702e1947a5SStefano Zampini   }
57128f4e0baSStefano Zampini   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
5722e1947a5SStefano Zampini   PetscFunctionReturn(0);
5732e1947a5SStefano Zampini }
574b4319ba4SBarry Smith 
575b4319ba4SBarry Smith #undef __FUNCT__
5763927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private"
5773927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
5783927de2eSStefano Zampini {
5793927de2eSStefano Zampini   Mat_IS          *matis = (Mat_IS*)(A->data);
5803927de2eSStefano Zampini   PetscInt        *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership;
581ecf5a873SStefano Zampini   const PetscInt  *global_indices_r,*global_indices_c;
5823927de2eSStefano Zampini   PetscInt        i,j,bs,rows,cols;
5833927de2eSStefano Zampini   PetscInt        lrows,lcols;
5843927de2eSStefano Zampini   PetscInt        local_rows,local_cols;
5853927de2eSStefano Zampini   PetscMPIInt     nsubdomains;
5863927de2eSStefano Zampini   PetscBool       isdense,issbaij;
5873927de2eSStefano Zampini   PetscErrorCode  ierr;
5883927de2eSStefano Zampini 
5893927de2eSStefano Zampini   PetscFunctionBegin;
5903927de2eSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr);
5913927de2eSStefano Zampini   ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr);
5923927de2eSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
5933927de2eSStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
5943927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
5953927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
596ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
597ecf5a873SStefano Zampini   if (A->rmap->mapping != A->cmap->mapping) {
5987230de76SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
599ecf5a873SStefano Zampini   } else {
600ecf5a873SStefano Zampini     global_indices_c = global_indices_r;
601ecf5a873SStefano Zampini   }
602ecf5a873SStefano Zampini 
6033927de2eSStefano Zampini   if (issbaij) {
6043927de2eSStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
6053927de2eSStefano Zampini   }
6063927de2eSStefano Zampini   /*
607ecf5a873SStefano Zampini      An SF reduce is needed to sum up properly on shared rows.
6083927de2eSStefano Zampini      Note that generally preallocation is not exact, since it overestimates nonzeros
6093927de2eSStefano Zampini   */
6103927de2eSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
6113927de2eSStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
6123927de2eSStefano Zampini   }
6133927de2eSStefano Zampini   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
6143927de2eSStefano Zampini   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
6153927de2eSStefano Zampini   /* All processes need to compute entire row ownership */
6163927de2eSStefano Zampini   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
6173927de2eSStefano Zampini   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
6183927de2eSStefano Zampini   for (i=0;i<nsubdomains;i++) {
6193927de2eSStefano Zampini     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
6203927de2eSStefano Zampini       row_ownership[j] = i;
6213927de2eSStefano Zampini     }
6223927de2eSStefano Zampini   }
6237230de76SStefano Zampini   ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
6243927de2eSStefano Zampini 
6253927de2eSStefano Zampini   /*
6263927de2eSStefano Zampini      my_dnz and my_onz contains exact contribution to preallocation from each local mat
6273927de2eSStefano Zampini      then, they will be summed up properly. This way, preallocation is always sufficient
6283927de2eSStefano Zampini   */
6293927de2eSStefano Zampini   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
6303927de2eSStefano Zampini   /* preallocation as a MATAIJ */
6313927de2eSStefano Zampini   if (isdense) { /* special case for dense local matrices */
6323927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
633ecf5a873SStefano Zampini       PetscInt index_row = global_indices_r[i];
6343927de2eSStefano Zampini       for (j=i;j<local_rows;j++) {
6353927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
636ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[j];
6373927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
6383927de2eSStefano Zampini           my_dnz[i] += 1;
6393927de2eSStefano Zampini         } else { /* offdiag block */
6403927de2eSStefano Zampini           my_onz[i] += 1;
6413927de2eSStefano Zampini         }
6423927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
6433927de2eSStefano Zampini         if (i != j) {
6443927de2eSStefano Zampini           owner = row_ownership[index_col];
6453927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
6463927de2eSStefano Zampini             my_dnz[j] += 1;
6473927de2eSStefano Zampini           } else {
6483927de2eSStefano Zampini             my_onz[j] += 1;
6493927de2eSStefano Zampini           }
6503927de2eSStefano Zampini         }
6513927de2eSStefano Zampini       }
6523927de2eSStefano Zampini     }
653bb1015c3SStefano Zampini   } else if (matis->A->ops->getrowij) {
654bb1015c3SStefano Zampini     const PetscInt *ii,*jj,*jptr;
655bb1015c3SStefano Zampini     PetscBool      done;
656bb1015c3SStefano Zampini     ierr = MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr);
657bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__);
658bb1015c3SStefano Zampini     jptr = jj;
659bb1015c3SStefano Zampini     for (i=0;i<local_rows;i++) {
660bb1015c3SStefano Zampini       PetscInt index_row = global_indices_r[i];
661bb1015c3SStefano Zampini       for (j=0;j<ii[i+1]-ii[i];j++,jptr++) {
662bb1015c3SStefano Zampini         PetscInt owner = row_ownership[index_row];
663bb1015c3SStefano Zampini         PetscInt index_col = global_indices_c[*jptr];
664bb1015c3SStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
665bb1015c3SStefano Zampini           my_dnz[i] += 1;
666bb1015c3SStefano Zampini         } else { /* offdiag block */
667bb1015c3SStefano Zampini           my_onz[i] += 1;
668bb1015c3SStefano Zampini         }
669bb1015c3SStefano Zampini         /* same as before, interchanging rows and cols */
670bb1015c3SStefano Zampini         if (issbaij && index_col != index_row) {
671bb1015c3SStefano Zampini           owner = row_ownership[index_col];
672bb1015c3SStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
673bb1015c3SStefano Zampini             my_dnz[*jptr] += 1;
674bb1015c3SStefano Zampini           } else {
675bb1015c3SStefano Zampini             my_onz[*jptr] += 1;
676bb1015c3SStefano Zampini           }
677bb1015c3SStefano Zampini         }
678bb1015c3SStefano Zampini       }
679bb1015c3SStefano Zampini     }
680bb1015c3SStefano Zampini     ierr = MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr);
681bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__);
682bb1015c3SStefano Zampini   } else { /* loop over rows and use MatGetRow */
6833927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
6843927de2eSStefano Zampini       const PetscInt *cols;
685ecf5a873SStefano Zampini       PetscInt       ncols,index_row = global_indices_r[i];
6863927de2eSStefano Zampini       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
6873927de2eSStefano Zampini       for (j=0;j<ncols;j++) {
6883927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
689ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[cols[j]];
6903927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
6913927de2eSStefano Zampini           my_dnz[i] += 1;
6923927de2eSStefano Zampini         } else { /* offdiag block */
6933927de2eSStefano Zampini           my_onz[i] += 1;
6943927de2eSStefano Zampini         }
6953927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
696d9a9e74cSStefano Zampini         if (issbaij && index_col != index_row) {
6973927de2eSStefano Zampini           owner = row_ownership[index_col];
6983927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
699d9a9e74cSStefano Zampini             my_dnz[cols[j]] += 1;
7003927de2eSStefano Zampini           } else {
701d9a9e74cSStefano Zampini             my_onz[cols[j]] += 1;
7023927de2eSStefano Zampini           }
7033927de2eSStefano Zampini         }
7043927de2eSStefano Zampini       }
7053927de2eSStefano Zampini       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
7063927de2eSStefano Zampini     }
7073927de2eSStefano Zampini   }
708ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
709ecf5a873SStefano Zampini   if (global_indices_c != global_indices_r) {
7107230de76SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
711ecf5a873SStefano Zampini   }
7123927de2eSStefano Zampini   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
713ecf5a873SStefano Zampini 
714ecf5a873SStefano Zampini   /* Reduce my_dnz and my_onz */
7153927de2eSStefano Zampini   if (maxreduce) {
7163927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
7173927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
718bb1015c3SStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
7193927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
7203927de2eSStefano Zampini   } else {
7213927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
7223927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
723bb1015c3SStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
7243927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
7253927de2eSStefano Zampini   }
7263927de2eSStefano Zampini   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
7273927de2eSStefano Zampini 
7283927de2eSStefano Zampini   /* Resize preallocation if overestimated */
7293927de2eSStefano Zampini   for (i=0;i<lrows;i++) {
7303927de2eSStefano Zampini     dnz[i] = PetscMin(dnz[i],lcols);
7313927de2eSStefano Zampini     onz[i] = PetscMin(onz[i],cols-lcols);
7323927de2eSStefano Zampini   }
7333927de2eSStefano Zampini   /* set preallocation */
7343927de2eSStefano Zampini   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
7353927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
7363927de2eSStefano Zampini     dnz[i] = dnz[i*bs]/bs;
7373927de2eSStefano Zampini     onz[i] = onz[i*bs]/bs;
7383927de2eSStefano Zampini   }
7393927de2eSStefano Zampini   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
7403927de2eSStefano Zampini   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
7413927de2eSStefano Zampini   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
7423927de2eSStefano Zampini   if (issbaij) {
7433927de2eSStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
7443927de2eSStefano Zampini   }
7453927de2eSStefano Zampini   PetscFunctionReturn(0);
7463927de2eSStefano Zampini }
7473927de2eSStefano Zampini 
7483927de2eSStefano Zampini #undef __FUNCT__
749b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
7507230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
751b7ce53b6SStefano Zampini {
752b7ce53b6SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
753d9a9e74cSStefano Zampini   Mat            local_mat;
754b7ce53b6SStefano Zampini   /* info on mat */
7553cfa4ea4SStefano Zampini   PetscInt       bs,rows,cols,lrows,lcols;
756b7ce53b6SStefano Zampini   PetscInt       local_rows,local_cols;
757686e3a49SStefano Zampini   PetscBool      isdense,issbaij,isseqaij;
7587c03b4e8SStefano Zampini   PetscMPIInt    nsubdomains;
759b7ce53b6SStefano Zampini   /* values insertion */
760b7ce53b6SStefano Zampini   PetscScalar    *array;
761b7ce53b6SStefano Zampini   /* work */
762b7ce53b6SStefano Zampini   PetscErrorCode ierr;
763b7ce53b6SStefano Zampini 
764b7ce53b6SStefano Zampini   PetscFunctionBegin;
765b7ce53b6SStefano Zampini   /* get info from mat */
7667c03b4e8SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
7677c03b4e8SStefano Zampini   if (nsubdomains == 1) {
7687c03b4e8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
7697c03b4e8SStefano Zampini       ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr);
7707c03b4e8SStefano Zampini     } else {
7717c03b4e8SStefano Zampini       ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
7727c03b4e8SStefano Zampini     }
7737c03b4e8SStefano Zampini     PetscFunctionReturn(0);
7747c03b4e8SStefano Zampini   }
775b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
776b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
7773cfa4ea4SStefano Zampini   ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr);
778b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
779b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
780686e3a49SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
781b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
782b7ce53b6SStefano Zampini 
783b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
784b7ce53b6SStefano Zampini     MatType     new_mat_type;
7853927de2eSStefano Zampini     PetscBool   issbaij_red;
786b7ce53b6SStefano Zampini 
787b7ce53b6SStefano Zampini     /* determining new matrix type */
788b2566f29SBarry Smith     ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
789b7ce53b6SStefano Zampini     if (issbaij_red) {
790b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
791b7ce53b6SStefano Zampini     } else {
792b7ce53b6SStefano Zampini       if (bs>1) {
793b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
794b7ce53b6SStefano Zampini       } else {
795b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
796b7ce53b6SStefano Zampini       }
797b7ce53b6SStefano Zampini     }
798b7ce53b6SStefano Zampini 
7993927de2eSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
8003cfa4ea4SStefano Zampini     ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr);
8013927de2eSStefano Zampini     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
8023927de2eSStefano Zampini     ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr);
8033927de2eSStefano Zampini     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
804b7ce53b6SStefano Zampini   } else {
8053cfa4ea4SStefano Zampini     PetscInt mbs,mrows,mcols,mlrows,mlcols;
806b7ce53b6SStefano Zampini     /* some checks */
807b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
808b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
8093cfa4ea4SStefano Zampini     ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr);
8106c4ed002SBarry Smith     if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
8116c4ed002SBarry Smith     if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
8126c4ed002SBarry Smith     if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows);
8136c4ed002SBarry Smith     if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols);
8146c4ed002SBarry Smith     if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
815b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
816b7ce53b6SStefano Zampini   }
817d9a9e74cSStefano Zampini 
818d9a9e74cSStefano Zampini   if (issbaij) {
819d9a9e74cSStefano Zampini     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr);
820d9a9e74cSStefano Zampini   } else {
821d9a9e74cSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
822d9a9e74cSStefano Zampini     local_mat = matis->A;
823d9a9e74cSStefano Zampini   }
824686e3a49SStefano Zampini 
825b7ce53b6SStefano Zampini   /* Set values */
826ecf5a873SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr);
827b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
828ecf5a873SStefano Zampini     PetscInt i,*dummy_rows,*dummy_cols;
829ecf5a873SStefano Zampini 
830ecf5a873SStefano Zampini     if (local_rows != local_cols) {
831ecf5a873SStefano Zampini       ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr);
832ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
833ecf5a873SStefano Zampini       for (i=0;i<local_cols;i++) dummy_cols[i] = i;
834ecf5a873SStefano Zampini     } else {
835ecf5a873SStefano Zampini       ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr);
836ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
837ecf5a873SStefano Zampini       dummy_cols = dummy_rows;
838ecf5a873SStefano Zampini     }
839b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
840d9a9e74cSStefano Zampini     ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr);
841ecf5a873SStefano Zampini     ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr);
842d9a9e74cSStefano Zampini     ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr);
843ecf5a873SStefano Zampini     if (dummy_rows != dummy_cols) {
844ecf5a873SStefano Zampini       ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr);
845ecf5a873SStefano Zampini     } else {
846ecf5a873SStefano Zampini       ierr = PetscFree(dummy_rows);CHKERRQ(ierr);
847ecf5a873SStefano Zampini     }
848686e3a49SStefano Zampini   } else if (isseqaij) {
849ecf5a873SStefano Zampini     PetscInt  i,nvtxs,*xadj,*adjncy;
850686e3a49SStefano Zampini     PetscBool done;
851686e3a49SStefano Zampini 
852d9a9e74cSStefano Zampini     ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
853bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__);
854d9a9e74cSStefano Zampini     ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr);
855686e3a49SStefano Zampini     for (i=0;i<nvtxs;i++) {
856ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr);
857686e3a49SStefano Zampini     }
858d9a9e74cSStefano Zampini     ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
859bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__);
860d9a9e74cSStefano Zampini     ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr);
861686e3a49SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
862ecf5a873SStefano Zampini     PetscInt i;
863c0962df8SStefano Zampini 
864686e3a49SStefano Zampini     for (i=0;i<local_rows;i++) {
865686e3a49SStefano Zampini       PetscInt       j;
866ecf5a873SStefano Zampini       const PetscInt *local_indices_cols;
867686e3a49SStefano Zampini 
868ecf5a873SStefano Zampini       ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
869ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr);
870ecf5a873SStefano Zampini       ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
871686e3a49SStefano Zampini     }
872b7ce53b6SStefano Zampini   }
873b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
874d9a9e74cSStefano Zampini   ierr = MatDestroy(&local_mat);CHKERRQ(ierr);
875b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
876b7ce53b6SStefano Zampini   if (isdense) {
877b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
878b7ce53b6SStefano Zampini   }
879b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
880b7ce53b6SStefano Zampini }
881b7ce53b6SStefano Zampini 
882b7ce53b6SStefano Zampini #undef __FUNCT__
883b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
884b7ce53b6SStefano Zampini /*@
885b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
886b7ce53b6SStefano Zampini 
887b7ce53b6SStefano Zampini   Input Parameter:
888b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
889b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
890b7ce53b6SStefano Zampini 
891b7ce53b6SStefano Zampini   Output Parameter:
892b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
893b7ce53b6SStefano Zampini 
894b7ce53b6SStefano Zampini   Level: developer
895b7ce53b6SStefano Zampini 
896eb82efa4SStefano Zampini   Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested.
897b7ce53b6SStefano Zampini 
898b7ce53b6SStefano Zampini .seealso: MATIS
899b7ce53b6SStefano Zampini @*/
900b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
901b7ce53b6SStefano Zampini {
902b7ce53b6SStefano Zampini   PetscErrorCode ierr;
903b7ce53b6SStefano Zampini 
904b7ce53b6SStefano Zampini   PetscFunctionBegin;
905b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
906b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
907b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
908b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
909b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
910b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
9116c4ed002SBarry Smith     if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix");
912b7ce53b6SStefano Zampini   }
913b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
914b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
915b7ce53b6SStefano Zampini }
916b7ce53b6SStefano Zampini 
917b7ce53b6SStefano Zampini #undef __FUNCT__
918ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
919ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
920ad6194a2SStefano Zampini {
921ad6194a2SStefano Zampini   PetscErrorCode ierr;
922ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
923ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
924ad6194a2SStefano Zampini   Mat            B,localmat;
925ad6194a2SStefano Zampini 
926ad6194a2SStefano Zampini   PetscFunctionBegin;
927ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
928ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
929ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
930e176bc59SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr);
931ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
932ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
933b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
934ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
935ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
936ad6194a2SStefano Zampini   *newmat = B;
937ad6194a2SStefano Zampini   PetscFunctionReturn(0);
938ad6194a2SStefano Zampini }
939ad6194a2SStefano Zampini 
940ad6194a2SStefano Zampini #undef __FUNCT__
94169796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
942a8116848SStefano Zampini static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
94369796d55SStefano Zampini {
94469796d55SStefano Zampini   PetscErrorCode ierr;
94569796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
94669796d55SStefano Zampini   PetscBool      local_sym;
94769796d55SStefano Zampini 
94869796d55SStefano Zampini   PetscFunctionBegin;
94969796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
950b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
95169796d55SStefano Zampini   PetscFunctionReturn(0);
95269796d55SStefano Zampini }
95369796d55SStefano Zampini 
95469796d55SStefano Zampini #undef __FUNCT__
95569796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
956a8116848SStefano Zampini static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
95769796d55SStefano Zampini {
95869796d55SStefano Zampini   PetscErrorCode ierr;
95969796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
96069796d55SStefano Zampini   PetscBool      local_sym;
96169796d55SStefano Zampini 
96269796d55SStefano Zampini   PetscFunctionBegin;
96369796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
964b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
96569796d55SStefano Zampini   PetscFunctionReturn(0);
96669796d55SStefano Zampini }
96769796d55SStefano Zampini 
96869796d55SStefano Zampini #undef __FUNCT__
969b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
970a8116848SStefano Zampini static PetscErrorCode MatDestroy_IS(Mat A)
971b4319ba4SBarry Smith {
972dfbe8321SBarry Smith   PetscErrorCode ierr;
973b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
974b4319ba4SBarry Smith 
975b4319ba4SBarry Smith   PetscFunctionBegin;
9766bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
977e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr);
978e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr);
9796bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
9806bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
9813fd1c9e7SStefano Zampini   ierr = VecDestroy(&b->counter);CHKERRQ(ierr);
982a8116848SStefano Zampini   ierr = ISDestroy(&b->getsub_ris);CHKERRQ(ierr);
983a8116848SStefano Zampini   ierr = ISDestroy(&b->getsub_cis);CHKERRQ(ierr);
984a8116848SStefano Zampini   if (b->sf != b->csf) {
985a8116848SStefano Zampini     ierr = PetscSFDestroy(&b->csf);CHKERRQ(ierr);
986a8116848SStefano Zampini     ierr = PetscFree2(b->csf_rootdata,b->csf_leafdata);CHKERRQ(ierr);
987a8116848SStefano Zampini   }
98828f4e0baSStefano Zampini   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
98928f4e0baSStefano Zampini   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
990bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
991dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
992bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
993b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
994b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
9952e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
996b4319ba4SBarry Smith   PetscFunctionReturn(0);
997b4319ba4SBarry Smith }
998b4319ba4SBarry Smith 
999b4319ba4SBarry Smith #undef __FUNCT__
1000b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
1001a8116848SStefano Zampini static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
1002b4319ba4SBarry Smith {
1003dfbe8321SBarry Smith   PetscErrorCode ierr;
1004b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
1005b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
1006b4319ba4SBarry Smith 
1007b4319ba4SBarry Smith   PetscFunctionBegin;
1008b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
1009e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1010e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1011b4319ba4SBarry Smith 
1012b4319ba4SBarry Smith   /* multiply the local matrix */
1013b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
1014b4319ba4SBarry Smith 
1015b4319ba4SBarry Smith   /* scatter product back into global memory */
10162dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
1017e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1018e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1019b4319ba4SBarry Smith   PetscFunctionReturn(0);
1020b4319ba4SBarry Smith }
1021b4319ba4SBarry Smith 
1022b4319ba4SBarry Smith #undef __FUNCT__
10232e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
1024a8116848SStefano Zampini static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
10252e74eeadSLisandro Dalcin {
1026650997f4SStefano Zampini   Vec            temp_vec;
10272e74eeadSLisandro Dalcin   PetscErrorCode ierr;
10282e74eeadSLisandro Dalcin 
10292e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
1030650997f4SStefano Zampini   if (v3 != v2) {
1031650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
1032650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
1033650997f4SStefano Zampini   } else {
1034650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
1035650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
1036650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
1037650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
1038650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
1039650997f4SStefano Zampini   }
10402e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
10412e74eeadSLisandro Dalcin }
10422e74eeadSLisandro Dalcin 
10432e74eeadSLisandro Dalcin #undef __FUNCT__
10442e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
1045a8116848SStefano Zampini static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x)
10462e74eeadSLisandro Dalcin {
10472e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
10482e74eeadSLisandro Dalcin   PetscErrorCode ierr;
10492e74eeadSLisandro Dalcin 
1050e176bc59SStefano Zampini   PetscFunctionBegin;
10512e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
1052e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1053e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
10542e74eeadSLisandro Dalcin 
10552e74eeadSLisandro Dalcin   /* multiply the local matrix */
1056e176bc59SStefano Zampini   ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr);
10572e74eeadSLisandro Dalcin 
10582e74eeadSLisandro Dalcin   /* scatter product back into global vector */
1059e176bc59SStefano Zampini   ierr = VecSet(x,0);CHKERRQ(ierr);
1060e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1061e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10622e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
10632e74eeadSLisandro Dalcin }
10642e74eeadSLisandro Dalcin 
10652e74eeadSLisandro Dalcin #undef __FUNCT__
10662e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
1067a8116848SStefano Zampini static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
10682e74eeadSLisandro Dalcin {
1069650997f4SStefano Zampini   Vec            temp_vec;
10702e74eeadSLisandro Dalcin   PetscErrorCode ierr;
10712e74eeadSLisandro Dalcin 
10722e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
1073650997f4SStefano Zampini   if (v3 != v2) {
1074650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
1075650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
1076650997f4SStefano Zampini   } else {
1077650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
1078650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
1079650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
1080650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
1081650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
1082650997f4SStefano Zampini   }
10832e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
10842e74eeadSLisandro Dalcin }
10852e74eeadSLisandro Dalcin 
10862e74eeadSLisandro Dalcin #undef __FUNCT__
1087b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
1088a8116848SStefano Zampini static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
1089b4319ba4SBarry Smith {
1090b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
1091dfbe8321SBarry Smith   PetscErrorCode ierr;
1092b4319ba4SBarry Smith   PetscViewer    sviewer;
1093b4319ba4SBarry Smith 
1094b4319ba4SBarry Smith   PetscFunctionBegin;
10953f08860eSBarry Smith   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
1096b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
10973f08860eSBarry Smith   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
10986e520ac8SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1099b4319ba4SBarry Smith   PetscFunctionReturn(0);
1100b4319ba4SBarry Smith }
1101b4319ba4SBarry Smith 
1102b4319ba4SBarry Smith #undef __FUNCT__
1103b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
1104a8116848SStefano Zampini static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
1105b4319ba4SBarry Smith {
1106dfbe8321SBarry Smith   PetscErrorCode ierr;
1107e176bc59SStefano Zampini   PetscInt       nr,rbs,nc,cbs;
1108b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1109b4319ba4SBarry Smith   IS             from,to;
1110e176bc59SStefano Zampini   Vec            cglobal,rglobal;
1111b4319ba4SBarry Smith 
1112b4319ba4SBarry Smith   PetscFunctionBegin;
1113784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
1114e176bc59SStefano Zampini   PetscCheckSameComm(A,1,cmapping,3);
11153bbff08aSStefano Zampini   /* Destroy any previous data */
111670cf5478SStefano Zampini   ierr = VecDestroy(&is->x);CHKERRQ(ierr);
111770cf5478SStefano Zampini   ierr = VecDestroy(&is->y);CHKERRQ(ierr);
11183fd1c9e7SStefano Zampini   ierr = VecDestroy(&is->counter);CHKERRQ(ierr);
1119e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr);
1120e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr);
11211c47cb0fSStefano Zampini   ierr = MatDestroy(&is->A);CHKERRQ(ierr);
112228f4e0baSStefano Zampini   ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
112328f4e0baSStefano Zampini   ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
11243bbff08aSStefano Zampini 
11253bbff08aSStefano Zampini   /* Setup Layout and set local to global maps */
1126fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1127fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1128fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
1129fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
1130b4319ba4SBarry Smith 
1131b4319ba4SBarry Smith   /* Create the local matrix A */
1132e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr);
1133e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr);
1134e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr);
1135e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr);
1136f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
1137e176bc59SStefano Zampini   ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr);
1138e176bc59SStefano Zampini   ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr);
1139e176bc59SStefano Zampini   ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr);
1140ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
1141ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
1142b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
1143b4319ba4SBarry Smith 
1144f26d0771SStefano Zampini   if (!is->islocalref) { /* setup scatters and local vectors for MatMult */
1145b4319ba4SBarry Smith     /* Create the local work vectors */
11463bbff08aSStefano Zampini     ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr);
1147b4319ba4SBarry Smith 
1148e176bc59SStefano Zampini     /* setup the global to local scatters */
1149e176bc59SStefano Zampini     ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr);
1150e176bc59SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr);
1151784ac674SJed Brown     ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
1152e176bc59SStefano Zampini     ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr);
1153e176bc59SStefano Zampini     if (rmapping != cmapping) {
1154e176bc59SStefano Zampini       ierr = ISDestroy(&to);CHKERRQ(ierr);
1155e176bc59SStefano Zampini       ierr = ISDestroy(&from);CHKERRQ(ierr);
1156e176bc59SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr);
1157e176bc59SStefano Zampini       ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr);
1158e176bc59SStefano Zampini       ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr);
1159e176bc59SStefano Zampini     } else {
1160e176bc59SStefano Zampini       ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr);
1161e176bc59SStefano Zampini       is->cctx = is->rctx;
1162e176bc59SStefano Zampini     }
11633fd1c9e7SStefano Zampini 
11643fd1c9e7SStefano Zampini     /* interface counter vector (local) */
11653fd1c9e7SStefano Zampini     ierr = VecDuplicate(is->y,&is->counter);CHKERRQ(ierr);
11663fd1c9e7SStefano Zampini     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
11673fd1c9e7SStefano Zampini     ierr = VecScatterBegin(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11683fd1c9e7SStefano Zampini     ierr = VecScatterEnd(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
11693fd1c9e7SStefano Zampini     ierr = VecScatterBegin(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11703fd1c9e7SStefano Zampini     ierr = VecScatterEnd(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11713fd1c9e7SStefano Zampini 
11723fd1c9e7SStefano Zampini     /* free workspace */
1173e176bc59SStefano Zampini     ierr = VecDestroy(&rglobal);CHKERRQ(ierr);
1174e176bc59SStefano Zampini     ierr = VecDestroy(&cglobal);CHKERRQ(ierr);
11756bf464f9SBarry Smith     ierr = ISDestroy(&to);CHKERRQ(ierr);
11766bf464f9SBarry Smith     ierr = ISDestroy(&from);CHKERRQ(ierr);
1177f26d0771SStefano Zampini   }
117848ff6bf3SStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
1179b4319ba4SBarry Smith   PetscFunctionReturn(0);
1180b4319ba4SBarry Smith }
1181b4319ba4SBarry Smith 
11822e74eeadSLisandro Dalcin #undef __FUNCT__
11832e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
1184a8116848SStefano Zampini static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
11852e74eeadSLisandro Dalcin {
11862e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
11872e74eeadSLisandro Dalcin   PetscErrorCode ierr;
118897563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
118997563a80SStefano Zampini   PetscInt       i,zm,zn;
119097563a80SStefano Zampini #endif
1191f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
11922e74eeadSLisandro Dalcin 
11932e74eeadSLisandro Dalcin   PetscFunctionBegin;
11942e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
1195f26d0771SStefano 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);
119697563a80SStefano Zampini   /* count negative indices */
119797563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
119897563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
11992e74eeadSLisandro Dalcin #endif
120097563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
120197563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
120297563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
120397563a80SStefano Zampini   /* count negative indices (should be the same as before) */
120497563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
120597563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
120697563a80SStefano 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");
120797563a80SStefano 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");
120897563a80SStefano Zampini #endif
12092e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
12102e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
12112e74eeadSLisandro Dalcin }
12122e74eeadSLisandro Dalcin 
1213b4319ba4SBarry Smith #undef __FUNCT__
121497563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS"
1215a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
121697563a80SStefano Zampini {
121797563a80SStefano Zampini   Mat_IS         *is = (Mat_IS*)mat->data;
121897563a80SStefano Zampini   PetscErrorCode ierr;
121997563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
122097563a80SStefano Zampini   PetscInt       i,zm,zn;
122197563a80SStefano Zampini #endif
1222f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
122397563a80SStefano Zampini 
122497563a80SStefano Zampini   PetscFunctionBegin;
122597563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
1226f26d0771SStefano 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);
122797563a80SStefano Zampini   /* count negative indices */
122897563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
122997563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
123097563a80SStefano Zampini #endif
123197563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
123297563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
123397563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
123497563a80SStefano Zampini   /* count negative indices (should be the same as before) */
123597563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
123697563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
123797563a80SStefano 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");
123897563a80SStefano 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");
123997563a80SStefano Zampini #endif
124097563a80SStefano Zampini   ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
124197563a80SStefano Zampini   PetscFunctionReturn(0);
124297563a80SStefano Zampini }
124397563a80SStefano Zampini 
124497563a80SStefano Zampini #undef __FUNCT__
1245b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
1246a8116848SStefano Zampini static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
1247b4319ba4SBarry Smith {
1248dfbe8321SBarry Smith   PetscErrorCode ierr;
1249b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1250b4319ba4SBarry Smith 
1251b4319ba4SBarry Smith   PetscFunctionBegin;
1252b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
1253b4319ba4SBarry Smith   PetscFunctionReturn(0);
1254b4319ba4SBarry Smith }
1255b4319ba4SBarry Smith 
1256b4319ba4SBarry Smith #undef __FUNCT__
1257f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
1258a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
1259f0006bf2SLisandro Dalcin {
1260f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
1261f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
1262f0006bf2SLisandro Dalcin 
1263f0006bf2SLisandro Dalcin   PetscFunctionBegin;
1264f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
1265f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
1266f0006bf2SLisandro Dalcin }
1267f0006bf2SLisandro Dalcin 
1268f0006bf2SLisandro Dalcin #undef __FUNCT__
1269f0ae7da4SStefano Zampini #define __FUNCT__ "MatISZeroRowsColumnsLocal_Private"
1270f0ae7da4SStefano Zampini static PetscErrorCode MatISZeroRowsColumnsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,PetscBool columns)
1271f0ae7da4SStefano Zampini {
1272f0ae7da4SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
1273f0ae7da4SStefano Zampini   PetscErrorCode ierr;
1274f0ae7da4SStefano Zampini 
1275f0ae7da4SStefano Zampini   PetscFunctionBegin;
1276f0ae7da4SStefano Zampini   if (!n) {
1277f0ae7da4SStefano Zampini     is->pure_neumann = PETSC_TRUE;
1278f0ae7da4SStefano Zampini   } else {
1279f0ae7da4SStefano Zampini     PetscInt i;
1280f0ae7da4SStefano Zampini     is->pure_neumann = PETSC_FALSE;
1281f0ae7da4SStefano Zampini 
1282f0ae7da4SStefano Zampini     if (columns) {
1283f0ae7da4SStefano Zampini       ierr = MatZeroRowsColumns(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
1284f0ae7da4SStefano Zampini     } else {
1285f0ae7da4SStefano Zampini       ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
1286f0ae7da4SStefano Zampini     }
1287f0ae7da4SStefano Zampini     if (diag != 0.) {
1288f0ae7da4SStefano Zampini       const PetscScalar *array;
1289f0ae7da4SStefano Zampini       ierr = VecGetArrayRead(is->counter,&array);CHKERRQ(ierr);
1290f0ae7da4SStefano Zampini       for (i=0; i<n; i++) {
1291f0ae7da4SStefano Zampini         ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
1292f0ae7da4SStefano Zampini       }
1293f0ae7da4SStefano Zampini       ierr = VecRestoreArrayRead(is->counter,&array);CHKERRQ(ierr);
1294f0ae7da4SStefano Zampini     }
1295f0ae7da4SStefano Zampini     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1296f0ae7da4SStefano Zampini     ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1297f0ae7da4SStefano Zampini   }
1298f0ae7da4SStefano Zampini   PetscFunctionReturn(0);
1299f0ae7da4SStefano Zampini }
1300f0ae7da4SStefano Zampini 
1301f0ae7da4SStefano Zampini #undef __FUNCT__
1302f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_Private_IS"
1303f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Private_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b,PetscBool columns)
13042e74eeadSLisandro Dalcin {
13056e520ac8SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
13066e520ac8SStefano Zampini   PetscInt       nr,nl,len,i;
13076e520ac8SStefano Zampini   PetscInt       *lrows;
13082e74eeadSLisandro Dalcin   PetscErrorCode ierr;
13092e74eeadSLisandro Dalcin 
13102e74eeadSLisandro Dalcin   PetscFunctionBegin;
1311f0ae7da4SStefano Zampini #if defined(PETSC_USE_DEBUG)
1312f0ae7da4SStefano Zampini   if (columns || diag != 0. || (x && b)) {
1313f0ae7da4SStefano Zampini     PetscBool cong;
1314f0ae7da4SStefano Zampini     ierr = PetscLayoutCompare(A->rmap,A->cmap,&cong);CHKERRQ(ierr);
1315f0ae7da4SStefano 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");
1316f0ae7da4SStefano 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");
1317f0ae7da4SStefano Zampini     if (!cong && x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"A->rmap and A->cmap need to be congruent");
1318f0ae7da4SStefano Zampini   }
1319f0ae7da4SStefano Zampini #endif
13206e520ac8SStefano Zampini   /* get locally owned rows */
1321f0ae7da4SStefano Zampini   ierr = PetscLayoutMapLocal_Private(A->rmap,n,rows,&len,&lrows,NULL);CHKERRQ(ierr);
13226e520ac8SStefano Zampini   /* fix right hand side if needed */
13236e520ac8SStefano Zampini   if (x && b) {
13246e520ac8SStefano Zampini     const PetscScalar *xx;
13256e520ac8SStefano Zampini     PetscScalar       *bb;
13266e520ac8SStefano Zampini 
13276e520ac8SStefano Zampini     ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr);
13286e520ac8SStefano Zampini     ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
13296e520ac8SStefano Zampini     for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]];
13306e520ac8SStefano Zampini     ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr);
13316e520ac8SStefano Zampini     ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr);
13322e74eeadSLisandro Dalcin   }
13336e520ac8SStefano Zampini   /* get rows associated to the local matrices */
13346e520ac8SStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
13356e520ac8SStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
13366e520ac8SStefano Zampini   }
13376e520ac8SStefano Zampini   ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr);
13386e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr);
13396e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
13406e520ac8SStefano Zampini   for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1;
13416e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
13426e520ac8SStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
13436e520ac8SStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
13446e520ac8SStefano Zampini   ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr);
13456e520ac8SStefano Zampini   for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i;
1346f0ae7da4SStefano Zampini   ierr = MatISZeroRowsColumnsLocal_Private(A,nr,lrows,diag,columns);CHKERRQ(ierr);
13476e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
13482e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
13492e74eeadSLisandro Dalcin }
13502e74eeadSLisandro Dalcin 
13512e74eeadSLisandro Dalcin #undef __FUNCT__
1352f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRows_IS"
1353f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1354b4319ba4SBarry Smith {
1355dfbe8321SBarry Smith   PetscErrorCode ierr;
1356b4319ba4SBarry Smith 
1357b4319ba4SBarry Smith   PetscFunctionBegin;
1358f0ae7da4SStefano Zampini   ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_FALSE);CHKERRQ(ierr);
1359f0ae7da4SStefano Zampini   PetscFunctionReturn(0);
1360f0ae7da4SStefano Zampini }
13612205254eSKarl Rupp 
1362f0ae7da4SStefano Zampini #undef __FUNCT__
1363f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_IS"
1364f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1365f0ae7da4SStefano Zampini {
1366f0ae7da4SStefano Zampini   PetscErrorCode ierr;
1367f0ae7da4SStefano Zampini 
1368f0ae7da4SStefano Zampini   PetscFunctionBegin;
1369f0ae7da4SStefano Zampini   ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_TRUE);CHKERRQ(ierr);
1370b4319ba4SBarry Smith   PetscFunctionReturn(0);
1371b4319ba4SBarry Smith }
1372b4319ba4SBarry Smith 
1373b4319ba4SBarry Smith #undef __FUNCT__
1374b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
1375a8116848SStefano Zampini static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
1376b4319ba4SBarry Smith {
1377b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1378dfbe8321SBarry Smith   PetscErrorCode ierr;
1379b4319ba4SBarry Smith 
1380b4319ba4SBarry Smith   PetscFunctionBegin;
1381b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
1382b4319ba4SBarry Smith   PetscFunctionReturn(0);
1383b4319ba4SBarry Smith }
1384b4319ba4SBarry Smith 
1385b4319ba4SBarry Smith #undef __FUNCT__
1386b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
1387a8116848SStefano Zampini static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
1388b4319ba4SBarry Smith {
1389b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1390dfbe8321SBarry Smith   PetscErrorCode ierr;
1391b4319ba4SBarry Smith 
1392b4319ba4SBarry Smith   PetscFunctionBegin;
1393b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
1394b4319ba4SBarry Smith   PetscFunctionReturn(0);
1395b4319ba4SBarry Smith }
1396b4319ba4SBarry Smith 
1397b4319ba4SBarry Smith #undef __FUNCT__
1398b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
1399a8116848SStefano Zampini static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
1400b4319ba4SBarry Smith {
1401b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
1402b4319ba4SBarry Smith 
1403b4319ba4SBarry Smith   PetscFunctionBegin;
1404b4319ba4SBarry Smith   *local = is->A;
1405b4319ba4SBarry Smith   PetscFunctionReturn(0);
1406b4319ba4SBarry Smith }
1407b4319ba4SBarry Smith 
1408b4319ba4SBarry Smith #undef __FUNCT__
1409b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
1410b4319ba4SBarry Smith /*@
1411b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
1412b4319ba4SBarry Smith 
1413b4319ba4SBarry Smith   Input Parameter:
1414b4319ba4SBarry Smith .  mat - the matrix
1415b4319ba4SBarry Smith 
1416b4319ba4SBarry Smith   Output Parameter:
1417eb82efa4SStefano Zampini .  local - the local matrix
1418b4319ba4SBarry Smith 
1419b4319ba4SBarry Smith   Level: advanced
1420b4319ba4SBarry Smith 
1421b4319ba4SBarry Smith   Notes:
1422b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
1423b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
1424b4319ba4SBarry Smith   of the MatSetValues() operation.
1425b4319ba4SBarry Smith 
1426b4319ba4SBarry Smith .seealso: MATIS
1427b4319ba4SBarry Smith @*/
14287087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
1429b4319ba4SBarry Smith {
14304ac538c5SBarry Smith   PetscErrorCode ierr;
1431b4319ba4SBarry Smith 
1432b4319ba4SBarry Smith   PetscFunctionBegin;
14330700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1434b4319ba4SBarry Smith   PetscValidPointer(local,2);
14354ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
1436b4319ba4SBarry Smith   PetscFunctionReturn(0);
1437b4319ba4SBarry Smith }
1438b4319ba4SBarry Smith 
14393b03a366Sstefano_zampini #undef __FUNCT__
14403b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
1441a8116848SStefano Zampini static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
14423b03a366Sstefano_zampini {
14433b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
14443b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
14453b03a366Sstefano_zampini   PetscErrorCode ierr;
14463b03a366Sstefano_zampini 
14473b03a366Sstefano_zampini   PetscFunctionBegin;
14484e4c7dbeSStefano Zampini   if (is->A) {
14493b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
14503b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
1451f0ae7da4SStefano 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);
14524e4c7dbeSStefano Zampini   }
14533b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
14543b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
14553b03a366Sstefano_zampini   is->A = local;
14563b03a366Sstefano_zampini   PetscFunctionReturn(0);
14573b03a366Sstefano_zampini }
14583b03a366Sstefano_zampini 
14593b03a366Sstefano_zampini #undef __FUNCT__
14603b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
14613b03a366Sstefano_zampini /*@
1462eb82efa4SStefano Zampini     MatISSetLocalMat - Replace the local matrix stored inside a MATIS object.
14633b03a366Sstefano_zampini 
14643b03a366Sstefano_zampini   Input Parameter:
14653b03a366Sstefano_zampini .  mat - the matrix
1466eb82efa4SStefano Zampini .  local - the local matrix
14673b03a366Sstefano_zampini 
14683b03a366Sstefano_zampini   Output Parameter:
14693b03a366Sstefano_zampini 
14703b03a366Sstefano_zampini   Level: advanced
14713b03a366Sstefano_zampini 
14723b03a366Sstefano_zampini   Notes:
14733b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
14743b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
14753b03a366Sstefano_zampini 
14763b03a366Sstefano_zampini .seealso: MATIS
14773b03a366Sstefano_zampini @*/
14783b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
14793b03a366Sstefano_zampini {
14803b03a366Sstefano_zampini   PetscErrorCode ierr;
14813b03a366Sstefano_zampini 
14823b03a366Sstefano_zampini   PetscFunctionBegin;
14833b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1484b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
14853b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
14863b03a366Sstefano_zampini   PetscFunctionReturn(0);
14873b03a366Sstefano_zampini }
14883b03a366Sstefano_zampini 
14896726f965SBarry Smith #undef __FUNCT__
14906726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
1491a8116848SStefano Zampini static PetscErrorCode MatZeroEntries_IS(Mat A)
14926726f965SBarry Smith {
14936726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
14946726f965SBarry Smith   PetscErrorCode ierr;
14956726f965SBarry Smith 
14966726f965SBarry Smith   PetscFunctionBegin;
14976726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
14986726f965SBarry Smith   PetscFunctionReturn(0);
14996726f965SBarry Smith }
15006726f965SBarry Smith 
15016726f965SBarry Smith #undef __FUNCT__
15022e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
1503a8116848SStefano Zampini static PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
15042e74eeadSLisandro Dalcin {
15052e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
15062e74eeadSLisandro Dalcin   PetscErrorCode ierr;
15072e74eeadSLisandro Dalcin 
15082e74eeadSLisandro Dalcin   PetscFunctionBegin;
15092e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
15102e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
15112e74eeadSLisandro Dalcin }
15122e74eeadSLisandro Dalcin 
15132e74eeadSLisandro Dalcin #undef __FUNCT__
15142e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
1515a8116848SStefano Zampini static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
15162e74eeadSLisandro Dalcin {
15172e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
15182e74eeadSLisandro Dalcin   PetscErrorCode ierr;
15192e74eeadSLisandro Dalcin 
15202e74eeadSLisandro Dalcin   PetscFunctionBegin;
15212e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
1522e176bc59SStefano Zampini   ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr);
15232e74eeadSLisandro Dalcin 
15242e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
15252e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
1526e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1527e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
15282e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
15292e74eeadSLisandro Dalcin }
15302e74eeadSLisandro Dalcin 
15312e74eeadSLisandro Dalcin #undef __FUNCT__
15326726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
1533a8116848SStefano Zampini static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
15346726f965SBarry Smith {
15356726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
15366726f965SBarry Smith   PetscErrorCode ierr;
15376726f965SBarry Smith 
15386726f965SBarry Smith   PetscFunctionBegin;
15394e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
15406726f965SBarry Smith   PetscFunctionReturn(0);
15416726f965SBarry Smith }
15426726f965SBarry Smith 
1543284134d9SBarry Smith #undef __FUNCT__
1544f26d0771SStefano Zampini #define __FUNCT__ "MatAXPY_IS"
1545f26d0771SStefano Zampini static PetscErrorCode MatAXPY_IS(Mat Y,PetscScalar a,Mat X,MatStructure str)
1546f26d0771SStefano Zampini {
1547f26d0771SStefano Zampini   Mat_IS         *y = (Mat_IS*)Y->data;
1548f26d0771SStefano Zampini   Mat_IS         *x;
1549f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1550f26d0771SStefano Zampini   PetscBool      ismatis;
1551f26d0771SStefano Zampini #endif
1552f26d0771SStefano Zampini   PetscErrorCode ierr;
1553f26d0771SStefano Zampini 
1554f26d0771SStefano Zampini   PetscFunctionBegin;
1555f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1556f26d0771SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)X,MATIS,&ismatis);CHKERRQ(ierr);
1557f26d0771SStefano Zampini   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS");
1558f26d0771SStefano Zampini #endif
1559f26d0771SStefano Zampini   x = (Mat_IS*)X->data;
1560f26d0771SStefano Zampini   ierr = MatAXPY(y->A,a,x->A,str);CHKERRQ(ierr);
1561f26d0771SStefano Zampini   PetscFunctionReturn(0);
1562f26d0771SStefano Zampini }
1563f26d0771SStefano Zampini 
1564f26d0771SStefano Zampini #undef __FUNCT__
1565f26d0771SStefano Zampini #define __FUNCT__ "MatGetLocalSubMatrix_IS"
1566f26d0771SStefano Zampini static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat)
1567f26d0771SStefano Zampini {
1568f26d0771SStefano Zampini   Mat                    lA;
1569f26d0771SStefano Zampini   Mat_IS                 *matis;
1570f26d0771SStefano Zampini   ISLocalToGlobalMapping rl2g,cl2g;
1571f26d0771SStefano Zampini   IS                     is;
1572f26d0771SStefano Zampini   const PetscInt         *rg,*rl;
1573f26d0771SStefano Zampini   PetscInt               nrg;
1574f26d0771SStefano Zampini   PetscInt               N,M,nrl,i,*idxs;
1575f26d0771SStefano Zampini   PetscErrorCode         ierr;
1576f26d0771SStefano Zampini 
1577f26d0771SStefano Zampini   PetscFunctionBegin;
1578f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr);
1579f26d0771SStefano Zampini   ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr);
1580f26d0771SStefano Zampini   ierr = ISGetIndices(row,&rl);CHKERRQ(ierr);
1581f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr);
1582f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1583f0ae7da4SStefano 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);
1584f26d0771SStefano Zampini #endif
1585f26d0771SStefano Zampini   ierr = PetscMalloc1(nrg,&idxs);CHKERRQ(ierr);
1586f26d0771SStefano Zampini   /* map from [0,nrl) to row */
1587f26d0771SStefano Zampini   for (i=0;i<nrl;i++) idxs[i] = rl[i];
1588f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1589f26d0771SStefano Zampini   for (i=nrl;i<nrg;i++) idxs[i] = nrg;
1590f26d0771SStefano Zampini #else
1591f26d0771SStefano Zampini   for (i=nrl;i<nrg;i++) idxs[i] = -1;
1592f26d0771SStefano Zampini #endif
1593f26d0771SStefano Zampini   ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr);
1594f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr);
1595f26d0771SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
1596f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
1597f26d0771SStefano Zampini   ierr = ISDestroy(&is);CHKERRQ(ierr);
1598f26d0771SStefano Zampini   /* compute new l2g map for columns */
1599f26d0771SStefano Zampini   if (col != row || A->rmap->mapping != A->cmap->mapping) {
1600f26d0771SStefano Zampini     const PetscInt *cg,*cl;
1601f26d0771SStefano Zampini     PetscInt       ncg;
1602f26d0771SStefano Zampini     PetscInt       ncl;
1603f26d0771SStefano Zampini 
1604f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr);
1605f26d0771SStefano Zampini     ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr);
1606f26d0771SStefano Zampini     ierr = ISGetIndices(col,&cl);CHKERRQ(ierr);
1607f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr);
1608f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1609f0ae7da4SStefano 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);
1610f26d0771SStefano Zampini #endif
1611f26d0771SStefano Zampini     ierr = PetscMalloc1(ncg,&idxs);CHKERRQ(ierr);
1612f26d0771SStefano Zampini     /* map from [0,ncl) to col */
1613f26d0771SStefano Zampini     for (i=0;i<ncl;i++) idxs[i] = cl[i];
1614f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1615f26d0771SStefano Zampini     for (i=ncl;i<ncg;i++) idxs[i] = ncg;
1616f26d0771SStefano Zampini #else
1617f26d0771SStefano Zampini     for (i=ncl;i<ncg;i++) idxs[i] = -1;
1618f26d0771SStefano Zampini #endif
1619f26d0771SStefano Zampini     ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr);
1620f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr);
1621f26d0771SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
1622f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
1623f26d0771SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
1624f26d0771SStefano Zampini   } else {
1625f26d0771SStefano Zampini     ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr);
1626f26d0771SStefano Zampini     cl2g = rl2g;
1627f26d0771SStefano Zampini   }
1628f26d0771SStefano Zampini   /* create the MATIS submatrix */
1629f26d0771SStefano Zampini   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1630f26d0771SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr);
1631f26d0771SStefano Zampini   ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1632f26d0771SStefano Zampini   ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr);
1633b0aa3428SStefano Zampini   matis = (Mat_IS*)((*submat)->data);
1634f26d0771SStefano Zampini   matis->islocalref = PETSC_TRUE;
1635f26d0771SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr);
1636f26d0771SStefano Zampini   ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr);
1637f26d0771SStefano Zampini   ierr = MatISSetLocalMat(*submat,lA);CHKERRQ(ierr);
1638f26d0771SStefano Zampini   ierr = MatSetUp(*submat);CHKERRQ(ierr);
1639f26d0771SStefano Zampini   ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1640f26d0771SStefano Zampini   ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1641f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
1642f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
1643f26d0771SStefano Zampini   /* remove unsupported ops */
1644f26d0771SStefano Zampini   ierr = PetscMemzero((*submat)->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1645f26d0771SStefano Zampini   (*submat)->ops->destroy               = MatDestroy_IS;
1646f26d0771SStefano Zampini   (*submat)->ops->setvalueslocal        = MatSetValuesLocal_SubMat_IS;
1647f26d0771SStefano Zampini   (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS;
1648f26d0771SStefano Zampini   (*submat)->ops->assemblybegin         = MatAssemblyBegin_IS;
1649f26d0771SStefano Zampini   (*submat)->ops->assemblyend           = MatAssemblyEnd_IS;
1650f26d0771SStefano Zampini   PetscFunctionReturn(0);
1651f26d0771SStefano Zampini }
1652f26d0771SStefano Zampini 
1653f26d0771SStefano Zampini #undef __FUNCT__
1654284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
1655284134d9SBarry Smith /*@
16563c212e90SHong Zhang     MatCreateIS - Creates a "process" unassembled matrix, assembled on each
1657284134d9SBarry Smith        process but not across processes.
1658284134d9SBarry Smith 
1659284134d9SBarry Smith    Input Parameters:
1660284134d9SBarry Smith +     comm    - MPI communicator that will share the matrix
1661e176bc59SStefano Zampini .     bs      - block size of the matrix
1662df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
1663e176bc59SStefano Zampini .     rmap    - local to global map for rows
1664e176bc59SStefano Zampini -     cmap    - local to global map for cols
1665284134d9SBarry Smith 
1666284134d9SBarry Smith    Output Parameter:
1667284134d9SBarry Smith .    A - the resulting matrix
1668284134d9SBarry Smith 
16698e6c10adSSatish Balay    Level: advanced
16708e6c10adSSatish Balay 
16713c212e90SHong Zhang    Notes: See MATIS for more details.
16726fdf41d1SStefano Zampini           m and n are NOT related to the size of the map; they represent the size of the local parts of the vectors
16736fdf41d1SStefano Zampini           used in MatMult operations. The sizes of rmap and cmap define the size of the local matrices.
16743c212e90SHong Zhang           If either rmap or cmap are NULL, then the matrix is assumed to be square.
1675284134d9SBarry Smith 
1676284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
1677284134d9SBarry Smith @*/
1678e176bc59SStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A)
1679284134d9SBarry Smith {
1680284134d9SBarry Smith   PetscErrorCode ierr;
1681284134d9SBarry Smith 
1682284134d9SBarry Smith   PetscFunctionBegin;
16836fdf41d1SStefano Zampini   if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mappings");
1684284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
16856fdf41d1SStefano Zampini   if (bs > 0) {
1686d3ee2243SStefano Zampini     ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
16876fdf41d1SStefano Zampini   }
1688284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1689284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
1690e176bc59SStefano Zampini   if (rmap && cmap) {
1691e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr);
1692e176bc59SStefano Zampini   } else if (!rmap) {
1693e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr);
1694e176bc59SStefano Zampini   } else {
1695e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr);
1696e176bc59SStefano Zampini   }
1697284134d9SBarry Smith   PetscFunctionReturn(0);
1698284134d9SBarry Smith }
1699284134d9SBarry Smith 
1700b4319ba4SBarry Smith /*MC
1701f26d0771SStefano Zampini    MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition methods (e.g. PCBDDC or KSPFETIDP).
1702b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
1703b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
1704b4319ba4SBarry Smith    product is handled "implicitly".
1705b4319ba4SBarry Smith 
1706b4319ba4SBarry Smith    Operations Provided:
17076726f965SBarry Smith +  MatMult()
17082e74eeadSLisandro Dalcin .  MatMultAdd()
17092e74eeadSLisandro Dalcin .  MatMultTranspose()
17102e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
17116726f965SBarry Smith .  MatZeroEntries()
17126726f965SBarry Smith .  MatSetOption()
17132e74eeadSLisandro Dalcin .  MatZeroRows()
17142e74eeadSLisandro Dalcin .  MatSetValues()
171597563a80SStefano Zampini .  MatSetValuesBlocked()
17166726f965SBarry Smith .  MatSetValuesLocal()
171797563a80SStefano Zampini .  MatSetValuesBlockedLocal()
17182e74eeadSLisandro Dalcin .  MatScale()
17192e74eeadSLisandro Dalcin .  MatGetDiagonal()
17202b404112SStefano Zampini .  MatMissingDiagonal()
17212b404112SStefano Zampini .  MatDuplicate()
17222b404112SStefano Zampini .  MatCopy()
1723f26d0771SStefano Zampini .  MatAXPY()
1724f26d0771SStefano Zampini .  MatGetSubMatrix()
1725f26d0771SStefano Zampini .  MatGetLocalSubMatrix()
1726d7f69cd0SStefano Zampini .  MatTranspose()
17276726f965SBarry Smith -  MatSetLocalToGlobalMapping()
1728b4319ba4SBarry Smith 
1729b4319ba4SBarry Smith    Options Database Keys:
1730b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
1731b4319ba4SBarry Smith 
1732b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
1733b4319ba4SBarry Smith 
1734b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
1735b4319ba4SBarry Smith 
1736b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
1737eb82efa4SStefano Zampini           MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation()
1738b4319ba4SBarry Smith 
1739b4319ba4SBarry Smith   Level: advanced
1740b4319ba4SBarry Smith 
1741f26d0771SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), MatCreateIS(), PCBDDC, KSPFETIDP
1742b4319ba4SBarry Smith 
1743b4319ba4SBarry Smith M*/
1744b4319ba4SBarry Smith 
1745b4319ba4SBarry Smith #undef __FUNCT__
1746b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
17478cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
1748b4319ba4SBarry Smith {
1749dfbe8321SBarry Smith   PetscErrorCode ierr;
1750b4319ba4SBarry Smith   Mat_IS         *b;
1751b4319ba4SBarry Smith 
1752b4319ba4SBarry Smith   PetscFunctionBegin;
1753b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1754b4319ba4SBarry Smith   A->data = (void*)b;
1755b4319ba4SBarry Smith 
1756e176bc59SStefano Zampini   /* matrix ops */
1757e176bc59SStefano Zampini   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1758b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
17592e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
17602e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
17612e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1762b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1763b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
17642e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
176598921651SStefano Zampini   A->ops->setvaluesblocked        = MatSetValuesBlocked_IS;
1766b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1767f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
17682e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1769f0ae7da4SStefano Zampini   A->ops->zerorowscolumns         = MatZeroRowsColumns_IS;
1770b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1771b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1772b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
17736726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
17742e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
17752e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
17766726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
177769796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
177869796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1779ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
17806bd84002SStefano Zampini   A->ops->missingdiagonal         = MatMissingDiagonal_IS;
17812b404112SStefano Zampini   A->ops->copy                    = MatCopy_IS;
1782659959c5SStefano Zampini   A->ops->getlocalsubmatrix       = MatGetLocalSubMatrix_IS;
1783a8116848SStefano Zampini   A->ops->getsubmatrix            = MatGetSubMatrix_IS;
1784f26d0771SStefano Zampini   A->ops->axpy                    = MatAXPY_IS;
17853fd1c9e7SStefano Zampini   A->ops->diagonalset             = MatDiagonalSet_IS;
17863fd1c9e7SStefano Zampini   A->ops->shift                   = MatShift_IS;
1787d7f69cd0SStefano Zampini   A->ops->transpose               = MatTranspose_IS;
17887fa8f2d3SStefano Zampini   A->ops->getinfo                 = MatGetInfo_IS;
1789b4319ba4SBarry Smith 
1790b7ce53b6SStefano Zampini   /* special MATIS functions */
1791bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1792bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1793b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
17942e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
179517667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1796b4319ba4SBarry Smith   PetscFunctionReturn(0);
1797b4319ba4SBarry Smith }
1798