xref: /petsc/src/mat/impls/is/matis.c (revision 7c03b4e8d47c574c21ea6318ab3f8833a74a88bb)
1be1d678aSKris Buschelman 
2b4319ba4SBarry Smith /*
3b4319ba4SBarry Smith     Creates a matrix class for using the Neumann-Neumann type preconditioners.
4b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
5b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
6b4319ba4SBarry Smith    product is handled "implicitly".
7b4319ba4SBarry Smith 
8b4319ba4SBarry Smith      We provide:
9b4319ba4SBarry Smith          MatMult()
10b4319ba4SBarry Smith 
11b4319ba4SBarry Smith     Currently this allows for only one subdomain per processor.
12b4319ba4SBarry Smith 
13b4319ba4SBarry Smith */
14b4319ba4SBarry Smith 
15c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h>      /*I "petscmat.h" I*/
162e1947a5SStefano Zampini #include <petscsf.h>
172e1947a5SStefano Zampini 
182e1947a5SStefano Zampini #undef __FUNCT__
192e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation"
20a88811baSStefano Zampini /*
21a88811baSStefano Zampini    MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix.
22a88811baSStefano Zampini 
23a88811baSStefano Zampini    Collective on MPI_Comm
24a88811baSStefano Zampini 
25a88811baSStefano Zampini    Input Parameters:
26a88811baSStefano Zampini +  B - the matrix
27a88811baSStefano Zampini .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
28a88811baSStefano Zampini            (same value is used for all local rows)
29a88811baSStefano Zampini .  d_nnz - array containing the number of nonzeros in the various rows of the
30a88811baSStefano Zampini            DIAGONAL portion of the local submatrix (possibly different for each row)
31a88811baSStefano Zampini            or NULL, if d_nz is used to specify the nonzero structure.
32a88811baSStefano Zampini            The size of this array is equal to the number of local rows, i.e 'm'.
33a88811baSStefano Zampini            For matrices that will be factored, you must leave room for (and set)
34a88811baSStefano Zampini            the diagonal entry even if it is zero.
35a88811baSStefano Zampini .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
36a88811baSStefano Zampini            submatrix (same value is used for all local rows).
37a88811baSStefano Zampini -  o_nnz - array containing the number of nonzeros in the various rows of the
38a88811baSStefano Zampini            OFF-DIAGONAL portion of the local submatrix (possibly different for
39a88811baSStefano Zampini            each row) or NULL, if o_nz is used to specify the nonzero
40a88811baSStefano Zampini            structure. The size of this array is equal to the number
41a88811baSStefano Zampini            of local rows, i.e 'm'.
42a88811baSStefano Zampini 
43a88811baSStefano Zampini    If the *_nnz parameter is given then the *_nz parameter is ignored
44a88811baSStefano Zampini 
45a88811baSStefano Zampini    Level: intermediate
46a88811baSStefano Zampini 
47a88811baSStefano Zampini    Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition
48a88811baSStefano Zampini           from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local
49a88811baSStefano Zampini           matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects.
50a88811baSStefano Zampini 
51a88811baSStefano Zampini .keywords: matrix
52a88811baSStefano Zampini 
53a88811baSStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat()
54a88811baSStefano Zampini @*/
552e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
562e1947a5SStefano Zampini {
572e1947a5SStefano Zampini   PetscErrorCode ierr;
582e1947a5SStefano Zampini 
592e1947a5SStefano Zampini   PetscFunctionBegin;
602e1947a5SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
612e1947a5SStefano Zampini   PetscValidType(B,1);
622e1947a5SStefano Zampini   ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
632e1947a5SStefano Zampini   PetscFunctionReturn(0);
642e1947a5SStefano Zampini }
652e1947a5SStefano Zampini 
662e1947a5SStefano Zampini #undef __FUNCT__
672e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS"
682e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
692e1947a5SStefano Zampini {
702e1947a5SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
712e1947a5SStefano Zampini   PetscSF        sf;
722e1947a5SStefano Zampini   PetscInt       bs,i,nroots,*rootdata,nleaves,*leafdata,nlocalcols;
732e1947a5SStefano Zampini   const PetscInt *gidxs;
742e1947a5SStefano Zampini   PetscErrorCode ierr;
752e1947a5SStefano Zampini 
762e1947a5SStefano Zampini   PetscFunctionBegin;
772e1947a5SStefano Zampini   if (!matis->A) {
782e1947a5SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
792e1947a5SStefano Zampini   }
802e1947a5SStefano Zampini   ierr = MatGetLocalSize(B,&nroots,NULL);CHKERRQ(ierr);
812e1947a5SStefano Zampini   ierr = MatGetSize(matis->A,&nleaves,&nlocalcols);CHKERRQ(ierr);
822e1947a5SStefano Zampini   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
832e1947a5SStefano Zampini   ierr = PetscCalloc2(nroots,&rootdata,nleaves,&leafdata);CHKERRQ(ierr);
842e1947a5SStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&sf);CHKERRQ(ierr);
852e1947a5SStefano Zampini   ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr);
862e1947a5SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&gidxs);CHKERRQ(ierr);
872e1947a5SStefano Zampini   ierr = PetscSFSetGraphLayout(sf,B->rmap,nleaves,NULL,PETSC_COPY_VALUES,gidxs);CHKERRQ(ierr);
882e1947a5SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&gidxs);CHKERRQ(ierr);
892e1947a5SStefano Zampini   if (!d_nnz) {
902e1947a5SStefano Zampini     for (i=0;i<nroots;i++) rootdata[i] += d_nz;
912e1947a5SStefano Zampini   } else {
922e1947a5SStefano Zampini     for (i=0;i<nroots;i++) rootdata[i] += d_nnz[i];
932e1947a5SStefano Zampini   }
942e1947a5SStefano Zampini   if (!o_nnz) {
952e1947a5SStefano Zampini     for (i=0;i<nroots;i++) rootdata[i] += o_nz;
962e1947a5SStefano Zampini   } else {
972e1947a5SStefano Zampini     for (i=0;i<nroots;i++) rootdata[i] += o_nnz[i];
982e1947a5SStefano Zampini   }
992e1947a5SStefano Zampini   ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
1002e1947a5SStefano Zampini   ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr);
1012e1947a5SStefano Zampini   for (i=0;i<nleaves;i++) {
1022e1947a5SStefano Zampini     leafdata[i] = PetscMin(leafdata[i],nlocalcols);
1032e1947a5SStefano Zampini   }
1042e1947a5SStefano Zampini   ierr = MatSeqAIJSetPreallocation(matis->A,0,leafdata);CHKERRQ(ierr);
1052e1947a5SStefano Zampini   for (i=0;i<nleaves/bs;i++) {
1062e1947a5SStefano Zampini     leafdata[i] = leafdata[i*bs]/bs;
1072e1947a5SStefano Zampini   }
1082e1947a5SStefano Zampini   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,leafdata);CHKERRQ(ierr);
1092e1947a5SStefano Zampini   for (i=0;i<nleaves/bs;i++) {
1102e1947a5SStefano Zampini     leafdata[i] = leafdata[i]-i;
1112e1947a5SStefano Zampini   }
1122e1947a5SStefano Zampini   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,leafdata);CHKERRQ(ierr);
1132e1947a5SStefano Zampini   ierr = PetscSFDestroy(&sf);CHKERRQ(ierr);
1142e1947a5SStefano Zampini   ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr);
1152e1947a5SStefano Zampini   PetscFunctionReturn(0);
1162e1947a5SStefano Zampini }
117b4319ba4SBarry Smith 
118b4319ba4SBarry Smith #undef __FUNCT__
119b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
120b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
121b7ce53b6SStefano Zampini {
122b7ce53b6SStefano Zampini   Mat_IS                 *matis = (Mat_IS*)(mat->data);
123b7ce53b6SStefano Zampini   /* info on mat */
124b7ce53b6SStefano Zampini   /* ISLocalToGlobalMapping rmapping,cmapping; */
125b7ce53b6SStefano Zampini   PetscInt               bs,rows,cols;
126b7ce53b6SStefano Zampini   PetscInt               lrows,lcols;
127b7ce53b6SStefano Zampini   PetscInt               local_rows,local_cols;
128b7ce53b6SStefano Zampini   PetscBool              isdense,issbaij,issbaij_red;
129*7c03b4e8SStefano Zampini   PetscMPIInt            nsubdomains;
130b7ce53b6SStefano Zampini   /* values insertion */
131b7ce53b6SStefano Zampini   PetscScalar            *array;
132b7ce53b6SStefano Zampini   PetscInt               *local_indices,*global_indices;
133b7ce53b6SStefano Zampini   /* work */
134b7ce53b6SStefano Zampini   PetscInt               i,j,index_row;
135b7ce53b6SStefano Zampini   PetscErrorCode         ierr;
136b7ce53b6SStefano Zampini 
137b7ce53b6SStefano Zampini   PetscFunctionBegin;
138b7ce53b6SStefano Zampini   /* MISSING CHECKS
139b7ce53b6SStefano Zampini     - rectangular case not covered (it is not allowed by MATIS)
140b7ce53b6SStefano Zampini   */
141b7ce53b6SStefano Zampini   /* get info from mat */
142*7c03b4e8SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
143*7c03b4e8SStefano Zampini   if (nsubdomains == 1) {
144*7c03b4e8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
145*7c03b4e8SStefano Zampini       ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr);
146*7c03b4e8SStefano Zampini     } else {
147*7c03b4e8SStefano Zampini       ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
148*7c03b4e8SStefano Zampini     }
149*7c03b4e8SStefano Zampini     PetscFunctionReturn(0);
150*7c03b4e8SStefano Zampini   }
151b7ce53b6SStefano Zampini   /* ierr = MatGetLocalToGlobalMapping(mat,&rmapping,&cmapping);CHKERRQ(ierr); */
152b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
153b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
154b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
155b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
156b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
157b7ce53b6SStefano Zampini 
158b7ce53b6SStefano Zampini   /* work */
159b7ce53b6SStefano Zampini   ierr = PetscMalloc1(local_rows,&local_indices);CHKERRQ(ierr);
160b7ce53b6SStefano Zampini   for (i=0;i<local_rows;i++) local_indices[i]=i;
161b7ce53b6SStefano Zampini   /* map indices of local mat to global values */
162b7ce53b6SStefano Zampini   ierr = PetscMalloc(PetscMax(local_cols,local_rows)*sizeof(*global_indices),&global_indices);CHKERRQ(ierr);
163b7ce53b6SStefano Zampini   /* ierr = ISLocalToGlobalMappingApply(rmapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); */
164b7ce53b6SStefano Zampini   ierr = ISLocalToGlobalMappingApply(matis->mapping,local_rows,local_indices,global_indices);CHKERRQ(ierr);
165b7ce53b6SStefano Zampini 
166b7ce53b6SStefano Zampini   if (issbaij) {
167b7ce53b6SStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
168b7ce53b6SStefano Zampini   }
169b7ce53b6SStefano Zampini 
170b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
171b7ce53b6SStefano Zampini     Mat         new_mat;
172b7ce53b6SStefano Zampini     MatType     new_mat_type;
173b7ce53b6SStefano Zampini     Vec         vec_dnz,vec_onz;
174b7ce53b6SStefano Zampini     PetscScalar *my_dnz,*my_onz;
175b7ce53b6SStefano Zampini     PetscInt    *dnz,*onz,*mat_ranges,*row_ownership;
176b7ce53b6SStefano Zampini     PetscInt    index_col,owner;
177b7ce53b6SStefano Zampini 
178b7ce53b6SStefano Zampini     /* determining new matrix type */
179b7ce53b6SStefano Zampini     ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
180b7ce53b6SStefano Zampini     if (issbaij_red) {
181b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
182b7ce53b6SStefano Zampini     } else {
183b7ce53b6SStefano Zampini       if (bs>1) {
184b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
185b7ce53b6SStefano Zampini       } else {
186b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
187b7ce53b6SStefano Zampini       }
188b7ce53b6SStefano Zampini     }
189b7ce53b6SStefano Zampini 
190b7ce53b6SStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),&new_mat);CHKERRQ(ierr);
191b7ce53b6SStefano Zampini     ierr = MatSetSizes(new_mat,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr);
192b7ce53b6SStefano Zampini     ierr = MatSetBlockSize(new_mat,bs);CHKERRQ(ierr);
193b7ce53b6SStefano Zampini     ierr = MatSetType(new_mat,new_mat_type);CHKERRQ(ierr);
194b7ce53b6SStefano Zampini     ierr = MatSetUp(new_mat);CHKERRQ(ierr);
195b7ce53b6SStefano Zampini     ierr = MatGetLocalSize(new_mat,&lrows,&lcols);CHKERRQ(ierr);
196b7ce53b6SStefano Zampini 
197b7ce53b6SStefano Zampini     /*
198b7ce53b6SStefano Zampini       preallocation
199b7ce53b6SStefano Zampini     */
200b7ce53b6SStefano Zampini 
201b7ce53b6SStefano Zampini     ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)new_mat),lrows,lcols,dnz,onz);CHKERRQ(ierr);
202b7ce53b6SStefano Zampini     /*
203b7ce53b6SStefano Zampini        Some vectors are needed to sum up properly on shared interface dofs.
204b7ce53b6SStefano Zampini        Preallocation macros cannot do the job.
205b7ce53b6SStefano Zampini        Note that preallocation is not exact, since it overestimates nonzeros
206b7ce53b6SStefano Zampini     */
2072a7a6963SBarry Smith     ierr = MatCreateVecs(new_mat,NULL,&vec_dnz);CHKERRQ(ierr);
208b7ce53b6SStefano Zampini     /* ierr = VecSetLocalToGlobalMapping(vec_dnz,rmapping);CHKERRQ(ierr); */
209b7ce53b6SStefano Zampini     ierr = VecSetLocalToGlobalMapping(vec_dnz,matis->mapping);CHKERRQ(ierr);
210b7ce53b6SStefano Zampini     ierr = VecDuplicate(vec_dnz,&vec_onz);CHKERRQ(ierr);
211b7ce53b6SStefano Zampini     /* All processes need to compute entire row ownership */
212b7ce53b6SStefano Zampini     ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
213b7ce53b6SStefano Zampini     ierr = MatGetOwnershipRanges(new_mat,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
214b7ce53b6SStefano Zampini     for (i=0;i<nsubdomains;i++) {
215b7ce53b6SStefano Zampini       for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
216b7ce53b6SStefano Zampini         row_ownership[j]=i;
217b7ce53b6SStefano Zampini       }
218b7ce53b6SStefano Zampini     }
219b7ce53b6SStefano Zampini 
220b7ce53b6SStefano Zampini     /*
221b7ce53b6SStefano Zampini        my_dnz and my_onz contains exact contribution to preallocation from each local mat
222b7ce53b6SStefano Zampini        then, they will be summed up properly. This way, preallocation is always sufficient
223b7ce53b6SStefano Zampini     */
224b7ce53b6SStefano Zampini     ierr = PetscMalloc1(local_rows,&my_dnz);CHKERRQ(ierr);
225b7ce53b6SStefano Zampini     ierr = PetscMalloc1(local_rows,&my_onz);CHKERRQ(ierr);
226b7ce53b6SStefano Zampini     ierr = PetscMemzero(my_dnz,local_rows*sizeof(*my_dnz));CHKERRQ(ierr);
227b7ce53b6SStefano Zampini     ierr = PetscMemzero(my_onz,local_rows*sizeof(*my_onz));CHKERRQ(ierr);
228b7ce53b6SStefano Zampini     /* preallocation as a MATAIJ */
229b7ce53b6SStefano Zampini     if (isdense) { /* special case for dense local matrices */
230b7ce53b6SStefano Zampini       for (i=0;i<local_rows;i++) {
231b7ce53b6SStefano Zampini         index_row = global_indices[i];
232b7ce53b6SStefano Zampini         for (j=i;j<local_rows;j++) {
233b7ce53b6SStefano Zampini           owner = row_ownership[index_row];
234b7ce53b6SStefano Zampini           index_col = global_indices[j];
235b7ce53b6SStefano Zampini           if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
236b7ce53b6SStefano Zampini             my_dnz[i] += 1.0;
237b7ce53b6SStefano Zampini           } else { /* offdiag block */
238b7ce53b6SStefano Zampini             my_onz[i] += 1.0;
239b7ce53b6SStefano Zampini           }
240b7ce53b6SStefano Zampini           /* same as before, interchanging rows and cols */
241b7ce53b6SStefano Zampini           if (i != j) {
242b7ce53b6SStefano Zampini             owner = row_ownership[index_col];
243b7ce53b6SStefano Zampini             if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
244b7ce53b6SStefano Zampini               my_dnz[j] += 1.0;
245b7ce53b6SStefano Zampini             } else {
246b7ce53b6SStefano Zampini               my_onz[j] += 1.0;
247b7ce53b6SStefano Zampini             }
248b7ce53b6SStefano Zampini           }
249b7ce53b6SStefano Zampini         }
250b7ce53b6SStefano Zampini       }
251b7ce53b6SStefano Zampini     } else {
252b7ce53b6SStefano Zampini       for (i=0;i<local_rows;i++) {
253b7ce53b6SStefano Zampini         PetscInt ncols;
254b7ce53b6SStefano Zampini         const PetscInt *cols;
255b7ce53b6SStefano Zampini         index_row = global_indices[i];
256b7ce53b6SStefano Zampini         ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
257b7ce53b6SStefano Zampini         for (j=0;j<ncols;j++) {
258b7ce53b6SStefano Zampini           owner = row_ownership[index_row];
259b7ce53b6SStefano Zampini           index_col = global_indices[cols[j]];
260b7ce53b6SStefano Zampini           if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
261b7ce53b6SStefano Zampini             my_dnz[i] += 1.0;
262b7ce53b6SStefano Zampini           } else { /* offdiag block */
263b7ce53b6SStefano Zampini             my_onz[i] += 1.0;
264b7ce53b6SStefano Zampini           }
265b7ce53b6SStefano Zampini           /* same as before, interchanging rows and cols */
266b7ce53b6SStefano Zampini           if (issbaij) {
267b7ce53b6SStefano Zampini             owner = row_ownership[index_col];
268b7ce53b6SStefano Zampini             if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
269b7ce53b6SStefano Zampini               my_dnz[j] += 1.0;
270b7ce53b6SStefano Zampini             } else {
271b7ce53b6SStefano Zampini               my_onz[j] += 1.0;
272b7ce53b6SStefano Zampini             }
273b7ce53b6SStefano Zampini           }
274b7ce53b6SStefano Zampini         }
275b7ce53b6SStefano Zampini         ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
276b7ce53b6SStefano Zampini       }
277b7ce53b6SStefano Zampini     }
278b7ce53b6SStefano Zampini     ierr = VecSet(vec_dnz,0.0);CHKERRQ(ierr);
279b7ce53b6SStefano Zampini     ierr = VecSet(vec_onz,0.0);CHKERRQ(ierr);
280b7ce53b6SStefano Zampini     if (local_rows) { /* multilevel guard */
281b7ce53b6SStefano Zampini       ierr = VecSetValuesLocal(vec_dnz,local_rows,local_indices,my_dnz,ADD_VALUES);CHKERRQ(ierr);
282b7ce53b6SStefano Zampini       ierr = VecSetValuesLocal(vec_onz,local_rows,local_indices,my_onz,ADD_VALUES);CHKERRQ(ierr);
283b7ce53b6SStefano Zampini     }
284b7ce53b6SStefano Zampini     ierr = VecAssemblyBegin(vec_dnz);CHKERRQ(ierr);
285b7ce53b6SStefano Zampini     ierr = VecAssemblyBegin(vec_onz);CHKERRQ(ierr);
286b7ce53b6SStefano Zampini     ierr = VecAssemblyEnd(vec_dnz);CHKERRQ(ierr);
287b7ce53b6SStefano Zampini     ierr = VecAssemblyEnd(vec_onz);CHKERRQ(ierr);
288b7ce53b6SStefano Zampini     ierr = PetscFree(my_dnz);CHKERRQ(ierr);
289b7ce53b6SStefano Zampini     ierr = PetscFree(my_onz);CHKERRQ(ierr);
290b7ce53b6SStefano Zampini     ierr = PetscFree(row_ownership);CHKERRQ(ierr);
291b7ce53b6SStefano Zampini 
292b7ce53b6SStefano Zampini     /* set computed preallocation in dnz and onz */
293b7ce53b6SStefano Zampini     ierr = VecGetArray(vec_dnz,&array);CHKERRQ(ierr);
294b7ce53b6SStefano Zampini     for (i=0; i<lrows; i++) dnz[i] = (PetscInt)PetscRealPart(array[i]);
295b7ce53b6SStefano Zampini     ierr = VecRestoreArray(vec_dnz,&array);CHKERRQ(ierr);
296b7ce53b6SStefano Zampini     ierr = VecGetArray(vec_onz,&array);CHKERRQ(ierr);
297b7ce53b6SStefano Zampini     for (i=0;i<lrows;i++) onz[i] = (PetscInt)PetscRealPart(array[i]);
298b7ce53b6SStefano Zampini     ierr = VecRestoreArray(vec_onz,&array);CHKERRQ(ierr);
299b7ce53b6SStefano Zampini     ierr = VecDestroy(&vec_dnz);CHKERRQ(ierr);
300b7ce53b6SStefano Zampini     ierr = VecDestroy(&vec_onz);CHKERRQ(ierr);
301b7ce53b6SStefano Zampini 
302b7ce53b6SStefano Zampini     /* Resize preallocation if overestimated */
303b7ce53b6SStefano Zampini     for (i=0;i<lrows;i++) {
304b7ce53b6SStefano Zampini       dnz[i] = PetscMin(dnz[i],lcols);
305b7ce53b6SStefano Zampini       onz[i] = PetscMin(onz[i],cols-lcols);
306b7ce53b6SStefano Zampini     }
307b7ce53b6SStefano Zampini     /* set preallocation */
308b7ce53b6SStefano Zampini     ierr = MatMPIAIJSetPreallocation(new_mat,0,dnz,0,onz);CHKERRQ(ierr);
309b7ce53b6SStefano Zampini     for (i=0;i<lrows/bs;i++) {
310b7ce53b6SStefano Zampini       dnz[i] = dnz[i*bs]/bs;
311b7ce53b6SStefano Zampini       onz[i] = onz[i*bs]/bs;
312b7ce53b6SStefano Zampini     }
313b7ce53b6SStefano Zampini     ierr = MatMPIBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr);
314b7ce53b6SStefano Zampini     for (i=0;i<lrows/bs;i++) {
315b7ce53b6SStefano Zampini       dnz[i] = dnz[i]-i;
316b7ce53b6SStefano Zampini     }
317b7ce53b6SStefano Zampini     ierr = MatMPISBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr);
318b7ce53b6SStefano Zampini     ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
319b7ce53b6SStefano Zampini     *M = new_mat;
320b7ce53b6SStefano Zampini   } else {
321b7ce53b6SStefano Zampini     PetscInt mbs,mrows,mcols;
322b7ce53b6SStefano Zampini     /* some checks */
323b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
324b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
325b7ce53b6SStefano Zampini     if (mrows != rows) {
326b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
327b7ce53b6SStefano Zampini     }
328b7ce53b6SStefano Zampini     if (mrows != rows) {
329b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
330b7ce53b6SStefano Zampini     }
331b7ce53b6SStefano Zampini     if (mbs != bs) {
332b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
333b7ce53b6SStefano Zampini     }
334b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
335b7ce53b6SStefano Zampini   }
336b7ce53b6SStefano Zampini   /* set local to global mappings */
337b7ce53b6SStefano Zampini   /* ierr = MatSetLocalToGlobalMapping(*M,rmapping,cmapping);CHKERRQ(ierr); */
338b7ce53b6SStefano Zampini   /* Set values */
339b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
340b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
341b7ce53b6SStefano Zampini     ierr = MatDenseGetArray(matis->A,&array);CHKERRQ(ierr);
342b7ce53b6SStefano Zampini     ierr = MatSetValues(*M,local_rows,global_indices,local_cols,global_indices,array,ADD_VALUES);CHKERRQ(ierr);
343b7ce53b6SStefano Zampini     ierr = MatDenseRestoreArray(matis->A,&array);CHKERRQ(ierr);
344b7ce53b6SStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
345b7ce53b6SStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
346b7ce53b6SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
347b7ce53b6SStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
348b7ce53b6SStefano Zampini     for (i=0;i<local_rows;i++) {
349b7ce53b6SStefano Zampini       ierr = MatGetRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
350b7ce53b6SStefano Zampini       /* ierr = MatSetValuesLocal(*M,1,&i,j,local_indices,array,ADD_VALUES);CHKERRQ(ierr); */
351b7ce53b6SStefano Zampini       ierr = ISLocalToGlobalMappingApply(matis->mapping,j,local_indices,global_indices);CHKERRQ(ierr);
352b7ce53b6SStefano Zampini       ierr = ISLocalToGlobalMappingApply(matis->mapping,1,&i,&index_row);CHKERRQ(ierr);
353b7ce53b6SStefano Zampini       ierr = MatSetValues(*M,1,&index_row,j,global_indices,array,ADD_VALUES);CHKERRQ(ierr);
354b7ce53b6SStefano Zampini       ierr = MatRestoreRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
355b7ce53b6SStefano Zampini     }
356b7ce53b6SStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
357b7ce53b6SStefano Zampini   }
358b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
359b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
360b7ce53b6SStefano Zampini   if (isdense) {
361b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
362b7ce53b6SStefano Zampini   }
363b7ce53b6SStefano Zampini   if (issbaij) {
364b7ce53b6SStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
365b7ce53b6SStefano Zampini   }
366b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
367b7ce53b6SStefano Zampini }
368b7ce53b6SStefano Zampini 
369b7ce53b6SStefano Zampini #undef __FUNCT__
370b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
371b7ce53b6SStefano Zampini /*@
372b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
373b7ce53b6SStefano Zampini 
374b7ce53b6SStefano Zampini   Input Parameter:
375b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
376b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
377b7ce53b6SStefano Zampini 
378b7ce53b6SStefano Zampini   Output Parameter:
379b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
380b7ce53b6SStefano Zampini 
381b7ce53b6SStefano Zampini   Level: developer
382b7ce53b6SStefano Zampini 
383b7ce53b6SStefano Zampini   Notes:
384b7ce53b6SStefano Zampini 
385b7ce53b6SStefano Zampini .seealso: MATIS
386b7ce53b6SStefano Zampini @*/
387b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
388b7ce53b6SStefano Zampini {
389b7ce53b6SStefano Zampini   PetscErrorCode ierr;
390b7ce53b6SStefano Zampini 
391b7ce53b6SStefano Zampini   PetscFunctionBegin;
392b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
393b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
394b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
395b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
396b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
397b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
398b7ce53b6SStefano Zampini   }
399b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
400b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
401b7ce53b6SStefano Zampini }
402b7ce53b6SStefano Zampini 
403b7ce53b6SStefano Zampini #undef __FUNCT__
404ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
405ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
406ad6194a2SStefano Zampini {
407ad6194a2SStefano Zampini   PetscErrorCode ierr;
408ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
409ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
410ad6194a2SStefano Zampini   Mat            B,localmat;
411ad6194a2SStefano Zampini 
412ad6194a2SStefano Zampini   PetscFunctionBegin;
413ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
414ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
415ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
416ad6194a2SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr);
417ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
418ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
419b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
420ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
421ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
422ad6194a2SStefano Zampini   *newmat = B;
423ad6194a2SStefano Zampini   PetscFunctionReturn(0);
424ad6194a2SStefano Zampini }
425ad6194a2SStefano Zampini 
426ad6194a2SStefano Zampini #undef __FUNCT__
42769796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
42869796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
42969796d55SStefano Zampini {
43069796d55SStefano Zampini   PetscErrorCode ierr;
43169796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
43269796d55SStefano Zampini   PetscBool      local_sym;
43369796d55SStefano Zampini 
43469796d55SStefano Zampini   PetscFunctionBegin;
43569796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
43669796d55SStefano Zampini   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
43769796d55SStefano Zampini   PetscFunctionReturn(0);
43869796d55SStefano Zampini }
43969796d55SStefano Zampini 
44069796d55SStefano Zampini #undef __FUNCT__
44169796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
44269796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
44369796d55SStefano Zampini {
44469796d55SStefano Zampini   PetscErrorCode ierr;
44569796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
44669796d55SStefano Zampini   PetscBool      local_sym;
44769796d55SStefano Zampini 
44869796d55SStefano Zampini   PetscFunctionBegin;
44969796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
45069796d55SStefano Zampini   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
45169796d55SStefano Zampini   PetscFunctionReturn(0);
45269796d55SStefano Zampini }
45369796d55SStefano Zampini 
45469796d55SStefano Zampini #undef __FUNCT__
455b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
456dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
457b4319ba4SBarry Smith {
458dfbe8321SBarry Smith   PetscErrorCode ierr;
459b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
460b4319ba4SBarry Smith 
461b4319ba4SBarry Smith   PetscFunctionBegin;
4626bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
4636bf464f9SBarry Smith   ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr);
4646bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
4656bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
4666bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr);
467bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
468dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
469bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
470b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
471b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
4722e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
473b4319ba4SBarry Smith   PetscFunctionReturn(0);
474b4319ba4SBarry Smith }
475b4319ba4SBarry Smith 
476b4319ba4SBarry Smith #undef __FUNCT__
477b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
478dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
479b4319ba4SBarry Smith {
480dfbe8321SBarry Smith   PetscErrorCode ierr;
481b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
482b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
483b4319ba4SBarry Smith 
484b4319ba4SBarry Smith   PetscFunctionBegin;
485b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
486ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
487ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
488b4319ba4SBarry Smith 
489b4319ba4SBarry Smith   /* multiply the local matrix */
490b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
491b4319ba4SBarry Smith 
492b4319ba4SBarry Smith   /* scatter product back into global memory */
4932dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
494ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
495ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
496b4319ba4SBarry Smith   PetscFunctionReturn(0);
497b4319ba4SBarry Smith }
498b4319ba4SBarry Smith 
499b4319ba4SBarry Smith #undef __FUNCT__
5002e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
501b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5022e74eeadSLisandro Dalcin {
503650997f4SStefano Zampini   Vec            temp_vec;
5042e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5052e74eeadSLisandro Dalcin 
5062e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
507650997f4SStefano Zampini   if (v3 != v2) {
508650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
509650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
510650997f4SStefano Zampini   } else {
511650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
512650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
513650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
514650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
515650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
516650997f4SStefano Zampini   }
5172e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5182e74eeadSLisandro Dalcin }
5192e74eeadSLisandro Dalcin 
5202e74eeadSLisandro Dalcin #undef __FUNCT__
5212e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
5222e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
5232e74eeadSLisandro Dalcin {
5242e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
5252e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5262e74eeadSLisandro Dalcin 
5272e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  y = A' * x */
5282e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
529ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
530ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5312e74eeadSLisandro Dalcin 
5322e74eeadSLisandro Dalcin   /* multiply the local matrix */
5332e74eeadSLisandro Dalcin   ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr);
5342e74eeadSLisandro Dalcin 
5352e74eeadSLisandro Dalcin   /* scatter product back into global vector */
5362e74eeadSLisandro Dalcin   ierr = VecSet(y,0);CHKERRQ(ierr);
537ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
538ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5392e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5402e74eeadSLisandro Dalcin }
5412e74eeadSLisandro Dalcin 
5422e74eeadSLisandro Dalcin #undef __FUNCT__
5432e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
5442e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5452e74eeadSLisandro Dalcin {
546650997f4SStefano Zampini   Vec            temp_vec;
5472e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5482e74eeadSLisandro Dalcin 
5492e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
550650997f4SStefano Zampini   if (v3 != v2) {
551650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
552650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
553650997f4SStefano Zampini   } else {
554650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
555650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
556650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
557650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
558650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
559650997f4SStefano Zampini   }
5602e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5612e74eeadSLisandro Dalcin }
5622e74eeadSLisandro Dalcin 
5632e74eeadSLisandro Dalcin #undef __FUNCT__
564b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
565dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
566b4319ba4SBarry Smith {
567b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
568dfbe8321SBarry Smith   PetscErrorCode ierr;
569b4319ba4SBarry Smith   PetscViewer    sviewer;
570b4319ba4SBarry Smith 
571b4319ba4SBarry Smith   PetscFunctionBegin;
572dd2fa690SBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
573b4319ba4SBarry Smith   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
574b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
575b4319ba4SBarry Smith   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
576b4319ba4SBarry Smith   PetscFunctionReturn(0);
577b4319ba4SBarry Smith }
578b4319ba4SBarry Smith 
579b4319ba4SBarry Smith #undef __FUNCT__
580b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
581784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
582b4319ba4SBarry Smith {
583dfbe8321SBarry Smith   PetscErrorCode ierr;
5844e4c7dbeSStefano Zampini   PetscInt       n,bs;
585b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
586b4319ba4SBarry Smith   IS             from,to;
587b4319ba4SBarry Smith   Vec            global;
588b4319ba4SBarry Smith 
589b4319ba4SBarry Smith   PetscFunctionBegin;
590784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
591ce94432eSBarry Smith   if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical");
59270cf5478SStefano Zampini   if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */
59370cf5478SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
59470cf5478SStefano Zampini     ierr = VecDestroy(&is->x);CHKERRQ(ierr);
59570cf5478SStefano Zampini     ierr = VecDestroy(&is->y);CHKERRQ(ierr);
59670cf5478SStefano Zampini     ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr);
5971c47cb0fSStefano Zampini     ierr = MatDestroy(&is->A);CHKERRQ(ierr);
59870cf5478SStefano Zampini   }
599784ac674SJed Brown   ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr);
6006bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
601784ac674SJed Brown   is->mapping = rmapping;
602fa7f1dd8SStefano Zampini /*
603fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
604fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
605fa7f1dd8SStefano Zampini */
606b4319ba4SBarry Smith 
607b4319ba4SBarry Smith   /* Create the local matrix A */
608784ac674SJed Brown   ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr);
6092e1947a5SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&bs);CHKERRQ(ierr);
610f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
6112e1947a5SStefano Zampini   if (bs > 1) {
6122e1947a5SStefano Zampini     ierr = MatSetType(is->A,MATSEQBAIJ);CHKERRQ(ierr);
6132e1947a5SStefano Zampini   } else {
6142e1947a5SStefano Zampini     ierr = MatSetType(is->A,MATSEQAIJ);CHKERRQ(ierr);
6152e1947a5SStefano Zampini   }
616f69a0ea3SMatthew Knepley   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
6174e4c7dbeSStefano Zampini   ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr);
618ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
619ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
620b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
621b4319ba4SBarry Smith 
622b4319ba4SBarry Smith   /* Create the local work vectors */
6234e4c7dbeSStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr);
6244e4c7dbeSStefano Zampini   ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr);
6254e4c7dbeSStefano Zampini   ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr);
626ff130e51SJed Brown   ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr);
627ff130e51SJed Brown   ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr);
6284e4c7dbeSStefano Zampini   ierr = VecSetFromOptions(is->x);CHKERRQ(ierr);
629b4319ba4SBarry Smith   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
630b4319ba4SBarry Smith 
631b4319ba4SBarry Smith   /* setup the global to local scatter */
632b4319ba4SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
633784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
6342a7a6963SBarry Smith   ierr = MatCreateVecs(A,&global,NULL);CHKERRQ(ierr);
635b4319ba4SBarry Smith   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
6366bf464f9SBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
6376bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
6386bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
639b4319ba4SBarry Smith   PetscFunctionReturn(0);
640b4319ba4SBarry Smith }
641b4319ba4SBarry Smith 
6422e74eeadSLisandro Dalcin #undef __FUNCT__
6432e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
6442e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
6452e74eeadSLisandro Dalcin {
6462e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
6472e74eeadSLisandro Dalcin   PetscInt       rows_l[2048],cols_l[2048];
6482e74eeadSLisandro Dalcin   PetscErrorCode ierr;
6492e74eeadSLisandro Dalcin 
6502e74eeadSLisandro Dalcin   PetscFunctionBegin;
6512e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
652f23aa3ddSBarry Smith   if (m > 2048 || n > 2048) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n);
6532e74eeadSLisandro Dalcin #endif
6542e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr);
6552e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr);
6562e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
6572e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
6582e74eeadSLisandro Dalcin }
6592e74eeadSLisandro Dalcin 
6602e74eeadSLisandro Dalcin #undef ISG2LMapSetUp
6612e74eeadSLisandro Dalcin #undef ISG2LMapApply
662b4319ba4SBarry Smith 
663b4319ba4SBarry Smith #undef __FUNCT__
664b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
66513f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
666b4319ba4SBarry Smith {
667dfbe8321SBarry Smith   PetscErrorCode ierr;
668b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
669b4319ba4SBarry Smith 
670b4319ba4SBarry Smith   PetscFunctionBegin;
671b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
672b4319ba4SBarry Smith   PetscFunctionReturn(0);
673b4319ba4SBarry Smith }
674b4319ba4SBarry Smith 
675b4319ba4SBarry Smith #undef __FUNCT__
676f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
677f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
678f0006bf2SLisandro Dalcin {
679f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
680f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
681f0006bf2SLisandro Dalcin 
682f0006bf2SLisandro Dalcin   PetscFunctionBegin;
683f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
684f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
685f0006bf2SLisandro Dalcin }
686f0006bf2SLisandro Dalcin 
687f0006bf2SLisandro Dalcin #undef __FUNCT__
6882e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
6892b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
6902e74eeadSLisandro Dalcin {
6912e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
6920298fd71SBarry Smith   PetscInt       n_l = 0, *rows_l = NULL;
6932e74eeadSLisandro Dalcin   PetscErrorCode ierr;
6942e74eeadSLisandro Dalcin 
6952e74eeadSLisandro Dalcin   PetscFunctionBegin;
696ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
6972e74eeadSLisandro Dalcin   if (n) {
698785e854fSJed Brown     ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr);
6992e74eeadSLisandro Dalcin     ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
7002e74eeadSLisandro Dalcin   }
7012b40b63fSBarry Smith   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
702c31cb41cSBarry Smith   ierr = PetscFree(rows_l);CHKERRQ(ierr);
7032e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7042e74eeadSLisandro Dalcin }
7052e74eeadSLisandro Dalcin 
7062e74eeadSLisandro Dalcin #undef __FUNCT__
707b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
7082b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
709b4319ba4SBarry Smith {
710b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
711dfbe8321SBarry Smith   PetscErrorCode ierr;
712f4df32b1SMatthew Knepley   PetscInt       i;
713b4319ba4SBarry Smith   PetscScalar    *array;
714b4319ba4SBarry Smith 
715b4319ba4SBarry Smith   PetscFunctionBegin;
716ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
717b4319ba4SBarry Smith   {
718b4319ba4SBarry Smith     /*
719b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
720b4319ba4SBarry Smith        work properly in the interface nodes.
721b4319ba4SBarry Smith     */
722b4319ba4SBarry Smith     Vec         counter;
723b4319ba4SBarry Smith     PetscScalar one=1.0, zero=0.0;
7242a7a6963SBarry Smith     ierr = MatCreateVecs(A,&counter,NULL);CHKERRQ(ierr);
7252dcb1b2aSMatthew Knepley     ierr = VecSet(counter,zero);CHKERRQ(ierr);
7262dcb1b2aSMatthew Knepley     ierr = VecSet(is->x,one);CHKERRQ(ierr);
727ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
728ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
729ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
730ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7316bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
732b4319ba4SBarry Smith   }
733958c9bccSBarry Smith   if (!n) {
734b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
735b4319ba4SBarry Smith   } else {
736b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
7372205254eSKarl Rupp 
738b4319ba4SBarry Smith     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
7392b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
740b4319ba4SBarry Smith     for (i=0; i<n; i++) {
741f4df32b1SMatthew Knepley       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
742b4319ba4SBarry Smith     }
743b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
744b4319ba4SBarry Smith     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
745b4319ba4SBarry Smith     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
746b4319ba4SBarry Smith   }
747b4319ba4SBarry Smith   PetscFunctionReturn(0);
748b4319ba4SBarry Smith }
749b4319ba4SBarry Smith 
750b4319ba4SBarry Smith #undef __FUNCT__
751b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
752dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
753b4319ba4SBarry Smith {
754b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
755dfbe8321SBarry Smith   PetscErrorCode ierr;
756b4319ba4SBarry Smith 
757b4319ba4SBarry Smith   PetscFunctionBegin;
758b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
759b4319ba4SBarry Smith   PetscFunctionReturn(0);
760b4319ba4SBarry Smith }
761b4319ba4SBarry Smith 
762b4319ba4SBarry Smith #undef __FUNCT__
763b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
764dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
765b4319ba4SBarry Smith {
766b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
767dfbe8321SBarry Smith   PetscErrorCode ierr;
768b4319ba4SBarry Smith 
769b4319ba4SBarry Smith   PetscFunctionBegin;
770b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
771b4319ba4SBarry Smith   PetscFunctionReturn(0);
772b4319ba4SBarry Smith }
773b4319ba4SBarry Smith 
774b4319ba4SBarry Smith #undef __FUNCT__
775b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
7767087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
777b4319ba4SBarry Smith {
778b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
779b4319ba4SBarry Smith 
780b4319ba4SBarry Smith   PetscFunctionBegin;
781b4319ba4SBarry Smith   *local = is->A;
782b4319ba4SBarry Smith   PetscFunctionReturn(0);
783b4319ba4SBarry Smith }
784b4319ba4SBarry Smith 
785b4319ba4SBarry Smith #undef __FUNCT__
786b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
787b4319ba4SBarry Smith /*@
788b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
789b4319ba4SBarry Smith 
790b4319ba4SBarry Smith   Input Parameter:
791b4319ba4SBarry Smith .  mat - the matrix
792b4319ba4SBarry Smith 
793b4319ba4SBarry Smith   Output Parameter:
794b4319ba4SBarry Smith .  local - the local matrix usually MATSEQAIJ
795b4319ba4SBarry Smith 
796b4319ba4SBarry Smith   Level: advanced
797b4319ba4SBarry Smith 
798b4319ba4SBarry Smith   Notes:
799b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
800b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
801b4319ba4SBarry Smith   of the MatSetValues() operation.
802b4319ba4SBarry Smith 
803b4319ba4SBarry Smith .seealso: MATIS
804b4319ba4SBarry Smith @*/
8057087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
806b4319ba4SBarry Smith {
8074ac538c5SBarry Smith   PetscErrorCode ierr;
808b4319ba4SBarry Smith 
809b4319ba4SBarry Smith   PetscFunctionBegin;
8100700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
811b4319ba4SBarry Smith   PetscValidPointer(local,2);
8124ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
813b4319ba4SBarry Smith   PetscFunctionReturn(0);
814b4319ba4SBarry Smith }
815b4319ba4SBarry Smith 
8163b03a366Sstefano_zampini #undef __FUNCT__
8173b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
8183b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
8193b03a366Sstefano_zampini {
8203b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
8213b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
8223b03a366Sstefano_zampini   PetscErrorCode ierr;
8233b03a366Sstefano_zampini 
8243b03a366Sstefano_zampini   PetscFunctionBegin;
8254e4c7dbeSStefano Zampini   if (is->A) {
8263b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
8273b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
828f23aa3ddSBarry 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);
8294e4c7dbeSStefano Zampini   }
8303b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
8313b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
8323b03a366Sstefano_zampini   is->A = local;
8333b03a366Sstefano_zampini   PetscFunctionReturn(0);
8343b03a366Sstefano_zampini }
8353b03a366Sstefano_zampini 
8363b03a366Sstefano_zampini #undef __FUNCT__
8373b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
8383b03a366Sstefano_zampini /*@
8393b03a366Sstefano_zampini     MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix.
8403b03a366Sstefano_zampini 
8413b03a366Sstefano_zampini   Input Parameter:
8423b03a366Sstefano_zampini .  mat - the matrix
8433b03a366Sstefano_zampini .  local - the local matrix usually MATSEQAIJ
8443b03a366Sstefano_zampini 
8453b03a366Sstefano_zampini   Output Parameter:
8463b03a366Sstefano_zampini 
8473b03a366Sstefano_zampini   Level: advanced
8483b03a366Sstefano_zampini 
8493b03a366Sstefano_zampini   Notes:
8503b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
8513b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
8523b03a366Sstefano_zampini 
8533b03a366Sstefano_zampini .seealso: MATIS
8543b03a366Sstefano_zampini @*/
8553b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
8563b03a366Sstefano_zampini {
8573b03a366Sstefano_zampini   PetscErrorCode ierr;
8583b03a366Sstefano_zampini 
8593b03a366Sstefano_zampini   PetscFunctionBegin;
8603b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
861b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
8623b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
8633b03a366Sstefano_zampini   PetscFunctionReturn(0);
8643b03a366Sstefano_zampini }
8653b03a366Sstefano_zampini 
8666726f965SBarry Smith #undef __FUNCT__
8676726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
8686726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
8696726f965SBarry Smith {
8706726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
8716726f965SBarry Smith   PetscErrorCode ierr;
8726726f965SBarry Smith 
8736726f965SBarry Smith   PetscFunctionBegin;
8746726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
8756726f965SBarry Smith   PetscFunctionReturn(0);
8766726f965SBarry Smith }
8776726f965SBarry Smith 
8786726f965SBarry Smith #undef __FUNCT__
8792e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
8802e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
8812e74eeadSLisandro Dalcin {
8822e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
8832e74eeadSLisandro Dalcin   PetscErrorCode ierr;
8842e74eeadSLisandro Dalcin 
8852e74eeadSLisandro Dalcin   PetscFunctionBegin;
8862e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
8872e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
8882e74eeadSLisandro Dalcin }
8892e74eeadSLisandro Dalcin 
8902e74eeadSLisandro Dalcin #undef __FUNCT__
8912e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
8922e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
8932e74eeadSLisandro Dalcin {
8942e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
8952e74eeadSLisandro Dalcin   PetscErrorCode ierr;
8962e74eeadSLisandro Dalcin 
8972e74eeadSLisandro Dalcin   PetscFunctionBegin;
8982e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
8992e74eeadSLisandro Dalcin   ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr);
9002e74eeadSLisandro Dalcin 
9012e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
9022e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
903ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
904ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9052e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9062e74eeadSLisandro Dalcin }
9072e74eeadSLisandro Dalcin 
9082e74eeadSLisandro Dalcin #undef __FUNCT__
9096726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
910ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
9116726f965SBarry Smith {
9126726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
9136726f965SBarry Smith   PetscErrorCode ierr;
9146726f965SBarry Smith 
9156726f965SBarry Smith   PetscFunctionBegin;
9164e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
9176726f965SBarry Smith   PetscFunctionReturn(0);
9186726f965SBarry Smith }
9196726f965SBarry Smith 
920284134d9SBarry Smith #undef __FUNCT__
921284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
922284134d9SBarry Smith /*@
923284134d9SBarry Smith     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
924284134d9SBarry Smith        process but not across processes.
925284134d9SBarry Smith 
926284134d9SBarry Smith    Input Parameters:
927284134d9SBarry Smith +     comm - MPI communicator that will share the matrix
9284e4c7dbeSStefano Zampini .     bs - local and global block size of the matrix
929df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
930284134d9SBarry Smith -     map - mapping that defines the global number for each local number
931284134d9SBarry Smith 
932284134d9SBarry Smith    Output Parameter:
933284134d9SBarry Smith .    A - the resulting matrix
934284134d9SBarry Smith 
9358e6c10adSSatish Balay    Level: advanced
9368e6c10adSSatish Balay 
937284134d9SBarry Smith    Notes: See MATIS for more details
9388cda6cd7SBarry Smith           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
9398cda6cd7SBarry Smith           by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
9408cda6cd7SBarry Smith           plus the ghost points to global indices.
941284134d9SBarry Smith 
942284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
943284134d9SBarry Smith @*/
9444e4c7dbeSStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
945284134d9SBarry Smith {
946284134d9SBarry Smith   PetscErrorCode ierr;
947284134d9SBarry Smith 
948284134d9SBarry Smith   PetscFunctionBegin;
949284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
950d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
951284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
952284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
9533b03a366Sstefano_zampini   ierr = MatSetUp(*A);CHKERRQ(ierr);
954784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr);
955284134d9SBarry Smith   PetscFunctionReturn(0);
956284134d9SBarry Smith }
957284134d9SBarry Smith 
958b4319ba4SBarry Smith /*MC
9596a818285SBarry Smith    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners, see PCBDDC.
960b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
961b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
962b4319ba4SBarry Smith    product is handled "implicitly".
963b4319ba4SBarry Smith 
964b4319ba4SBarry Smith    Operations Provided:
9656726f965SBarry Smith +  MatMult()
9662e74eeadSLisandro Dalcin .  MatMultAdd()
9672e74eeadSLisandro Dalcin .  MatMultTranspose()
9682e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
9696726f965SBarry Smith .  MatZeroEntries()
9706726f965SBarry Smith .  MatSetOption()
9712e74eeadSLisandro Dalcin .  MatZeroRows()
9726726f965SBarry Smith .  MatZeroRowsLocal()
9732e74eeadSLisandro Dalcin .  MatSetValues()
9746726f965SBarry Smith .  MatSetValuesLocal()
9752e74eeadSLisandro Dalcin .  MatScale()
9762e74eeadSLisandro Dalcin .  MatGetDiagonal()
9776726f965SBarry Smith -  MatSetLocalToGlobalMapping()
978b4319ba4SBarry Smith 
979b4319ba4SBarry Smith    Options Database Keys:
980b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
981b4319ba4SBarry Smith 
982b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
983b4319ba4SBarry Smith 
984b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
985b4319ba4SBarry Smith 
986b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
987b4319ba4SBarry Smith           MatISGetLocalMat()
988b4319ba4SBarry Smith 
989b4319ba4SBarry Smith   Level: advanced
990b4319ba4SBarry Smith 
9916a818285SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), PCBDDC
992b4319ba4SBarry Smith 
993b4319ba4SBarry Smith M*/
994b4319ba4SBarry Smith 
995b4319ba4SBarry Smith #undef __FUNCT__
996b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
9978cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
998b4319ba4SBarry Smith {
999dfbe8321SBarry Smith   PetscErrorCode ierr;
1000b4319ba4SBarry Smith   Mat_IS         *b;
1001b4319ba4SBarry Smith 
1002b4319ba4SBarry Smith   PetscFunctionBegin;
1003b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1004b4319ba4SBarry Smith   A->data = (void*)b;
1005b4319ba4SBarry Smith   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1006b4319ba4SBarry Smith 
1007b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
10082e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
10092e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
10102e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1011b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1012b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
10132e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
1014b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1015f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
10162e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1017b4319ba4SBarry Smith   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
1018b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1019b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1020b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
10216726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
10222e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
10232e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
10246726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
102569796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
102669796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1027ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
1028b4319ba4SBarry Smith 
102926283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
103026283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1031b4319ba4SBarry Smith 
1032b7ce53b6SStefano Zampini   /* zeros pointer */
1033b4319ba4SBarry Smith   b->A       = 0;
1034b4319ba4SBarry Smith   b->ctx     = 0;
1035b4319ba4SBarry Smith   b->x       = 0;
1036b4319ba4SBarry Smith   b->y       = 0;
1037b09366fdSStefano Zampini   b->mapping = 0;
1038b7ce53b6SStefano Zampini 
1039b7ce53b6SStefano Zampini   /* special MATIS functions */
1040bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1041bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1042b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
10432e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
104417667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1045b4319ba4SBarry Smith   PetscFunctionReturn(0);
1046b4319ba4SBarry Smith }
1047