xref: /petsc/src/mat/impls/is/matis.c (revision bb1015c3d26b4b687e2446ad76998fe3eecfe9bd)
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);
157230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat,PetscInt,const PetscInt[],PetscScalar);
16a72627d2SStefano Zampini 
1728f4e0baSStefano Zampini #undef __FUNCT__
18f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesLocal_SubMat_IS"
19f26d0771SStefano Zampini static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
20f26d0771SStefano Zampini {
21f26d0771SStefano Zampini   PetscErrorCode ierr;
22f26d0771SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
23f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
24f26d0771SStefano Zampini 
25f26d0771SStefano Zampini   PetscFunctionBegin;
26f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
27f26d0771SStefano 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);
28f26d0771SStefano Zampini #endif
29f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr);
30f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr);
31f26d0771SStefano Zampini   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
32f26d0771SStefano Zampini   PetscFunctionReturn(0);
33f26d0771SStefano Zampini }
34f26d0771SStefano Zampini 
35f26d0771SStefano Zampini #undef __FUNCT__
36f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesBlockedLocal_SubMat_IS"
37f26d0771SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
38f26d0771SStefano Zampini {
39f26d0771SStefano Zampini   PetscErrorCode ierr;
40f26d0771SStefano Zampini   Mat_IS         *is = (Mat_IS*)A->data;
41f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
42f26d0771SStefano Zampini 
43f26d0771SStefano Zampini   PetscFunctionBegin;
44f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
45f26d0771SStefano 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);
46f26d0771SStefano Zampini #endif
47f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr);
48f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr);
49f26d0771SStefano Zampini   ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
50f26d0771SStefano Zampini   PetscFunctionReturn(0);
51f26d0771SStefano Zampini }
52f26d0771SStefano Zampini 
53f26d0771SStefano Zampini #undef __FUNCT__
54a8116848SStefano Zampini #define __FUNCT__ "PetscLayoutMapLocal_Private"
55a8116848SStefano Zampini static PetscErrorCode PetscLayoutMapLocal_Private(PetscLayout map,PetscInt N,const PetscInt idxs[],const PetscInt gidxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs)
56a8116848SStefano Zampini {
57a8116848SStefano Zampini   PetscInt      *owners = map->range;
58a8116848SStefano Zampini   PetscInt       n      = map->n;
59a8116848SStefano Zampini   PetscSF        sf;
60a8116848SStefano Zampini   PetscInt      *lidxs,*work = NULL;
61a8116848SStefano Zampini   PetscSFNode   *ridxs;
62a8116848SStefano Zampini   PetscMPIInt    rank;
63a8116848SStefano Zampini   PetscInt       r, p = 0, len = 0;
64a8116848SStefano Zampini   PetscErrorCode ierr;
65a8116848SStefano Zampini 
66a8116848SStefano Zampini   PetscFunctionBegin;
67a8116848SStefano Zampini   /* Create SF where leaves are input idxs and roots are owned idxs (code adapted from MatZeroRowsMapLocal_Private) */
68a8116848SStefano Zampini   ierr = MPI_Comm_rank(map->comm,&rank);CHKERRQ(ierr);
69a8116848SStefano Zampini   ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr);
70a8116848SStefano Zampini   for (r = 0; r < n; ++r) lidxs[r] = -1;
71a8116848SStefano Zampini   ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr);
72a8116848SStefano Zampini   for (r = 0; r < N; ++r) {
73a8116848SStefano Zampini     const PetscInt idx = idxs[r];
74a8116848SStefano 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);
75a8116848SStefano Zampini     if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */
76a8116848SStefano Zampini       ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr);
77a8116848SStefano Zampini     }
78a8116848SStefano Zampini     ridxs[r].rank = p;
79a8116848SStefano Zampini     ridxs[r].index = idxs[r] - owners[p];
80a8116848SStefano Zampini   }
81a8116848SStefano Zampini   ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr);
82a8116848SStefano Zampini   ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr);
83a8116848SStefano Zampini   ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr);
84a8116848SStefano Zampini   ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr);
85a8116848SStefano Zampini   if ((gidxs && ogidxs) || ogidxs) {
86a8116848SStefano Zampini     ierr = PetscMalloc1(n,&work);CHKERRQ(ierr);
87a8116848SStefano Zampini     if (!gidxs) { /* if gidxs is not provided, global indices are inferred */
88a8116848SStefano Zampini       PetscInt cum = 0,start,*work2;
89a8116848SStefano Zampini       ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr);
90a8116848SStefano Zampini       for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++;
91a8116848SStefano Zampini       ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRQ(ierr);
92a8116848SStefano Zampini       start -= cum;
93a8116848SStefano Zampini       cum = 0;
94a8116848SStefano Zampini       for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++;
95a8116848SStefano Zampini       ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr);
96a8116848SStefano Zampini       ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr);
97a8116848SStefano Zampini       ierr = PetscFree(work2);CHKERRQ(ierr);
98a8116848SStefano Zampini     } else {
99a8116848SStefano Zampini       ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)gidxs,work,MPIU_REPLACE);CHKERRQ(ierr);
100a8116848SStefano Zampini       ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)gidxs,work,MPIU_REPLACE);CHKERRQ(ierr);
101a8116848SStefano Zampini     }
102a8116848SStefano Zampini   }
103a8116848SStefano Zampini   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
104a8116848SStefano Zampini   /* Compress and put in indices */
105a8116848SStefano Zampini   for (r = 0; r < n; ++r)
106a8116848SStefano Zampini     if (lidxs[r] >= 0) {
107a8116848SStefano Zampini       if (work) work[len] = work[r];
108a8116848SStefano Zampini       lidxs[len++] = r;
109a8116848SStefano Zampini     }
110a8116848SStefano Zampini   if (on) *on = len;
111a8116848SStefano Zampini   if (oidxs) *oidxs = lidxs;
112a8116848SStefano Zampini   if (ogidxs) *ogidxs = work;
113a8116848SStefano Zampini   PetscFunctionReturn(0);
114a8116848SStefano Zampini }
115a8116848SStefano Zampini 
116a8116848SStefano Zampini #undef __FUNCT__
117a8116848SStefano Zampini #define __FUNCT__ "MatGetSubMatrix_IS"
118a8116848SStefano Zampini static PetscErrorCode MatGetSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat)
119a8116848SStefano Zampini {
120a8116848SStefano Zampini   Mat               locmat,newlocmat;
121a8116848SStefano Zampini   Mat_IS            *newmatis;
122a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG)
123a8116848SStefano Zampini   Vec               rtest,ltest;
124a8116848SStefano Zampini   const PetscScalar *array;
125a8116848SStefano Zampini #endif
126a8116848SStefano Zampini   const PetscInt    *idxs;
127a8116848SStefano Zampini   PetscInt          i,m,n;
128a8116848SStefano Zampini   PetscErrorCode    ierr;
129a8116848SStefano Zampini 
130a8116848SStefano Zampini   PetscFunctionBegin;
131a8116848SStefano Zampini   if (scall == MAT_REUSE_MATRIX) {
132a8116848SStefano Zampini     PetscBool ismatis;
133a8116848SStefano Zampini 
134a8116848SStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr);
135a8116848SStefano Zampini     if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type");
136a8116848SStefano Zampini     newmatis = (Mat_IS*)(*newmat)->data;
137a8116848SStefano Zampini     if (!newmatis->getsub_ris) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS");
138a8116848SStefano Zampini     if (!newmatis->getsub_cis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS");
139a8116848SStefano Zampini   }
140a8116848SStefano Zampini   /* irow and icol may not have duplicate entries */
141a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG)
142a8116848SStefano Zampini   ierr = MatCreateVecs(mat,&ltest,&rtest);CHKERRQ(ierr);
143a8116848SStefano Zampini   ierr = ISGetLocalSize(irow,&n);CHKERRQ(ierr);
144a8116848SStefano Zampini   ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr);
145a8116848SStefano Zampini   for (i=0;i<n;i++) {
146a8116848SStefano Zampini     ierr = VecSetValue(rtest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr);
147a8116848SStefano Zampini   }
148a8116848SStefano Zampini   ierr = VecAssemblyBegin(rtest);CHKERRQ(ierr);
149a8116848SStefano Zampini   ierr = VecAssemblyEnd(rtest);CHKERRQ(ierr);
150a8116848SStefano Zampini   ierr = VecGetLocalSize(rtest,&n);CHKERRQ(ierr);
151a8116848SStefano Zampini   ierr = VecGetOwnershipRange(rtest,&m,NULL);CHKERRQ(ierr);
152a8116848SStefano Zampini   ierr = VecGetArrayRead(rtest,&array);CHKERRQ(ierr);
153fd479f66SMatthew 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]));
154a8116848SStefano Zampini   ierr = VecRestoreArrayRead(rtest,&array);CHKERRQ(ierr);
155a8116848SStefano Zampini   ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr);
156a8116848SStefano Zampini   ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr);
157a8116848SStefano Zampini   ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr);
158a8116848SStefano Zampini   for (i=0;i<n;i++) {
159a8116848SStefano Zampini     ierr = VecSetValue(ltest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr);
160a8116848SStefano Zampini   }
161a8116848SStefano Zampini   ierr = VecAssemblyBegin(ltest);CHKERRQ(ierr);
162a8116848SStefano Zampini   ierr = VecAssemblyEnd(ltest);CHKERRQ(ierr);
163a8116848SStefano Zampini   ierr = VecGetLocalSize(ltest,&n);CHKERRQ(ierr);
164a8116848SStefano Zampini   ierr = VecGetOwnershipRange(ltest,&m,NULL);CHKERRQ(ierr);
165a8116848SStefano Zampini   ierr = VecGetArrayRead(ltest,&array);CHKERRQ(ierr);
166fd479f66SMatthew 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]));
167a8116848SStefano Zampini   ierr = VecRestoreArrayRead(ltest,&array);CHKERRQ(ierr);
168a8116848SStefano Zampini   ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr);
169a8116848SStefano Zampini   ierr = VecDestroy(&rtest);CHKERRQ(ierr);
170a8116848SStefano Zampini   ierr = VecDestroy(&ltest);CHKERRQ(ierr);
171a8116848SStefano Zampini #endif
172a8116848SStefano Zampini   if (scall == MAT_INITIAL_MATRIX) {
173a8116848SStefano Zampini     Mat_IS                 *matis = (Mat_IS*)mat->data;
174a8116848SStefano Zampini     ISLocalToGlobalMapping rl2g;
175a8116848SStefano Zampini     IS                     is;
176a8116848SStefano Zampini     PetscInt               *lidxs,*lgidxs,*newgidxs;
1776dd40735SStefano Zampini     PetscInt               ll,newloc;
178a8116848SStefano Zampini     MPI_Comm               comm;
179a8116848SStefano Zampini 
180a8116848SStefano Zampini     ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
181a8116848SStefano Zampini     ierr = ISGetLocalSize(irow,&m);CHKERRQ(ierr);
182a8116848SStefano Zampini     ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr);
183a8116848SStefano Zampini     ierr = MatCreate(comm,newmat);CHKERRQ(ierr);
184a8116848SStefano Zampini     ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr);
185a8116848SStefano Zampini     ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
186a8116848SStefano Zampini     /* communicate irow to their owners in the layout */
187a8116848SStefano Zampini     ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr);
188a8116848SStefano Zampini     ierr = PetscLayoutMapLocal_Private(mat->rmap,m,idxs,NULL,&ll,&lidxs,&lgidxs);CHKERRQ(ierr);
189a8116848SStefano Zampini     ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr);
190a8116848SStefano Zampini     if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
191a8116848SStefano Zampini       ierr = MatISComputeSF_Private(mat);CHKERRQ(ierr);
192a8116848SStefano Zampini     }
193a8116848SStefano Zampini     ierr = PetscMemzero(matis->sf_rootdata,matis->sf_nroots*sizeof(PetscInt));CHKERRQ(ierr);
194a8116848SStefano Zampini     for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1;
195a8116848SStefano Zampini     ierr = PetscFree(lidxs);CHKERRQ(ierr);
196a8116848SStefano Zampini     ierr = PetscFree(lgidxs);CHKERRQ(ierr);
197a8116848SStefano Zampini     ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
198a8116848SStefano Zampini     ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
199a8116848SStefano Zampini     for (i=0,newloc=0;i<matis->sf_nleaves;i++) if (matis->sf_leafdata[i]) newloc++;
200a8116848SStefano Zampini     ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr);
201a8116848SStefano Zampini     ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr);
202a8116848SStefano Zampini     for (i=0,newloc=0;i<matis->sf_nleaves;i++)
203a8116848SStefano Zampini       if (matis->sf_leafdata[i]) {
204a8116848SStefano Zampini         lidxs[newloc] = i;
205a8116848SStefano Zampini         newgidxs[newloc++] = matis->sf_leafdata[i]-1;
206a8116848SStefano Zampini       }
207a8116848SStefano Zampini     ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
208a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
209a8116848SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
210a8116848SStefano Zampini     /* local is to extract local submatrix */
211a8116848SStefano Zampini     newmatis = (Mat_IS*)(*newmat)->data;
212a8116848SStefano Zampini     ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris);CHKERRQ(ierr);
213a8116848SStefano Zampini     if (mat->congruentlayouts == -1) { /* first time we compare rows and cols layouts */
214a8116848SStefano Zampini       PetscBool cong;
215a8116848SStefano Zampini       ierr = PetscLayoutCompare(mat->rmap,mat->cmap,&cong);CHKERRQ(ierr);
216a8116848SStefano Zampini       if (cong) mat->congruentlayouts = 1;
217a8116848SStefano Zampini       else      mat->congruentlayouts = 0;
218a8116848SStefano Zampini     }
219a8116848SStefano Zampini     if (mat->congruentlayouts && irow == icol) {
220a8116848SStefano Zampini       ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g);CHKERRQ(ierr);
221a8116848SStefano Zampini       ierr = PetscObjectReference((PetscObject)newmatis->getsub_ris);CHKERRQ(ierr);
222a8116848SStefano Zampini       newmatis->getsub_cis = newmatis->getsub_ris;
223a8116848SStefano Zampini     } else {
224a8116848SStefano Zampini       ISLocalToGlobalMapping cl2g;
225a8116848SStefano Zampini 
226a8116848SStefano Zampini       /* communicate icol to their owners in the layout */
227a8116848SStefano Zampini       ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr);
228a8116848SStefano Zampini       ierr = PetscLayoutMapLocal_Private(mat->cmap,n,idxs,NULL,&ll,&lidxs,&lgidxs);CHKERRQ(ierr);
229a8116848SStefano Zampini       ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr);
230a8116848SStefano Zampini       ierr = PetscMemzero(matis->csf_rootdata,matis->csf_nroots*sizeof(PetscInt));CHKERRQ(ierr);
231a8116848SStefano Zampini       for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1;
232a8116848SStefano Zampini       ierr = PetscFree(lidxs);CHKERRQ(ierr);
233a8116848SStefano Zampini       ierr = PetscFree(lgidxs);CHKERRQ(ierr);
234a8116848SStefano Zampini       ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr);
235a8116848SStefano Zampini       ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr);
236a8116848SStefano Zampini       for (i=0,newloc=0;i<matis->csf_nleaves;i++) if (matis->csf_leafdata[i]) newloc++;
237a8116848SStefano Zampini       ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr);
238a8116848SStefano Zampini       ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr);
239a8116848SStefano Zampini       for (i=0,newloc=0;i<matis->csf_nleaves;i++)
240a8116848SStefano Zampini         if (matis->csf_leafdata[i]) {
241a8116848SStefano Zampini           lidxs[newloc] = i;
242a8116848SStefano Zampini           newgidxs[newloc++] = matis->csf_leafdata[i]-1;
243a8116848SStefano Zampini         }
244a8116848SStefano Zampini       ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
245a8116848SStefano Zampini       ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
246a8116848SStefano Zampini       ierr = ISDestroy(&is);CHKERRQ(ierr);
247a8116848SStefano Zampini       /* local is to extract local submatrix */
248a8116848SStefano Zampini       ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis);CHKERRQ(ierr);
249a8116848SStefano Zampini       ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr);
250a8116848SStefano Zampini       ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
251a8116848SStefano Zampini     }
252a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
253a8116848SStefano Zampini   } else {
254a8116848SStefano Zampini     ierr = MatISGetLocalMat(*newmat,&newlocmat);CHKERRQ(ierr);
255a8116848SStefano Zampini   }
256a8116848SStefano Zampini   ierr = MatISGetLocalMat(mat,&locmat);CHKERRQ(ierr);
257a8116848SStefano Zampini   newmatis = (Mat_IS*)(*newmat)->data;
258a8116848SStefano Zampini   ierr = MatGetSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat);CHKERRQ(ierr);
259a8116848SStefano Zampini   if (scall == MAT_INITIAL_MATRIX) {
260a8116848SStefano Zampini     ierr = MatISSetLocalMat(*newmat,newlocmat);CHKERRQ(ierr);
261a8116848SStefano Zampini     ierr = MatDestroy(&newlocmat);CHKERRQ(ierr);
262a8116848SStefano Zampini   }
263a8116848SStefano Zampini   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
264a8116848SStefano Zampini   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
265a8116848SStefano Zampini   PetscFunctionReturn(0);
266a8116848SStefano Zampini }
267a8116848SStefano Zampini 
268a8116848SStefano Zampini #undef __FUNCT__
2692b404112SStefano Zampini #define __FUNCT__ "MatCopy_IS"
270a8116848SStefano Zampini static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str)
2712b404112SStefano Zampini {
2722b404112SStefano Zampini   Mat_IS         *a = (Mat_IS*)A->data,*b;
2732b404112SStefano Zampini   PetscBool      ismatis;
2742b404112SStefano Zampini   PetscErrorCode ierr;
2752b404112SStefano Zampini 
2762b404112SStefano Zampini   PetscFunctionBegin;
2772b404112SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr);
2782b404112SStefano Zampini   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented");
2792b404112SStefano Zampini   b = (Mat_IS*)B->data;
2802b404112SStefano Zampini   ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr);
2812b404112SStefano Zampini   PetscFunctionReturn(0);
2822b404112SStefano Zampini }
2832b404112SStefano Zampini 
2842b404112SStefano Zampini #undef __FUNCT__
2856bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS"
286a8116848SStefano Zampini static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool  *missing,PetscInt *d)
2876bd84002SStefano Zampini {
288527b2640SStefano Zampini   Vec               v;
289527b2640SStefano Zampini   const PetscScalar *array;
290527b2640SStefano Zampini   PetscInt          i,n;
2916bd84002SStefano Zampini   PetscErrorCode    ierr;
2926bd84002SStefano Zampini 
2936bd84002SStefano Zampini   PetscFunctionBegin;
294527b2640SStefano Zampini   *missing = PETSC_FALSE;
295527b2640SStefano Zampini   ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr);
296527b2640SStefano Zampini   ierr = MatGetDiagonal(A,v);CHKERRQ(ierr);
297527b2640SStefano Zampini   ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr);
298527b2640SStefano Zampini   ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr);
299527b2640SStefano Zampini   for (i=0;i<n;i++) if (array[i] == 0.) break;
300527b2640SStefano Zampini   ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr);
301527b2640SStefano Zampini   ierr = VecDestroy(&v);CHKERRQ(ierr);
302527b2640SStefano Zampini   if (i != n) *missing = PETSC_TRUE;
303527b2640SStefano Zampini   if (d) {
304527b2640SStefano Zampini     *d = -1;
305527b2640SStefano Zampini     if (*missing) {
306527b2640SStefano Zampini       PetscInt rstart;
307527b2640SStefano Zampini       ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
308527b2640SStefano Zampini       *d = i+rstart;
309527b2640SStefano Zampini     }
310527b2640SStefano Zampini   }
3116bd84002SStefano Zampini   PetscFunctionReturn(0);
3126bd84002SStefano Zampini }
3136bd84002SStefano Zampini 
3146bd84002SStefano Zampini #undef __FUNCT__
31528f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private"
316a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B)
31728f4e0baSStefano Zampini {
31828f4e0baSStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
31928f4e0baSStefano Zampini   const PetscInt *gidxs;
32028f4e0baSStefano Zampini   PetscErrorCode ierr;
32128f4e0baSStefano Zampini 
32228f4e0baSStefano Zampini   PetscFunctionBegin;
32328f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr);
32428f4e0baSStefano Zampini   ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr);
32528f4e0baSStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr);
3263bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
32728f4e0baSStefano Zampini   ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
3283bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
32928f4e0baSStefano Zampini   ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr);
330a8116848SStefano Zampini   if (B->rmap->mapping != B->cmap->mapping) { /* setup SF for columns */
331a8116848SStefano Zampini     ierr = MatGetSize(matis->A,NULL,&matis->csf_nleaves);CHKERRQ(ierr);
332a8116848SStefano Zampini     ierr = MatGetLocalSize(B,NULL,&matis->csf_nroots);CHKERRQ(ierr);
333a8116848SStefano Zampini     ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf);CHKERRQ(ierr);
334a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr);
335a8116848SStefano Zampini     ierr = PetscSFSetGraphLayout(matis->csf,B->cmap,matis->csf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
336a8116848SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr);
337a8116848SStefano Zampini     ierr = PetscMalloc2(matis->csf_nroots,&matis->csf_rootdata,matis->csf_nleaves,&matis->csf_leafdata);CHKERRQ(ierr);
338a8116848SStefano Zampini   } else {
339a8116848SStefano Zampini     matis->csf = matis->sf;
340a8116848SStefano Zampini     matis->csf_nleaves = matis->sf_nleaves;
341a8116848SStefano Zampini     matis->csf_nroots = matis->sf_nroots;
342a8116848SStefano Zampini     matis->csf_leafdata = matis->sf_leafdata;
343a8116848SStefano Zampini     matis->csf_rootdata = matis->sf_rootdata;
344a8116848SStefano Zampini   }
34528f4e0baSStefano Zampini   PetscFunctionReturn(0);
34628f4e0baSStefano Zampini }
3472e1947a5SStefano Zampini 
3482e1947a5SStefano Zampini #undef __FUNCT__
3492e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation"
350eb82efa4SStefano Zampini /*@
351a88811baSStefano Zampini    MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix.
352a88811baSStefano Zampini 
353a88811baSStefano Zampini    Collective on MPI_Comm
354a88811baSStefano Zampini 
355a88811baSStefano Zampini    Input Parameters:
356a88811baSStefano Zampini +  B - the matrix
357a88811baSStefano Zampini .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
358a88811baSStefano Zampini            (same value is used for all local rows)
359a88811baSStefano Zampini .  d_nnz - array containing the number of nonzeros in the various rows of the
360a88811baSStefano Zampini            DIAGONAL portion of the local submatrix (possibly different for each row)
361a88811baSStefano Zampini            or NULL, if d_nz is used to specify the nonzero structure.
362a88811baSStefano Zampini            The size of this array is equal to the number of local rows, i.e 'm'.
363a88811baSStefano Zampini            For matrices that will be factored, you must leave room for (and set)
364a88811baSStefano Zampini            the diagonal entry even if it is zero.
365a88811baSStefano Zampini .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
366a88811baSStefano Zampini            submatrix (same value is used for all local rows).
367a88811baSStefano Zampini -  o_nnz - array containing the number of nonzeros in the various rows of the
368a88811baSStefano Zampini            OFF-DIAGONAL portion of the local submatrix (possibly different for
369a88811baSStefano Zampini            each row) or NULL, if o_nz is used to specify the nonzero
370a88811baSStefano Zampini            structure. The size of this array is equal to the number
371a88811baSStefano Zampini            of local rows, i.e 'm'.
372a88811baSStefano Zampini 
373a88811baSStefano Zampini    If the *_nnz parameter is given then the *_nz parameter is ignored
374a88811baSStefano Zampini 
375a88811baSStefano Zampini    Level: intermediate
376a88811baSStefano Zampini 
377a88811baSStefano Zampini    Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition
378a88811baSStefano Zampini           from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local
379a88811baSStefano Zampini           matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects.
380a88811baSStefano Zampini 
381a88811baSStefano Zampini .keywords: matrix
382a88811baSStefano Zampini 
3833c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS
384a88811baSStefano Zampini @*/
3852e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
3862e1947a5SStefano Zampini {
3872e1947a5SStefano Zampini   PetscErrorCode ierr;
3882e1947a5SStefano Zampini 
3892e1947a5SStefano Zampini   PetscFunctionBegin;
3902e1947a5SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
3912e1947a5SStefano Zampini   PetscValidType(B,1);
3922e1947a5SStefano Zampini   ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
3932e1947a5SStefano Zampini   PetscFunctionReturn(0);
3942e1947a5SStefano Zampini }
3952e1947a5SStefano Zampini 
3962e1947a5SStefano Zampini #undef __FUNCT__
3972e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS"
3987230de76SStefano Zampini static PetscErrorCode  MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
3992e1947a5SStefano Zampini {
4002e1947a5SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
40128f4e0baSStefano Zampini   PetscInt       bs,i,nlocalcols;
4022e1947a5SStefano Zampini   PetscErrorCode ierr;
4032e1947a5SStefano Zampini 
4042e1947a5SStefano Zampini   PetscFunctionBegin;
4056c4ed002SBarry Smith   if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
40628f4e0baSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
40728f4e0baSStefano Zampini     ierr = MatISComputeSF_Private(B);CHKERRQ(ierr);
40828f4e0baSStefano Zampini   }
4092e1947a5SStefano Zampini   if (!d_nnz) {
41028f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz;
4112e1947a5SStefano Zampini   } else {
41228f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i];
4132e1947a5SStefano Zampini   }
4142e1947a5SStefano Zampini   if (!o_nnz) {
41528f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz;
4162e1947a5SStefano Zampini   } else {
41728f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i];
4182e1947a5SStefano Zampini   }
41928f4e0baSStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
42028f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr);
42128f4e0baSStefano Zampini   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
42228f4e0baSStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
42328f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves;i++) {
42428f4e0baSStefano Zampini     matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols);
4252e1947a5SStefano Zampini   }
42628f4e0baSStefano Zampini   ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr);
42728f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
42828f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs;
4292e1947a5SStefano Zampini   }
43028f4e0baSStefano Zampini   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
43128f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
43228f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i]-i;
4332e1947a5SStefano Zampini   }
43428f4e0baSStefano Zampini   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
4352e1947a5SStefano Zampini   PetscFunctionReturn(0);
4362e1947a5SStefano Zampini }
437b4319ba4SBarry Smith 
438b4319ba4SBarry Smith #undef __FUNCT__
4393927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private"
4403927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
4413927de2eSStefano Zampini {
4423927de2eSStefano Zampini   Mat_IS          *matis = (Mat_IS*)(A->data);
4433927de2eSStefano Zampini   PetscInt        *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership;
444ecf5a873SStefano Zampini   const PetscInt  *global_indices_r,*global_indices_c;
4453927de2eSStefano Zampini   PetscInt        i,j,bs,rows,cols;
4463927de2eSStefano Zampini   PetscInt        lrows,lcols;
4473927de2eSStefano Zampini   PetscInt        local_rows,local_cols;
4483927de2eSStefano Zampini   PetscMPIInt     nsubdomains;
4493927de2eSStefano Zampini   PetscBool       isdense,issbaij;
4503927de2eSStefano Zampini   PetscErrorCode  ierr;
4513927de2eSStefano Zampini 
4523927de2eSStefano Zampini   PetscFunctionBegin;
4533927de2eSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr);
4543927de2eSStefano Zampini   ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr);
4553927de2eSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
4563927de2eSStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
4573927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
4583927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
459ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
460ecf5a873SStefano Zampini   if (A->rmap->mapping != A->cmap->mapping) {
4617230de76SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
462ecf5a873SStefano Zampini   } else {
463ecf5a873SStefano Zampini     global_indices_c = global_indices_r;
464ecf5a873SStefano Zampini   }
465ecf5a873SStefano Zampini 
4663927de2eSStefano Zampini   if (issbaij) {
4673927de2eSStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
4683927de2eSStefano Zampini   }
4693927de2eSStefano Zampini   /*
470ecf5a873SStefano Zampini      An SF reduce is needed to sum up properly on shared rows.
4713927de2eSStefano Zampini      Note that generally preallocation is not exact, since it overestimates nonzeros
4723927de2eSStefano Zampini   */
4733927de2eSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
4743927de2eSStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
4753927de2eSStefano Zampini   }
4763927de2eSStefano Zampini   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
4773927de2eSStefano Zampini   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
4783927de2eSStefano Zampini   /* All processes need to compute entire row ownership */
4793927de2eSStefano Zampini   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
4803927de2eSStefano Zampini   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
4813927de2eSStefano Zampini   for (i=0;i<nsubdomains;i++) {
4823927de2eSStefano Zampini     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
4833927de2eSStefano Zampini       row_ownership[j] = i;
4843927de2eSStefano Zampini     }
4853927de2eSStefano Zampini   }
4867230de76SStefano Zampini   ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
4873927de2eSStefano Zampini 
4883927de2eSStefano Zampini   /*
4893927de2eSStefano Zampini      my_dnz and my_onz contains exact contribution to preallocation from each local mat
4903927de2eSStefano Zampini      then, they will be summed up properly. This way, preallocation is always sufficient
4913927de2eSStefano Zampini   */
4923927de2eSStefano Zampini   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
4933927de2eSStefano Zampini   /* preallocation as a MATAIJ */
4943927de2eSStefano Zampini   if (isdense) { /* special case for dense local matrices */
4953927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
496ecf5a873SStefano Zampini       PetscInt index_row = global_indices_r[i];
4973927de2eSStefano Zampini       for (j=i;j<local_rows;j++) {
4983927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
499ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[j];
5003927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
5013927de2eSStefano Zampini           my_dnz[i] += 1;
5023927de2eSStefano Zampini         } else { /* offdiag block */
5033927de2eSStefano Zampini           my_onz[i] += 1;
5043927de2eSStefano Zampini         }
5053927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
5063927de2eSStefano Zampini         if (i != j) {
5073927de2eSStefano Zampini           owner = row_ownership[index_col];
5083927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
5093927de2eSStefano Zampini             my_dnz[j] += 1;
5103927de2eSStefano Zampini           } else {
5113927de2eSStefano Zampini             my_onz[j] += 1;
5123927de2eSStefano Zampini           }
5133927de2eSStefano Zampini         }
5143927de2eSStefano Zampini       }
5153927de2eSStefano Zampini     }
516*bb1015c3SStefano Zampini   } else if (matis->A->ops->getrowij) {
517*bb1015c3SStefano Zampini     const PetscInt *ii,*jj,*jptr;
518*bb1015c3SStefano Zampini     PetscBool      done;
519*bb1015c3SStefano Zampini     ierr = MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr);
520*bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__);
521*bb1015c3SStefano Zampini     jptr = jj;
522*bb1015c3SStefano Zampini     for (i=0;i<local_rows;i++) {
523*bb1015c3SStefano Zampini       PetscInt index_row = global_indices_r[i];
524*bb1015c3SStefano Zampini       for (j=0;j<ii[i+1]-ii[i];j++,jptr++) {
525*bb1015c3SStefano Zampini         PetscInt owner = row_ownership[index_row];
526*bb1015c3SStefano Zampini         PetscInt index_col = global_indices_c[*jptr];
527*bb1015c3SStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
528*bb1015c3SStefano Zampini           my_dnz[i] += 1;
529*bb1015c3SStefano Zampini         } else { /* offdiag block */
530*bb1015c3SStefano Zampini           my_onz[i] += 1;
531*bb1015c3SStefano Zampini         }
532*bb1015c3SStefano Zampini         /* same as before, interchanging rows and cols */
533*bb1015c3SStefano Zampini         if (issbaij && index_col != index_row) {
534*bb1015c3SStefano Zampini           owner = row_ownership[index_col];
535*bb1015c3SStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
536*bb1015c3SStefano Zampini             my_dnz[*jptr] += 1;
537*bb1015c3SStefano Zampini           } else {
538*bb1015c3SStefano Zampini             my_onz[*jptr] += 1;
539*bb1015c3SStefano Zampini           }
540*bb1015c3SStefano Zampini         }
541*bb1015c3SStefano Zampini       }
542*bb1015c3SStefano Zampini     }
543*bb1015c3SStefano Zampini     ierr = MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr);
544*bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__);
545*bb1015c3SStefano Zampini   } else { /* loop over rows and use MatGetRow */
5463927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
5473927de2eSStefano Zampini       const PetscInt *cols;
548ecf5a873SStefano Zampini       PetscInt       ncols,index_row = global_indices_r[i];
5493927de2eSStefano Zampini       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
5503927de2eSStefano Zampini       for (j=0;j<ncols;j++) {
5513927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
552ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[cols[j]];
5533927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
5543927de2eSStefano Zampini           my_dnz[i] += 1;
5553927de2eSStefano Zampini         } else { /* offdiag block */
5563927de2eSStefano Zampini           my_onz[i] += 1;
5573927de2eSStefano Zampini         }
5583927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
559d9a9e74cSStefano Zampini         if (issbaij && index_col != index_row) {
5603927de2eSStefano Zampini           owner = row_ownership[index_col];
5613927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
562d9a9e74cSStefano Zampini             my_dnz[cols[j]] += 1;
5633927de2eSStefano Zampini           } else {
564d9a9e74cSStefano Zampini             my_onz[cols[j]] += 1;
5653927de2eSStefano Zampini           }
5663927de2eSStefano Zampini         }
5673927de2eSStefano Zampini       }
5683927de2eSStefano Zampini       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
5693927de2eSStefano Zampini     }
5703927de2eSStefano Zampini   }
571ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
572ecf5a873SStefano Zampini   if (global_indices_c != global_indices_r) {
5737230de76SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
574ecf5a873SStefano Zampini   }
5753927de2eSStefano Zampini   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
576ecf5a873SStefano Zampini 
577ecf5a873SStefano Zampini   /* Reduce my_dnz and my_onz */
5783927de2eSStefano Zampini   if (maxreduce) {
5793927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
5803927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
581*bb1015c3SStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
5823927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
5833927de2eSStefano Zampini   } else {
5843927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
5853927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
586*bb1015c3SStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
5873927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
5883927de2eSStefano Zampini   }
5893927de2eSStefano Zampini   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
5903927de2eSStefano Zampini 
5913927de2eSStefano Zampini   /* Resize preallocation if overestimated */
5923927de2eSStefano Zampini   for (i=0;i<lrows;i++) {
5933927de2eSStefano Zampini     dnz[i] = PetscMin(dnz[i],lcols);
5943927de2eSStefano Zampini     onz[i] = PetscMin(onz[i],cols-lcols);
5953927de2eSStefano Zampini   }
5963927de2eSStefano Zampini   /* set preallocation */
5973927de2eSStefano Zampini   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
5983927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
5993927de2eSStefano Zampini     dnz[i] = dnz[i*bs]/bs;
6003927de2eSStefano Zampini     onz[i] = onz[i*bs]/bs;
6013927de2eSStefano Zampini   }
6023927de2eSStefano Zampini   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
6033927de2eSStefano Zampini   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
6043927de2eSStefano Zampini   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
6053927de2eSStefano Zampini   if (issbaij) {
6063927de2eSStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
6073927de2eSStefano Zampini   }
6083927de2eSStefano Zampini   PetscFunctionReturn(0);
6093927de2eSStefano Zampini }
6103927de2eSStefano Zampini 
6113927de2eSStefano Zampini #undef __FUNCT__
612b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
6137230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
614b7ce53b6SStefano Zampini {
615b7ce53b6SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
616d9a9e74cSStefano Zampini   Mat            local_mat;
617b7ce53b6SStefano Zampini   /* info on mat */
6183cfa4ea4SStefano Zampini   PetscInt       bs,rows,cols,lrows,lcols;
619b7ce53b6SStefano Zampini   PetscInt       local_rows,local_cols;
620686e3a49SStefano Zampini   PetscBool      isdense,issbaij,isseqaij;
6217c03b4e8SStefano Zampini   PetscMPIInt    nsubdomains;
622b7ce53b6SStefano Zampini   /* values insertion */
623b7ce53b6SStefano Zampini   PetscScalar    *array;
624b7ce53b6SStefano Zampini   /* work */
625b7ce53b6SStefano Zampini   PetscErrorCode ierr;
626b7ce53b6SStefano Zampini 
627b7ce53b6SStefano Zampini   PetscFunctionBegin;
628b7ce53b6SStefano Zampini   /* get info from mat */
6297c03b4e8SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
6307c03b4e8SStefano Zampini   if (nsubdomains == 1) {
6317c03b4e8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
6327c03b4e8SStefano Zampini       ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr);
6337c03b4e8SStefano Zampini     } else {
6347c03b4e8SStefano Zampini       ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
6357c03b4e8SStefano Zampini     }
6367c03b4e8SStefano Zampini     PetscFunctionReturn(0);
6377c03b4e8SStefano Zampini   }
638b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
639b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
6403cfa4ea4SStefano Zampini   ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr);
641b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
642b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
643686e3a49SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
644b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
645b7ce53b6SStefano Zampini 
646b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
647b7ce53b6SStefano Zampini     MatType     new_mat_type;
6483927de2eSStefano Zampini     PetscBool   issbaij_red;
649b7ce53b6SStefano Zampini 
650b7ce53b6SStefano Zampini     /* determining new matrix type */
651b2566f29SBarry Smith     ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
652b7ce53b6SStefano Zampini     if (issbaij_red) {
653b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
654b7ce53b6SStefano Zampini     } else {
655b7ce53b6SStefano Zampini       if (bs>1) {
656b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
657b7ce53b6SStefano Zampini       } else {
658b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
659b7ce53b6SStefano Zampini       }
660b7ce53b6SStefano Zampini     }
661b7ce53b6SStefano Zampini 
6623927de2eSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
6633cfa4ea4SStefano Zampini     ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr);
6643927de2eSStefano Zampini     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
6653927de2eSStefano Zampini     ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr);
6663927de2eSStefano Zampini     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
667b7ce53b6SStefano Zampini   } else {
6683cfa4ea4SStefano Zampini     PetscInt mbs,mrows,mcols,mlrows,mlcols;
669b7ce53b6SStefano Zampini     /* some checks */
670b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
671b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
6723cfa4ea4SStefano Zampini     ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr);
6736c4ed002SBarry Smith     if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
6746c4ed002SBarry Smith     if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
6756c4ed002SBarry Smith     if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows);
6766c4ed002SBarry Smith     if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols);
6776c4ed002SBarry Smith     if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
678b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
679b7ce53b6SStefano Zampini   }
680d9a9e74cSStefano Zampini 
681d9a9e74cSStefano Zampini   if (issbaij) {
682d9a9e74cSStefano Zampini     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr);
683d9a9e74cSStefano Zampini   } else {
684d9a9e74cSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
685d9a9e74cSStefano Zampini     local_mat = matis->A;
686d9a9e74cSStefano Zampini   }
687686e3a49SStefano Zampini 
688b7ce53b6SStefano Zampini   /* Set values */
689ecf5a873SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr);
690b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
691ecf5a873SStefano Zampini     PetscInt i,*dummy_rows,*dummy_cols;
692ecf5a873SStefano Zampini 
693ecf5a873SStefano Zampini     if (local_rows != local_cols) {
694ecf5a873SStefano Zampini       ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr);
695ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
696ecf5a873SStefano Zampini       for (i=0;i<local_cols;i++) dummy_cols[i] = i;
697ecf5a873SStefano Zampini     } else {
698ecf5a873SStefano Zampini       ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr);
699ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
700ecf5a873SStefano Zampini       dummy_cols = dummy_rows;
701ecf5a873SStefano Zampini     }
702b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
703d9a9e74cSStefano Zampini     ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr);
704ecf5a873SStefano Zampini     ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr);
705d9a9e74cSStefano Zampini     ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr);
706ecf5a873SStefano Zampini     if (dummy_rows != dummy_cols) {
707ecf5a873SStefano Zampini       ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr);
708ecf5a873SStefano Zampini     } else {
709ecf5a873SStefano Zampini       ierr = PetscFree(dummy_rows);CHKERRQ(ierr);
710ecf5a873SStefano Zampini     }
711686e3a49SStefano Zampini   } else if (isseqaij) {
712ecf5a873SStefano Zampini     PetscInt  i,nvtxs,*xadj,*adjncy;
713686e3a49SStefano Zampini     PetscBool done;
714686e3a49SStefano Zampini 
715d9a9e74cSStefano Zampini     ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
716*bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__);
717d9a9e74cSStefano Zampini     ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr);
718686e3a49SStefano Zampini     for (i=0;i<nvtxs;i++) {
719ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr);
720686e3a49SStefano Zampini     }
721d9a9e74cSStefano Zampini     ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
722*bb1015c3SStefano Zampini     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__);
723d9a9e74cSStefano Zampini     ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr);
724686e3a49SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
725ecf5a873SStefano Zampini     PetscInt i;
726c0962df8SStefano Zampini 
727686e3a49SStefano Zampini     for (i=0;i<local_rows;i++) {
728686e3a49SStefano Zampini       PetscInt       j;
729ecf5a873SStefano Zampini       const PetscInt *local_indices_cols;
730686e3a49SStefano Zampini 
731ecf5a873SStefano Zampini       ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
732ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr);
733ecf5a873SStefano Zampini       ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
734686e3a49SStefano Zampini     }
735b7ce53b6SStefano Zampini   }
736b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
737d9a9e74cSStefano Zampini   ierr = MatDestroy(&local_mat);CHKERRQ(ierr);
738b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
739b7ce53b6SStefano Zampini   if (isdense) {
740b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
741b7ce53b6SStefano Zampini   }
742b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
743b7ce53b6SStefano Zampini }
744b7ce53b6SStefano Zampini 
745b7ce53b6SStefano Zampini #undef __FUNCT__
746b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
747b7ce53b6SStefano Zampini /*@
748b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
749b7ce53b6SStefano Zampini 
750b7ce53b6SStefano Zampini   Input Parameter:
751b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
752b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
753b7ce53b6SStefano Zampini 
754b7ce53b6SStefano Zampini   Output Parameter:
755b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
756b7ce53b6SStefano Zampini 
757b7ce53b6SStefano Zampini   Level: developer
758b7ce53b6SStefano Zampini 
759eb82efa4SStefano Zampini   Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested.
760b7ce53b6SStefano Zampini 
761b7ce53b6SStefano Zampini .seealso: MATIS
762b7ce53b6SStefano Zampini @*/
763b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
764b7ce53b6SStefano Zampini {
765b7ce53b6SStefano Zampini   PetscErrorCode ierr;
766b7ce53b6SStefano Zampini 
767b7ce53b6SStefano Zampini   PetscFunctionBegin;
768b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
769b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
770b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
771b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
772b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
773b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
7746c4ed002SBarry Smith     if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix");
775b7ce53b6SStefano Zampini   }
776b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
777b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
778b7ce53b6SStefano Zampini }
779b7ce53b6SStefano Zampini 
780b7ce53b6SStefano Zampini #undef __FUNCT__
781ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
782ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
783ad6194a2SStefano Zampini {
784ad6194a2SStefano Zampini   PetscErrorCode ierr;
785ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
786ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
787ad6194a2SStefano Zampini   Mat            B,localmat;
788ad6194a2SStefano Zampini 
789ad6194a2SStefano Zampini   PetscFunctionBegin;
790ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
791ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
792ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
793e176bc59SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr);
794ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
795ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
796b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
797ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
798ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
799ad6194a2SStefano Zampini   *newmat = B;
800ad6194a2SStefano Zampini   PetscFunctionReturn(0);
801ad6194a2SStefano Zampini }
802ad6194a2SStefano Zampini 
803ad6194a2SStefano Zampini #undef __FUNCT__
80469796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
805a8116848SStefano Zampini static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
80669796d55SStefano Zampini {
80769796d55SStefano Zampini   PetscErrorCode ierr;
80869796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
80969796d55SStefano Zampini   PetscBool      local_sym;
81069796d55SStefano Zampini 
81169796d55SStefano Zampini   PetscFunctionBegin;
81269796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
813b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
81469796d55SStefano Zampini   PetscFunctionReturn(0);
81569796d55SStefano Zampini }
81669796d55SStefano Zampini 
81769796d55SStefano Zampini #undef __FUNCT__
81869796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
819a8116848SStefano Zampini static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
82069796d55SStefano Zampini {
82169796d55SStefano Zampini   PetscErrorCode ierr;
82269796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
82369796d55SStefano Zampini   PetscBool      local_sym;
82469796d55SStefano Zampini 
82569796d55SStefano Zampini   PetscFunctionBegin;
82669796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
827b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
82869796d55SStefano Zampini   PetscFunctionReturn(0);
82969796d55SStefano Zampini }
83069796d55SStefano Zampini 
83169796d55SStefano Zampini #undef __FUNCT__
832b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
833a8116848SStefano Zampini static PetscErrorCode MatDestroy_IS(Mat A)
834b4319ba4SBarry Smith {
835dfbe8321SBarry Smith   PetscErrorCode ierr;
836b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
837b4319ba4SBarry Smith 
838b4319ba4SBarry Smith   PetscFunctionBegin;
8396bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
840e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr);
841e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr);
8426bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
8436bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
844a8116848SStefano Zampini   ierr = ISDestroy(&b->getsub_ris);CHKERRQ(ierr);
845a8116848SStefano Zampini   ierr = ISDestroy(&b->getsub_cis);CHKERRQ(ierr);
846a8116848SStefano Zampini   if (b->sf != b->csf) {
847a8116848SStefano Zampini     ierr = PetscSFDestroy(&b->csf);CHKERRQ(ierr);
848a8116848SStefano Zampini     ierr = PetscFree2(b->csf_rootdata,b->csf_leafdata);CHKERRQ(ierr);
849a8116848SStefano Zampini   }
85028f4e0baSStefano Zampini   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
85128f4e0baSStefano Zampini   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
852bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
853dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
854bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
855b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
856b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
8572e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
858b4319ba4SBarry Smith   PetscFunctionReturn(0);
859b4319ba4SBarry Smith }
860b4319ba4SBarry Smith 
861b4319ba4SBarry Smith #undef __FUNCT__
862b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
863a8116848SStefano Zampini static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
864b4319ba4SBarry Smith {
865dfbe8321SBarry Smith   PetscErrorCode ierr;
866b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
867b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
868b4319ba4SBarry Smith 
869b4319ba4SBarry Smith   PetscFunctionBegin;
870b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
871e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
872e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
873b4319ba4SBarry Smith 
874b4319ba4SBarry Smith   /* multiply the local matrix */
875b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
876b4319ba4SBarry Smith 
877b4319ba4SBarry Smith   /* scatter product back into global memory */
8782dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
879e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
880e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
881b4319ba4SBarry Smith   PetscFunctionReturn(0);
882b4319ba4SBarry Smith }
883b4319ba4SBarry Smith 
884b4319ba4SBarry Smith #undef __FUNCT__
8852e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
886a8116848SStefano Zampini static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
8872e74eeadSLisandro Dalcin {
888650997f4SStefano Zampini   Vec            temp_vec;
8892e74eeadSLisandro Dalcin   PetscErrorCode ierr;
8902e74eeadSLisandro Dalcin 
8912e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
892650997f4SStefano Zampini   if (v3 != v2) {
893650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
894650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
895650997f4SStefano Zampini   } else {
896650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
897650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
898650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
899650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
900650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
901650997f4SStefano Zampini   }
9022e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9032e74eeadSLisandro Dalcin }
9042e74eeadSLisandro Dalcin 
9052e74eeadSLisandro Dalcin #undef __FUNCT__
9062e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
907a8116848SStefano Zampini static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x)
9082e74eeadSLisandro Dalcin {
9092e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9102e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9112e74eeadSLisandro Dalcin 
912e176bc59SStefano Zampini   PetscFunctionBegin;
9132e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
914e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
915e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
9162e74eeadSLisandro Dalcin 
9172e74eeadSLisandro Dalcin   /* multiply the local matrix */
918e176bc59SStefano Zampini   ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr);
9192e74eeadSLisandro Dalcin 
9202e74eeadSLisandro Dalcin   /* scatter product back into global vector */
921e176bc59SStefano Zampini   ierr = VecSet(x,0);CHKERRQ(ierr);
922e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
923e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9242e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9252e74eeadSLisandro Dalcin }
9262e74eeadSLisandro Dalcin 
9272e74eeadSLisandro Dalcin #undef __FUNCT__
9282e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
929a8116848SStefano Zampini static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
9302e74eeadSLisandro Dalcin {
931650997f4SStefano Zampini   Vec            temp_vec;
9322e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9332e74eeadSLisandro Dalcin 
9342e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
935650997f4SStefano Zampini   if (v3 != v2) {
936650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
937650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
938650997f4SStefano Zampini   } else {
939650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
940650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
941650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
942650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
943650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
944650997f4SStefano Zampini   }
9452e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9462e74eeadSLisandro Dalcin }
9472e74eeadSLisandro Dalcin 
9482e74eeadSLisandro Dalcin #undef __FUNCT__
949b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
950a8116848SStefano Zampini static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
951b4319ba4SBarry Smith {
952b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
953dfbe8321SBarry Smith   PetscErrorCode ierr;
954b4319ba4SBarry Smith   PetscViewer    sviewer;
955b4319ba4SBarry Smith 
956b4319ba4SBarry Smith   PetscFunctionBegin;
9573f08860eSBarry Smith   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
958b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
9593f08860eSBarry Smith   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
9606e520ac8SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
961b4319ba4SBarry Smith   PetscFunctionReturn(0);
962b4319ba4SBarry Smith }
963b4319ba4SBarry Smith 
964b4319ba4SBarry Smith #undef __FUNCT__
965b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
966a8116848SStefano Zampini static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
967b4319ba4SBarry Smith {
968dfbe8321SBarry Smith   PetscErrorCode ierr;
969e176bc59SStefano Zampini   PetscInt       nr,rbs,nc,cbs;
970b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
971b4319ba4SBarry Smith   IS             from,to;
972e176bc59SStefano Zampini   Vec            cglobal,rglobal;
973b4319ba4SBarry Smith 
974b4319ba4SBarry Smith   PetscFunctionBegin;
975784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
976e176bc59SStefano Zampini   PetscCheckSameComm(A,1,cmapping,3);
9773bbff08aSStefano Zampini   /* Destroy any previous data */
97870cf5478SStefano Zampini   ierr = VecDestroy(&is->x);CHKERRQ(ierr);
97970cf5478SStefano Zampini   ierr = VecDestroy(&is->y);CHKERRQ(ierr);
980e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr);
981e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr);
9821c47cb0fSStefano Zampini   ierr = MatDestroy(&is->A);CHKERRQ(ierr);
98328f4e0baSStefano Zampini   ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
98428f4e0baSStefano Zampini   ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
9853bbff08aSStefano Zampini 
9863bbff08aSStefano Zampini   /* Setup Layout and set local to global maps */
987fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
988fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
989fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
990fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
991b4319ba4SBarry Smith 
992b4319ba4SBarry Smith   /* Create the local matrix A */
993e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr);
994e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr);
995e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr);
996e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr);
997f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
998e176bc59SStefano Zampini   ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr);
999e176bc59SStefano Zampini   ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr);
1000e176bc59SStefano Zampini   ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr);
1001ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
1002ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
1003b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
1004b4319ba4SBarry Smith 
1005f26d0771SStefano Zampini   if (!is->islocalref) { /* setup scatters and local vectors for MatMult */
1006b4319ba4SBarry Smith     /* Create the local work vectors */
10073bbff08aSStefano Zampini     ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr);
1008b4319ba4SBarry Smith 
1009e176bc59SStefano Zampini     /* setup the global to local scatters */
1010e176bc59SStefano Zampini     ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr);
1011e176bc59SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr);
1012784ac674SJed Brown     ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
1013e176bc59SStefano Zampini     ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr);
1014e176bc59SStefano Zampini     if (rmapping != cmapping) {
1015e176bc59SStefano Zampini       ierr = ISDestroy(&to);CHKERRQ(ierr);
1016e176bc59SStefano Zampini       ierr = ISDestroy(&from);CHKERRQ(ierr);
1017e176bc59SStefano Zampini       ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr);
1018e176bc59SStefano Zampini       ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr);
1019e176bc59SStefano Zampini       ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr);
1020e176bc59SStefano Zampini     } else {
1021e176bc59SStefano Zampini       ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr);
1022e176bc59SStefano Zampini       is->cctx = is->rctx;
1023e176bc59SStefano Zampini     }
1024e176bc59SStefano Zampini     ierr = VecDestroy(&rglobal);CHKERRQ(ierr);
1025e176bc59SStefano Zampini     ierr = VecDestroy(&cglobal);CHKERRQ(ierr);
10266bf464f9SBarry Smith     ierr = ISDestroy(&to);CHKERRQ(ierr);
10276bf464f9SBarry Smith     ierr = ISDestroy(&from);CHKERRQ(ierr);
1028f26d0771SStefano Zampini   }
102948ff6bf3SStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
1030b4319ba4SBarry Smith   PetscFunctionReturn(0);
1031b4319ba4SBarry Smith }
1032b4319ba4SBarry Smith 
10332e74eeadSLisandro Dalcin #undef __FUNCT__
10342e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
1035a8116848SStefano Zampini static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
10362e74eeadSLisandro Dalcin {
10372e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
10382e74eeadSLisandro Dalcin   PetscErrorCode ierr;
103997563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
104097563a80SStefano Zampini   PetscInt       i,zm,zn;
104197563a80SStefano Zampini #endif
1042f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
10432e74eeadSLisandro Dalcin 
10442e74eeadSLisandro Dalcin   PetscFunctionBegin;
10452e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
1046f26d0771SStefano 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);
104797563a80SStefano Zampini   /* count negative indices */
104897563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
104997563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
10502e74eeadSLisandro Dalcin #endif
105197563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
105297563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
105397563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
105497563a80SStefano Zampini   /* count negative indices (should be the same as before) */
105597563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
105697563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
105797563a80SStefano 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");
105897563a80SStefano 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");
105997563a80SStefano Zampini #endif
10602e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
10612e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
10622e74eeadSLisandro Dalcin }
10632e74eeadSLisandro Dalcin 
1064b4319ba4SBarry Smith #undef __FUNCT__
106597563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS"
1066a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
106797563a80SStefano Zampini {
106897563a80SStefano Zampini   Mat_IS         *is = (Mat_IS*)mat->data;
106997563a80SStefano Zampini   PetscErrorCode ierr;
107097563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
107197563a80SStefano Zampini   PetscInt       i,zm,zn;
107297563a80SStefano Zampini #endif
1073f26d0771SStefano Zampini   PetscInt       rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION];
107497563a80SStefano Zampini 
107597563a80SStefano Zampini   PetscFunctionBegin;
107697563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
1077f26d0771SStefano 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);
107897563a80SStefano Zampini   /* count negative indices */
107997563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
108097563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
108197563a80SStefano Zampini #endif
108297563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
108397563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
108497563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
108597563a80SStefano Zampini   /* count negative indices (should be the same as before) */
108697563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
108797563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
108897563a80SStefano 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");
108997563a80SStefano 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");
109097563a80SStefano Zampini #endif
109197563a80SStefano Zampini   ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
109297563a80SStefano Zampini   PetscFunctionReturn(0);
109397563a80SStefano Zampini }
109497563a80SStefano Zampini 
109597563a80SStefano Zampini #undef __FUNCT__
1096b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
1097a8116848SStefano Zampini static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
1098b4319ba4SBarry Smith {
1099dfbe8321SBarry Smith   PetscErrorCode ierr;
1100b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1101b4319ba4SBarry Smith 
1102b4319ba4SBarry Smith   PetscFunctionBegin;
1103b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
1104b4319ba4SBarry Smith   PetscFunctionReturn(0);
1105b4319ba4SBarry Smith }
1106b4319ba4SBarry Smith 
1107b4319ba4SBarry Smith #undef __FUNCT__
1108f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
1109a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
1110f0006bf2SLisandro Dalcin {
1111f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
1112f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
1113f0006bf2SLisandro Dalcin 
1114f0006bf2SLisandro Dalcin   PetscFunctionBegin;
1115f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
1116f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
1117f0006bf2SLisandro Dalcin }
1118f0006bf2SLisandro Dalcin 
1119f0006bf2SLisandro Dalcin #undef __FUNCT__
11202e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
1121a8116848SStefano Zampini static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
11222e74eeadSLisandro Dalcin {
11236e520ac8SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
11246e520ac8SStefano Zampini   PetscInt       nr,nl,len,i;
11256e520ac8SStefano Zampini   PetscInt       *lrows;
11262e74eeadSLisandro Dalcin   PetscErrorCode ierr;
11272e74eeadSLisandro Dalcin 
11282e74eeadSLisandro Dalcin   PetscFunctionBegin;
11296e520ac8SStefano Zampini   /* get locally owned rows */
11306e520ac8SStefano Zampini   ierr = MatZeroRowsMapLocal_Private(A,n,rows,&len,&lrows);CHKERRQ(ierr);
11316e520ac8SStefano Zampini   /* fix right hand side if needed */
11326e520ac8SStefano Zampini   if (x && b) {
11336e520ac8SStefano Zampini     const PetscScalar *xx;
11346e520ac8SStefano Zampini     PetscScalar       *bb;
11356e520ac8SStefano Zampini 
11366e520ac8SStefano Zampini     ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr);
11376e520ac8SStefano Zampini     ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
11386e520ac8SStefano Zampini     for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]];
11396e520ac8SStefano Zampini     ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr);
11406e520ac8SStefano Zampini     ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr);
11412e74eeadSLisandro Dalcin   }
11426e520ac8SStefano Zampini   /* get rows associated to the local matrices */
11436e520ac8SStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
11446e520ac8SStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
11456e520ac8SStefano Zampini   }
11466e520ac8SStefano Zampini   ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr);
11476e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr);
11486e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
11496e520ac8SStefano Zampini   for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1;
11506e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
11516e520ac8SStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
11526e520ac8SStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
11536e520ac8SStefano Zampini   ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr);
11546e520ac8SStefano Zampini   for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i;
11557230de76SStefano Zampini   ierr = MatISZeroRowsLocal_Private(A,nr,lrows,diag);CHKERRQ(ierr);
11566e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
11572e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
11582e74eeadSLisandro Dalcin }
11592e74eeadSLisandro Dalcin 
11602e74eeadSLisandro Dalcin #undef __FUNCT__
11617230de76SStefano Zampini #define __FUNCT__ "MatISZeroRowsLocal_Private"
11627230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag)
1163b4319ba4SBarry Smith {
1164b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1165dfbe8321SBarry Smith   PetscErrorCode ierr;
1166b4319ba4SBarry Smith 
1167b4319ba4SBarry Smith   PetscFunctionBegin;
11686e520ac8SStefano Zampini   if (diag != 0.) {
1169b4319ba4SBarry Smith     /*
1170b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
1171b4319ba4SBarry Smith        work properly in the interface nodes.
1172b4319ba4SBarry Smith     */
1173b4319ba4SBarry Smith     Vec counter;
1174e176bc59SStefano Zampini     ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr);
1175e176bc59SStefano Zampini     ierr = VecSet(counter,0.);CHKERRQ(ierr);
1176e176bc59SStefano Zampini     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
1177e176bc59SStefano Zampini     ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1178e176bc59SStefano Zampini     ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1179e176bc59SStefano Zampini     ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1180e176bc59SStefano Zampini     ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
11816bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
1182b4319ba4SBarry Smith   }
11832b404112SStefano Zampini 
1184958c9bccSBarry Smith   if (!n) {
1185b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
1186b4319ba4SBarry Smith   } else {
11877230de76SStefano Zampini     PetscInt i;
1188b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
11892205254eSKarl Rupp 
11902b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
11916e520ac8SStefano Zampini     if (diag != 0.) {
11926e520ac8SStefano Zampini       const PetscScalar *array;
11936e520ac8SStefano Zampini       ierr = VecGetArrayRead(is->y,&array);CHKERRQ(ierr);
1194b4319ba4SBarry Smith       for (i=0; i<n; i++) {
1195f4df32b1SMatthew Knepley         ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
1196b4319ba4SBarry Smith       }
11976e520ac8SStefano Zampini       ierr = VecRestoreArrayRead(is->y,&array);CHKERRQ(ierr);
11986e520ac8SStefano Zampini     }
1199b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1200b4319ba4SBarry Smith     ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1201b4319ba4SBarry Smith   }
1202b4319ba4SBarry Smith   PetscFunctionReturn(0);
1203b4319ba4SBarry Smith }
1204b4319ba4SBarry Smith 
1205b4319ba4SBarry Smith #undef __FUNCT__
1206b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
1207a8116848SStefano Zampini static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
1208b4319ba4SBarry Smith {
1209b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1210dfbe8321SBarry Smith   PetscErrorCode ierr;
1211b4319ba4SBarry Smith 
1212b4319ba4SBarry Smith   PetscFunctionBegin;
1213b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
1214b4319ba4SBarry Smith   PetscFunctionReturn(0);
1215b4319ba4SBarry Smith }
1216b4319ba4SBarry Smith 
1217b4319ba4SBarry Smith #undef __FUNCT__
1218b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
1219a8116848SStefano Zampini static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
1220b4319ba4SBarry Smith {
1221b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
1222dfbe8321SBarry Smith   PetscErrorCode ierr;
1223b4319ba4SBarry Smith 
1224b4319ba4SBarry Smith   PetscFunctionBegin;
1225b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
1226b4319ba4SBarry Smith   PetscFunctionReturn(0);
1227b4319ba4SBarry Smith }
1228b4319ba4SBarry Smith 
1229b4319ba4SBarry Smith #undef __FUNCT__
1230b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
1231a8116848SStefano Zampini static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
1232b4319ba4SBarry Smith {
1233b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
1234b4319ba4SBarry Smith 
1235b4319ba4SBarry Smith   PetscFunctionBegin;
1236b4319ba4SBarry Smith   *local = is->A;
1237b4319ba4SBarry Smith   PetscFunctionReturn(0);
1238b4319ba4SBarry Smith }
1239b4319ba4SBarry Smith 
1240b4319ba4SBarry Smith #undef __FUNCT__
1241b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
1242b4319ba4SBarry Smith /*@
1243b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
1244b4319ba4SBarry Smith 
1245b4319ba4SBarry Smith   Input Parameter:
1246b4319ba4SBarry Smith .  mat - the matrix
1247b4319ba4SBarry Smith 
1248b4319ba4SBarry Smith   Output Parameter:
1249eb82efa4SStefano Zampini .  local - the local matrix
1250b4319ba4SBarry Smith 
1251b4319ba4SBarry Smith   Level: advanced
1252b4319ba4SBarry Smith 
1253b4319ba4SBarry Smith   Notes:
1254b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
1255b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
1256b4319ba4SBarry Smith   of the MatSetValues() operation.
1257b4319ba4SBarry Smith 
1258b4319ba4SBarry Smith .seealso: MATIS
1259b4319ba4SBarry Smith @*/
12607087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
1261b4319ba4SBarry Smith {
12624ac538c5SBarry Smith   PetscErrorCode ierr;
1263b4319ba4SBarry Smith 
1264b4319ba4SBarry Smith   PetscFunctionBegin;
12650700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1266b4319ba4SBarry Smith   PetscValidPointer(local,2);
12674ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
1268b4319ba4SBarry Smith   PetscFunctionReturn(0);
1269b4319ba4SBarry Smith }
1270b4319ba4SBarry Smith 
12713b03a366Sstefano_zampini #undef __FUNCT__
12723b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
1273a8116848SStefano Zampini static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
12743b03a366Sstefano_zampini {
12753b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
12763b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
12773b03a366Sstefano_zampini   PetscErrorCode ierr;
12783b03a366Sstefano_zampini 
12793b03a366Sstefano_zampini   PetscFunctionBegin;
12804e4c7dbeSStefano Zampini   if (is->A) {
12813b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
12823b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
1283f23aa3ddSBarry Smith     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)\n",orows,ocols,nrows,ncols);
12844e4c7dbeSStefano Zampini   }
12853b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
12863b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
12873b03a366Sstefano_zampini   is->A = local;
12883b03a366Sstefano_zampini   PetscFunctionReturn(0);
12893b03a366Sstefano_zampini }
12903b03a366Sstefano_zampini 
12913b03a366Sstefano_zampini #undef __FUNCT__
12923b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
12933b03a366Sstefano_zampini /*@
1294eb82efa4SStefano Zampini     MatISSetLocalMat - Replace the local matrix stored inside a MATIS object.
12953b03a366Sstefano_zampini 
12963b03a366Sstefano_zampini   Input Parameter:
12973b03a366Sstefano_zampini .  mat - the matrix
1298eb82efa4SStefano Zampini .  local - the local matrix
12993b03a366Sstefano_zampini 
13003b03a366Sstefano_zampini   Output Parameter:
13013b03a366Sstefano_zampini 
13023b03a366Sstefano_zampini   Level: advanced
13033b03a366Sstefano_zampini 
13043b03a366Sstefano_zampini   Notes:
13053b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
13063b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
13073b03a366Sstefano_zampini 
13083b03a366Sstefano_zampini .seealso: MATIS
13093b03a366Sstefano_zampini @*/
13103b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
13113b03a366Sstefano_zampini {
13123b03a366Sstefano_zampini   PetscErrorCode ierr;
13133b03a366Sstefano_zampini 
13143b03a366Sstefano_zampini   PetscFunctionBegin;
13153b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
1316b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
13173b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
13183b03a366Sstefano_zampini   PetscFunctionReturn(0);
13193b03a366Sstefano_zampini }
13203b03a366Sstefano_zampini 
13216726f965SBarry Smith #undef __FUNCT__
13226726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
1323a8116848SStefano Zampini static PetscErrorCode MatZeroEntries_IS(Mat A)
13246726f965SBarry Smith {
13256726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
13266726f965SBarry Smith   PetscErrorCode ierr;
13276726f965SBarry Smith 
13286726f965SBarry Smith   PetscFunctionBegin;
13296726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
13306726f965SBarry Smith   PetscFunctionReturn(0);
13316726f965SBarry Smith }
13326726f965SBarry Smith 
13336726f965SBarry Smith #undef __FUNCT__
13342e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
1335a8116848SStefano Zampini static PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
13362e74eeadSLisandro Dalcin {
13372e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
13382e74eeadSLisandro Dalcin   PetscErrorCode ierr;
13392e74eeadSLisandro Dalcin 
13402e74eeadSLisandro Dalcin   PetscFunctionBegin;
13412e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
13422e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
13432e74eeadSLisandro Dalcin }
13442e74eeadSLisandro Dalcin 
13452e74eeadSLisandro Dalcin #undef __FUNCT__
13462e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
1347a8116848SStefano Zampini static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
13482e74eeadSLisandro Dalcin {
13492e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
13502e74eeadSLisandro Dalcin   PetscErrorCode ierr;
13512e74eeadSLisandro Dalcin 
13522e74eeadSLisandro Dalcin   PetscFunctionBegin;
13532e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
1354e176bc59SStefano Zampini   ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr);
13552e74eeadSLisandro Dalcin 
13562e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
13572e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
1358e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1359e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
13602e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
13612e74eeadSLisandro Dalcin }
13622e74eeadSLisandro Dalcin 
13632e74eeadSLisandro Dalcin #undef __FUNCT__
13646726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
1365a8116848SStefano Zampini static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
13666726f965SBarry Smith {
13676726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
13686726f965SBarry Smith   PetscErrorCode ierr;
13696726f965SBarry Smith 
13706726f965SBarry Smith   PetscFunctionBegin;
13714e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
13726726f965SBarry Smith   PetscFunctionReturn(0);
13736726f965SBarry Smith }
13746726f965SBarry Smith 
1375284134d9SBarry Smith #undef __FUNCT__
1376f26d0771SStefano Zampini #define __FUNCT__ "MatAXPY_IS"
1377f26d0771SStefano Zampini static PetscErrorCode MatAXPY_IS(Mat Y,PetscScalar a,Mat X,MatStructure str)
1378f26d0771SStefano Zampini {
1379f26d0771SStefano Zampini   Mat_IS         *y = (Mat_IS*)Y->data;
1380f26d0771SStefano Zampini   Mat_IS         *x;
1381f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1382f26d0771SStefano Zampini   PetscBool      ismatis;
1383f26d0771SStefano Zampini #endif
1384f26d0771SStefano Zampini   PetscErrorCode ierr;
1385f26d0771SStefano Zampini 
1386f26d0771SStefano Zampini   PetscFunctionBegin;
1387f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1388f26d0771SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)X,MATIS,&ismatis);CHKERRQ(ierr);
1389f26d0771SStefano Zampini   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS");
1390f26d0771SStefano Zampini #endif
1391f26d0771SStefano Zampini   x = (Mat_IS*)X->data;
1392f26d0771SStefano Zampini   ierr = MatAXPY(y->A,a,x->A,str);CHKERRQ(ierr);
1393f26d0771SStefano Zampini   PetscFunctionReturn(0);
1394f26d0771SStefano Zampini }
1395f26d0771SStefano Zampini 
1396f26d0771SStefano Zampini #undef __FUNCT__
1397f26d0771SStefano Zampini #define __FUNCT__ "MatGetLocalSubMatrix_IS"
1398f26d0771SStefano Zampini static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat)
1399f26d0771SStefano Zampini {
1400f26d0771SStefano Zampini   Mat                    lA;
1401f26d0771SStefano Zampini   Mat_IS                 *matis;
1402f26d0771SStefano Zampini   ISLocalToGlobalMapping rl2g,cl2g;
1403f26d0771SStefano Zampini   IS                     is;
1404f26d0771SStefano Zampini   const PetscInt         *rg,*rl;
1405f26d0771SStefano Zampini   PetscInt               nrg;
1406f26d0771SStefano Zampini   PetscInt               N,M,nrl,i,*idxs;
1407f26d0771SStefano Zampini   PetscErrorCode         ierr;
1408f26d0771SStefano Zampini 
1409f26d0771SStefano Zampini   PetscFunctionBegin;
1410f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr);
1411f26d0771SStefano Zampini   ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr);
1412f26d0771SStefano Zampini   ierr = ISGetIndices(row,&rl);CHKERRQ(ierr);
1413f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr);
1414f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1415f26d0771SStefano 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\n",i,rl[i],nrg);
1416f26d0771SStefano Zampini #endif
1417f26d0771SStefano Zampini   ierr = PetscMalloc1(nrg,&idxs);CHKERRQ(ierr);
1418f26d0771SStefano Zampini   /* map from [0,nrl) to row */
1419f26d0771SStefano Zampini   for (i=0;i<nrl;i++) idxs[i] = rl[i];
1420f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1421f26d0771SStefano Zampini   for (i=nrl;i<nrg;i++) idxs[i] = nrg;
1422f26d0771SStefano Zampini #else
1423f26d0771SStefano Zampini   for (i=nrl;i<nrg;i++) idxs[i] = -1;
1424f26d0771SStefano Zampini #endif
1425f26d0771SStefano Zampini   ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr);
1426f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr);
1427f26d0771SStefano Zampini   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
1428f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
1429f26d0771SStefano Zampini   ierr = ISDestroy(&is);CHKERRQ(ierr);
1430f26d0771SStefano Zampini   /* compute new l2g map for columns */
1431f26d0771SStefano Zampini   if (col != row || A->rmap->mapping != A->cmap->mapping) {
1432f26d0771SStefano Zampini     const PetscInt *cg,*cl;
1433f26d0771SStefano Zampini     PetscInt       ncg;
1434f26d0771SStefano Zampini     PetscInt       ncl;
1435f26d0771SStefano Zampini 
1436f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr);
1437f26d0771SStefano Zampini     ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr);
1438f26d0771SStefano Zampini     ierr = ISGetIndices(col,&cl);CHKERRQ(ierr);
1439f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr);
1440f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1441f26d0771SStefano 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\n",i,cl[i],ncg);
1442f26d0771SStefano Zampini #endif
1443f26d0771SStefano Zampini     ierr = PetscMalloc1(ncg,&idxs);CHKERRQ(ierr);
1444f26d0771SStefano Zampini     /* map from [0,ncl) to col */
1445f26d0771SStefano Zampini     for (i=0;i<ncl;i++) idxs[i] = cl[i];
1446f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG)
1447f26d0771SStefano Zampini     for (i=ncl;i<ncg;i++) idxs[i] = ncg;
1448f26d0771SStefano Zampini #else
1449f26d0771SStefano Zampini     for (i=ncl;i<ncg;i++) idxs[i] = -1;
1450f26d0771SStefano Zampini #endif
1451f26d0771SStefano Zampini     ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr);
1452f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr);
1453f26d0771SStefano Zampini     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
1454f26d0771SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
1455f26d0771SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
1456f26d0771SStefano Zampini   } else {
1457f26d0771SStefano Zampini     ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr);
1458f26d0771SStefano Zampini     cl2g = rl2g;
1459f26d0771SStefano Zampini   }
1460f26d0771SStefano Zampini   /* create the MATIS submatrix */
1461f26d0771SStefano Zampini   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1462f26d0771SStefano Zampini   ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr);
1463f26d0771SStefano Zampini   ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
1464f26d0771SStefano Zampini   ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr);
1465b0aa3428SStefano Zampini   matis = (Mat_IS*)((*submat)->data);
1466f26d0771SStefano Zampini   matis->islocalref = PETSC_TRUE;
1467f26d0771SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr);
1468f26d0771SStefano Zampini   ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr);
1469f26d0771SStefano Zampini   ierr = MatISSetLocalMat(*submat,lA);CHKERRQ(ierr);
1470f26d0771SStefano Zampini   ierr = MatSetUp(*submat);CHKERRQ(ierr);
1471f26d0771SStefano Zampini   ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1472f26d0771SStefano Zampini   ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1473f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
1474f26d0771SStefano Zampini   ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
1475f26d0771SStefano Zampini   /* remove unsupported ops */
1476f26d0771SStefano Zampini   ierr = PetscMemzero((*submat)->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1477f26d0771SStefano Zampini   (*submat)->ops->destroy               = MatDestroy_IS;
1478f26d0771SStefano Zampini   (*submat)->ops->setvalueslocal        = MatSetValuesLocal_SubMat_IS;
1479f26d0771SStefano Zampini   (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS;
1480f26d0771SStefano Zampini   (*submat)->ops->assemblybegin         = MatAssemblyBegin_IS;
1481f26d0771SStefano Zampini   (*submat)->ops->assemblyend           = MatAssemblyEnd_IS;
1482f26d0771SStefano Zampini   PetscFunctionReturn(0);
1483f26d0771SStefano Zampini }
1484f26d0771SStefano Zampini 
1485f26d0771SStefano Zampini #undef __FUNCT__
1486284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
1487284134d9SBarry Smith /*@
14883c212e90SHong Zhang     MatCreateIS - Creates a "process" unassembled matrix, assembled on each
1489284134d9SBarry Smith        process but not across processes.
1490284134d9SBarry Smith 
1491284134d9SBarry Smith    Input Parameters:
1492284134d9SBarry Smith +     comm    - MPI communicator that will share the matrix
1493e176bc59SStefano Zampini .     bs      - block size of the matrix
1494df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
1495e176bc59SStefano Zampini .     rmap    - local to global map for rows
1496e176bc59SStefano Zampini -     cmap    - local to global map for cols
1497284134d9SBarry Smith 
1498284134d9SBarry Smith    Output Parameter:
1499284134d9SBarry Smith .    A - the resulting matrix
1500284134d9SBarry Smith 
15018e6c10adSSatish Balay    Level: advanced
15028e6c10adSSatish Balay 
15033c212e90SHong Zhang    Notes: See MATIS for more details.
15043c212e90SHong Zhang           m and n are NOT related to the size of the map; they are the size of the part of the vector owned
1505e176bc59SStefano Zampini           by that process. The sizes of rmap and cmap define the size of the local matrices.
15063c212e90SHong Zhang           If either rmap or cmap are NULL, then the matrix is assumed to be square.
1507284134d9SBarry Smith 
1508284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
1509284134d9SBarry Smith @*/
1510e176bc59SStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A)
1511284134d9SBarry Smith {
1512284134d9SBarry Smith   PetscErrorCode ierr;
1513284134d9SBarry Smith 
1514284134d9SBarry Smith   PetscFunctionBegin;
1515e176bc59SStefano Zampini   if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping");
1516284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1517d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
1518284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1519284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
1520e176bc59SStefano Zampini   if (rmap && cmap) {
1521e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr);
1522e176bc59SStefano Zampini   } else if (!rmap) {
1523e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr);
1524e176bc59SStefano Zampini   } else {
1525e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr);
1526e176bc59SStefano Zampini   }
1527284134d9SBarry Smith   PetscFunctionReturn(0);
1528284134d9SBarry Smith }
1529284134d9SBarry Smith 
1530b4319ba4SBarry Smith /*MC
1531f26d0771SStefano Zampini    MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition methods (e.g. PCBDDC or KSPFETIDP).
1532b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
1533b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
1534b4319ba4SBarry Smith    product is handled "implicitly".
1535b4319ba4SBarry Smith 
1536b4319ba4SBarry Smith    Operations Provided:
15376726f965SBarry Smith +  MatMult()
15382e74eeadSLisandro Dalcin .  MatMultAdd()
15392e74eeadSLisandro Dalcin .  MatMultTranspose()
15402e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
15416726f965SBarry Smith .  MatZeroEntries()
15426726f965SBarry Smith .  MatSetOption()
15432e74eeadSLisandro Dalcin .  MatZeroRows()
15442e74eeadSLisandro Dalcin .  MatSetValues()
154597563a80SStefano Zampini .  MatSetValuesBlocked()
15466726f965SBarry Smith .  MatSetValuesLocal()
154797563a80SStefano Zampini .  MatSetValuesBlockedLocal()
15482e74eeadSLisandro Dalcin .  MatScale()
15492e74eeadSLisandro Dalcin .  MatGetDiagonal()
15502b404112SStefano Zampini .  MatMissingDiagonal()
15512b404112SStefano Zampini .  MatDuplicate()
15522b404112SStefano Zampini .  MatCopy()
1553f26d0771SStefano Zampini .  MatAXPY()
1554f26d0771SStefano Zampini .  MatGetSubMatrix()
1555f26d0771SStefano Zampini .  MatGetLocalSubMatrix()
15566726f965SBarry Smith -  MatSetLocalToGlobalMapping()
1557b4319ba4SBarry Smith 
1558b4319ba4SBarry Smith    Options Database Keys:
1559b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
1560b4319ba4SBarry Smith 
1561b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
1562b4319ba4SBarry Smith 
1563b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
1564b4319ba4SBarry Smith 
1565b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
1566eb82efa4SStefano Zampini           MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation()
1567b4319ba4SBarry Smith 
1568b4319ba4SBarry Smith   Level: advanced
1569b4319ba4SBarry Smith 
1570f26d0771SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), MatCreateIS(), PCBDDC, KSPFETIDP
1571b4319ba4SBarry Smith 
1572b4319ba4SBarry Smith M*/
1573b4319ba4SBarry Smith 
1574b4319ba4SBarry Smith #undef __FUNCT__
1575b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
15768cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
1577b4319ba4SBarry Smith {
1578dfbe8321SBarry Smith   PetscErrorCode ierr;
1579b4319ba4SBarry Smith   Mat_IS         *b;
1580b4319ba4SBarry Smith 
1581b4319ba4SBarry Smith   PetscFunctionBegin;
1582b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1583b4319ba4SBarry Smith   A->data = (void*)b;
1584b4319ba4SBarry Smith 
1585e176bc59SStefano Zampini   /* matrix ops */
1586e176bc59SStefano Zampini   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1587b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
15882e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
15892e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
15902e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1591b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1592b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
15932e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
159498921651SStefano Zampini   A->ops->setvaluesblocked        = MatSetValuesBlocked_IS;
1595b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1596f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
15972e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1598b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1599b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1600b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
16016726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
16022e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
16032e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
16046726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
160569796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
160669796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1607ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
16086bd84002SStefano Zampini   A->ops->missingdiagonal         = MatMissingDiagonal_IS;
16092b404112SStefano Zampini   A->ops->copy                    = MatCopy_IS;
1610659959c5SStefano Zampini   A->ops->getlocalsubmatrix       = MatGetLocalSubMatrix_IS;
1611a8116848SStefano Zampini   A->ops->getsubmatrix            = MatGetSubMatrix_IS;
1612f26d0771SStefano Zampini   A->ops->axpy                    = MatAXPY_IS;
1613b4319ba4SBarry Smith 
1614b7ce53b6SStefano Zampini   /* special MATIS functions */
1615bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1616bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1617b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
16182e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
161917667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1620b4319ba4SBarry Smith   PetscFunctionReturn(0);
1621b4319ba4SBarry Smith }
1622