xref: /petsc/src/mat/impls/is/matis.c (revision 3927de2ed0f0aec79c45dc8add25bfe718efd5fa)
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*/
1628f4e0baSStefano Zampini 
17a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat);
18a72627d2SStefano Zampini 
1928f4e0baSStefano Zampini #undef __FUNCT__
2028f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private"
21a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B)
2228f4e0baSStefano Zampini {
2328f4e0baSStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
2428f4e0baSStefano Zampini   const PetscInt *gidxs;
2528f4e0baSStefano Zampini   PetscErrorCode ierr;
2628f4e0baSStefano Zampini 
2728f4e0baSStefano Zampini   PetscFunctionBegin;
2828f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr);
2928f4e0baSStefano Zampini   ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr);
3028f4e0baSStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr);
3128f4e0baSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&gidxs);CHKERRQ(ierr);
3228f4e0baSStefano Zampini   /* PETSC_OWN_POINTER refers to ilocal which is NULL */
3328f4e0baSStefano Zampini   ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
3428f4e0baSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&gidxs);CHKERRQ(ierr);
3528f4e0baSStefano Zampini   ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr);
3628f4e0baSStefano Zampini   PetscFunctionReturn(0);
3728f4e0baSStefano Zampini }
382e1947a5SStefano Zampini 
392e1947a5SStefano Zampini #undef __FUNCT__
402e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation"
41a88811baSStefano Zampini /*
42a88811baSStefano Zampini    MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix.
43a88811baSStefano Zampini 
44a88811baSStefano Zampini    Collective on MPI_Comm
45a88811baSStefano Zampini 
46a88811baSStefano Zampini    Input Parameters:
47a88811baSStefano Zampini +  B - the matrix
48a88811baSStefano Zampini .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
49a88811baSStefano Zampini            (same value is used for all local rows)
50a88811baSStefano Zampini .  d_nnz - array containing the number of nonzeros in the various rows of the
51a88811baSStefano Zampini            DIAGONAL portion of the local submatrix (possibly different for each row)
52a88811baSStefano Zampini            or NULL, if d_nz is used to specify the nonzero structure.
53a88811baSStefano Zampini            The size of this array is equal to the number of local rows, i.e 'm'.
54a88811baSStefano Zampini            For matrices that will be factored, you must leave room for (and set)
55a88811baSStefano Zampini            the diagonal entry even if it is zero.
56a88811baSStefano Zampini .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
57a88811baSStefano Zampini            submatrix (same value is used for all local rows).
58a88811baSStefano Zampini -  o_nnz - array containing the number of nonzeros in the various rows of the
59a88811baSStefano Zampini            OFF-DIAGONAL portion of the local submatrix (possibly different for
60a88811baSStefano Zampini            each row) or NULL, if o_nz is used to specify the nonzero
61a88811baSStefano Zampini            structure. The size of this array is equal to the number
62a88811baSStefano Zampini            of local rows, i.e 'm'.
63a88811baSStefano Zampini 
64a88811baSStefano Zampini    If the *_nnz parameter is given then the *_nz parameter is ignored
65a88811baSStefano Zampini 
66a88811baSStefano Zampini    Level: intermediate
67a88811baSStefano Zampini 
68a88811baSStefano Zampini    Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition
69a88811baSStefano Zampini           from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local
70a88811baSStefano Zampini           matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects.
71a88811baSStefano Zampini 
72a88811baSStefano Zampini .keywords: matrix
73a88811baSStefano Zampini 
74a88811baSStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat()
75a88811baSStefano Zampini @*/
762e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
772e1947a5SStefano Zampini {
782e1947a5SStefano Zampini   PetscErrorCode ierr;
792e1947a5SStefano Zampini 
802e1947a5SStefano Zampini   PetscFunctionBegin;
812e1947a5SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
822e1947a5SStefano Zampini   PetscValidType(B,1);
832e1947a5SStefano Zampini   ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
842e1947a5SStefano Zampini   PetscFunctionReturn(0);
852e1947a5SStefano Zampini }
862e1947a5SStefano Zampini 
872e1947a5SStefano Zampini #undef __FUNCT__
882e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS"
892e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
902e1947a5SStefano Zampini {
912e1947a5SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
9228f4e0baSStefano Zampini   PetscInt       bs,i,nlocalcols;
932e1947a5SStefano Zampini   PetscErrorCode ierr;
942e1947a5SStefano Zampini 
952e1947a5SStefano Zampini   PetscFunctionBegin;
962e1947a5SStefano Zampini   if (!matis->A) {
972e1947a5SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
982e1947a5SStefano Zampini   }
9928f4e0baSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
10028f4e0baSStefano Zampini     ierr = MatISComputeSF_Private(B);CHKERRQ(ierr);
10128f4e0baSStefano Zampini   }
1022e1947a5SStefano Zampini   if (!d_nnz) {
10328f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz;
1042e1947a5SStefano Zampini   } else {
10528f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i];
1062e1947a5SStefano Zampini   }
1072e1947a5SStefano Zampini   if (!o_nnz) {
10828f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz;
1092e1947a5SStefano Zampini   } else {
11028f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i];
1112e1947a5SStefano Zampini   }
11228f4e0baSStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
11328f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr);
11428f4e0baSStefano Zampini   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
11528f4e0baSStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
11628f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves;i++) {
11728f4e0baSStefano Zampini     matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols);
1182e1947a5SStefano Zampini   }
11928f4e0baSStefano Zampini   ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr);
12028f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
12128f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs;
1222e1947a5SStefano Zampini   }
12328f4e0baSStefano Zampini   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
12428f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
12528f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i]-i;
1262e1947a5SStefano Zampini   }
12728f4e0baSStefano Zampini   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
1282e1947a5SStefano Zampini   PetscFunctionReturn(0);
1292e1947a5SStefano Zampini }
130b4319ba4SBarry Smith 
131b4319ba4SBarry Smith #undef __FUNCT__
132*3927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private"
133*3927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
134*3927de2eSStefano Zampini {
135*3927de2eSStefano Zampini   Mat_IS          *matis = (Mat_IS*)(A->data);
136*3927de2eSStefano Zampini   PetscInt        *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership;
137*3927de2eSStefano Zampini   const PetscInt* global_indices;
138*3927de2eSStefano Zampini   PetscInt        i,j,bs,rows,cols;
139*3927de2eSStefano Zampini   PetscInt        lrows,lcols;
140*3927de2eSStefano Zampini   PetscInt        local_rows,local_cols;
141*3927de2eSStefano Zampini   PetscMPIInt     nsubdomains;
142*3927de2eSStefano Zampini   PetscBool       isdense,issbaij;
143*3927de2eSStefano Zampini   PetscErrorCode  ierr;
144*3927de2eSStefano Zampini 
145*3927de2eSStefano Zampini   PetscFunctionBegin;
146*3927de2eSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr);
147*3927de2eSStefano Zampini   ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr);
148*3927de2eSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
149*3927de2eSStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
150*3927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
151*3927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
152*3927de2eSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&global_indices);CHKERRQ(ierr);
153*3927de2eSStefano Zampini   if (issbaij) {
154*3927de2eSStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
155*3927de2eSStefano Zampini   }
156*3927de2eSStefano Zampini   /*
157*3927de2eSStefano Zampini      An SF reduce is needed to sum up properly on shared interface dofs.
158*3927de2eSStefano Zampini      Note that generally preallocation is not exact, since it overestimates nonzeros
159*3927de2eSStefano Zampini   */
160*3927de2eSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
161*3927de2eSStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
162*3927de2eSStefano Zampini   }
163*3927de2eSStefano Zampini   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
164*3927de2eSStefano Zampini   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
165*3927de2eSStefano Zampini   /* All processes need to compute entire row ownership */
166*3927de2eSStefano Zampini   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
167*3927de2eSStefano Zampini   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
168*3927de2eSStefano Zampini   for (i=0;i<nsubdomains;i++) {
169*3927de2eSStefano Zampini     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
170*3927de2eSStefano Zampini       row_ownership[j]=i;
171*3927de2eSStefano Zampini     }
172*3927de2eSStefano Zampini   }
173*3927de2eSStefano Zampini 
174*3927de2eSStefano Zampini   /*
175*3927de2eSStefano Zampini      my_dnz and my_onz contains exact contribution to preallocation from each local mat
176*3927de2eSStefano Zampini      then, they will be summed up properly. This way, preallocation is always sufficient
177*3927de2eSStefano Zampini   */
178*3927de2eSStefano Zampini   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
179*3927de2eSStefano Zampini   /* preallocation as a MATAIJ */
180*3927de2eSStefano Zampini   if (isdense) { /* special case for dense local matrices */
181*3927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
182*3927de2eSStefano Zampini       PetscInt index_row = global_indices[i];
183*3927de2eSStefano Zampini       for (j=i;j<local_rows;j++) {
184*3927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
185*3927de2eSStefano Zampini         PetscInt index_col = global_indices[j];
186*3927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
187*3927de2eSStefano Zampini           my_dnz[i] += 1;
188*3927de2eSStefano Zampini         } else { /* offdiag block */
189*3927de2eSStefano Zampini           my_onz[i] += 1;
190*3927de2eSStefano Zampini         }
191*3927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
192*3927de2eSStefano Zampini         if (i != j) {
193*3927de2eSStefano Zampini           owner = row_ownership[index_col];
194*3927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
195*3927de2eSStefano Zampini             my_dnz[j] += 1;
196*3927de2eSStefano Zampini           } else {
197*3927de2eSStefano Zampini             my_onz[j] += 1;
198*3927de2eSStefano Zampini           }
199*3927de2eSStefano Zampini         }
200*3927de2eSStefano Zampini       }
201*3927de2eSStefano Zampini     }
202*3927de2eSStefano Zampini   } else {
203*3927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
204*3927de2eSStefano Zampini       const PetscInt *cols;
205*3927de2eSStefano Zampini       PetscInt       ncols,index_row = global_indices[i];
206*3927de2eSStefano Zampini       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
207*3927de2eSStefano Zampini       for (j=0;j<ncols;j++) {
208*3927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
209*3927de2eSStefano Zampini         PetscInt index_col = global_indices[cols[j]];
210*3927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
211*3927de2eSStefano Zampini           my_dnz[i] += 1;
212*3927de2eSStefano Zampini         } else { /* offdiag block */
213*3927de2eSStefano Zampini           my_onz[i] += 1;
214*3927de2eSStefano Zampini         }
215*3927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
216*3927de2eSStefano Zampini         if (issbaij) {
217*3927de2eSStefano Zampini           owner = row_ownership[index_col];
218*3927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
219*3927de2eSStefano Zampini             my_dnz[j] += 1;
220*3927de2eSStefano Zampini           } else {
221*3927de2eSStefano Zampini             my_onz[j] += 1;
222*3927de2eSStefano Zampini           }
223*3927de2eSStefano Zampini         }
224*3927de2eSStefano Zampini       }
225*3927de2eSStefano Zampini       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
226*3927de2eSStefano Zampini     }
227*3927de2eSStefano Zampini   }
228*3927de2eSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&global_indices);CHKERRQ(ierr);
229*3927de2eSStefano Zampini   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
230*3927de2eSStefano Zampini   if (maxreduce) {
231*3927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
232*3927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
233*3927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
234*3927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
235*3927de2eSStefano Zampini   } else {
236*3927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
237*3927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
238*3927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
239*3927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
240*3927de2eSStefano Zampini   }
241*3927de2eSStefano Zampini   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
242*3927de2eSStefano Zampini 
243*3927de2eSStefano Zampini   /* Resize preallocation if overestimated */
244*3927de2eSStefano Zampini   for (i=0;i<lrows;i++) {
245*3927de2eSStefano Zampini     dnz[i] = PetscMin(dnz[i],lcols);
246*3927de2eSStefano Zampini     onz[i] = PetscMin(onz[i],cols-lcols);
247*3927de2eSStefano Zampini   }
248*3927de2eSStefano Zampini   /* set preallocation */
249*3927de2eSStefano Zampini   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
250*3927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
251*3927de2eSStefano Zampini     dnz[i] = dnz[i*bs]/bs;
252*3927de2eSStefano Zampini     onz[i] = onz[i*bs]/bs;
253*3927de2eSStefano Zampini   }
254*3927de2eSStefano Zampini   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
255*3927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
256*3927de2eSStefano Zampini     dnz[i] = dnz[i]-i;
257*3927de2eSStefano Zampini   }
258*3927de2eSStefano Zampini   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
259*3927de2eSStefano Zampini   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
260*3927de2eSStefano Zampini   if (issbaij) {
261*3927de2eSStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
262*3927de2eSStefano Zampini   }
263*3927de2eSStefano Zampini   PetscFunctionReturn(0);
264*3927de2eSStefano Zampini }
265*3927de2eSStefano Zampini 
266*3927de2eSStefano Zampini #undef __FUNCT__
267b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
268b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
269b7ce53b6SStefano Zampini {
270b7ce53b6SStefano Zampini   Mat_IS                 *matis = (Mat_IS*)(mat->data);
271b7ce53b6SStefano Zampini   /* info on mat */
272b7ce53b6SStefano Zampini   /* ISLocalToGlobalMapping rmapping,cmapping; */
273b7ce53b6SStefano Zampini   PetscInt               bs,rows,cols;
274b7ce53b6SStefano Zampini   PetscInt               local_rows,local_cols;
275*3927de2eSStefano Zampini   PetscBool              isdense,issbaij;
2767c03b4e8SStefano Zampini   PetscMPIInt            nsubdomains;
277b7ce53b6SStefano Zampini   /* values insertion */
278b7ce53b6SStefano Zampini   PetscScalar            *array;
279b7ce53b6SStefano Zampini   PetscInt               *local_indices,*global_indices;
280b7ce53b6SStefano Zampini   /* work */
281b7ce53b6SStefano Zampini   PetscInt               i,j,index_row;
282b7ce53b6SStefano Zampini   PetscErrorCode         ierr;
283b7ce53b6SStefano Zampini 
284b7ce53b6SStefano Zampini   PetscFunctionBegin;
285b7ce53b6SStefano Zampini   /* MISSING CHECKS
286b7ce53b6SStefano Zampini     - rectangular case not covered (it is not allowed by MATIS)
287b7ce53b6SStefano Zampini   */
288b7ce53b6SStefano Zampini   /* get info from mat */
2897c03b4e8SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
2907c03b4e8SStefano Zampini   if (nsubdomains == 1) {
2917c03b4e8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
2927c03b4e8SStefano Zampini       ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr);
2937c03b4e8SStefano Zampini     } else {
2947c03b4e8SStefano Zampini       ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
2957c03b4e8SStefano Zampini     }
2967c03b4e8SStefano Zampini     PetscFunctionReturn(0);
2977c03b4e8SStefano Zampini   }
298b7ce53b6SStefano Zampini   /* ierr = MatGetLocalToGlobalMapping(mat,&rmapping,&cmapping);CHKERRQ(ierr); */
299b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
300b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
301b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
302b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
303b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
304b7ce53b6SStefano Zampini 
305b7ce53b6SStefano Zampini   /* work */
306b7ce53b6SStefano Zampini   ierr = PetscMalloc1(local_rows,&local_indices);CHKERRQ(ierr);
307b7ce53b6SStefano Zampini   for (i=0;i<local_rows;i++) local_indices[i]=i;
308b7ce53b6SStefano Zampini   /* map indices of local mat to global values */
309854ce69bSBarry Smith   ierr = PetscMalloc1(PetscMax(local_cols,local_rows),&global_indices);CHKERRQ(ierr);
310b7ce53b6SStefano Zampini   /* ierr = ISLocalToGlobalMappingApply(rmapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); */
311b7ce53b6SStefano Zampini   ierr = ISLocalToGlobalMappingApply(matis->mapping,local_rows,local_indices,global_indices);CHKERRQ(ierr);
312b7ce53b6SStefano Zampini 
313b7ce53b6SStefano Zampini   if (issbaij) {
314b7ce53b6SStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
315b7ce53b6SStefano Zampini   }
316b7ce53b6SStefano Zampini 
317b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
318b7ce53b6SStefano Zampini     MatType     new_mat_type;
319*3927de2eSStefano Zampini     PetscBool   issbaij_red;
320b7ce53b6SStefano Zampini 
321b7ce53b6SStefano Zampini     /* determining new matrix type */
322b7ce53b6SStefano Zampini     ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
323b7ce53b6SStefano Zampini     if (issbaij_red) {
324b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
325b7ce53b6SStefano Zampini     } else {
326b7ce53b6SStefano Zampini       if (bs>1) {
327b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
328b7ce53b6SStefano Zampini       } else {
329b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
330b7ce53b6SStefano Zampini       }
331b7ce53b6SStefano Zampini     }
332b7ce53b6SStefano Zampini 
333*3927de2eSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
334*3927de2eSStefano Zampini     ierr = MatSetSizes(*M,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr);
335*3927de2eSStefano Zampini     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
336*3927de2eSStefano Zampini     ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr);
337*3927de2eSStefano Zampini     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
338b7ce53b6SStefano Zampini   } else {
339b7ce53b6SStefano Zampini     PetscInt mbs,mrows,mcols;
340b7ce53b6SStefano Zampini     /* some checks */
341b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
342b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
343b7ce53b6SStefano Zampini     if (mrows != rows) {
344b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
345b7ce53b6SStefano Zampini     }
346b7ce53b6SStefano Zampini     if (mrows != rows) {
347b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
348b7ce53b6SStefano Zampini     }
349b7ce53b6SStefano Zampini     if (mbs != bs) {
350b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
351b7ce53b6SStefano Zampini     }
352b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
353b7ce53b6SStefano Zampini   }
354b7ce53b6SStefano Zampini   /* set local to global mappings */
355b7ce53b6SStefano Zampini   /* ierr = MatSetLocalToGlobalMapping(*M,rmapping,cmapping);CHKERRQ(ierr); */
356b7ce53b6SStefano Zampini   /* Set values */
357b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
358b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
359b7ce53b6SStefano Zampini     ierr = MatDenseGetArray(matis->A,&array);CHKERRQ(ierr);
360b7ce53b6SStefano Zampini     ierr = MatSetValues(*M,local_rows,global_indices,local_cols,global_indices,array,ADD_VALUES);CHKERRQ(ierr);
361b7ce53b6SStefano Zampini     ierr = MatDenseRestoreArray(matis->A,&array);CHKERRQ(ierr);
362b7ce53b6SStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
363b7ce53b6SStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
364b7ce53b6SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
365b7ce53b6SStefano Zampini     ierr = PetscFree(local_indices);CHKERRQ(ierr);
366b7ce53b6SStefano Zampini     for (i=0;i<local_rows;i++) {
367b7ce53b6SStefano Zampini       ierr = MatGetRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
368b7ce53b6SStefano Zampini       /* ierr = MatSetValuesLocal(*M,1,&i,j,local_indices,array,ADD_VALUES);CHKERRQ(ierr); */
369b7ce53b6SStefano Zampini       ierr = ISLocalToGlobalMappingApply(matis->mapping,j,local_indices,global_indices);CHKERRQ(ierr);
370b7ce53b6SStefano Zampini       ierr = ISLocalToGlobalMappingApply(matis->mapping,1,&i,&index_row);CHKERRQ(ierr);
371b7ce53b6SStefano Zampini       ierr = MatSetValues(*M,1,&index_row,j,global_indices,array,ADD_VALUES);CHKERRQ(ierr);
372b7ce53b6SStefano Zampini       ierr = MatRestoreRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
373b7ce53b6SStefano Zampini     }
374b7ce53b6SStefano Zampini     ierr = PetscFree(global_indices);CHKERRQ(ierr);
375b7ce53b6SStefano Zampini   }
376b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
377b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
378b7ce53b6SStefano Zampini   if (isdense) {
379b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
380b7ce53b6SStefano Zampini   }
381b7ce53b6SStefano Zampini   if (issbaij) {
382b7ce53b6SStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
383b7ce53b6SStefano Zampini   }
384b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
385b7ce53b6SStefano Zampini }
386b7ce53b6SStefano Zampini 
387b7ce53b6SStefano Zampini #undef __FUNCT__
388b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
389b7ce53b6SStefano Zampini /*@
390b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
391b7ce53b6SStefano Zampini 
392b7ce53b6SStefano Zampini   Input Parameter:
393b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
394b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
395b7ce53b6SStefano Zampini 
396b7ce53b6SStefano Zampini   Output Parameter:
397b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
398b7ce53b6SStefano Zampini 
399b7ce53b6SStefano Zampini   Level: developer
400b7ce53b6SStefano Zampini 
401b7ce53b6SStefano Zampini   Notes:
402b7ce53b6SStefano Zampini 
403b7ce53b6SStefano Zampini .seealso: MATIS
404b7ce53b6SStefano Zampini @*/
405b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
406b7ce53b6SStefano Zampini {
407b7ce53b6SStefano Zampini   PetscErrorCode ierr;
408b7ce53b6SStefano Zampini 
409b7ce53b6SStefano Zampini   PetscFunctionBegin;
410b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
411b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
412b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
413b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
414b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
415b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
416b7ce53b6SStefano Zampini   }
417b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
418b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
419b7ce53b6SStefano Zampini }
420b7ce53b6SStefano Zampini 
421b7ce53b6SStefano Zampini #undef __FUNCT__
422ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
423ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
424ad6194a2SStefano Zampini {
425ad6194a2SStefano Zampini   PetscErrorCode ierr;
426ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
427ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
428ad6194a2SStefano Zampini   Mat            B,localmat;
429ad6194a2SStefano Zampini 
430ad6194a2SStefano Zampini   PetscFunctionBegin;
431ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
432ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
433ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
434ad6194a2SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr);
435ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
436ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
437b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
438ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
439ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
440ad6194a2SStefano Zampini   *newmat = B;
441ad6194a2SStefano Zampini   PetscFunctionReturn(0);
442ad6194a2SStefano Zampini }
443ad6194a2SStefano Zampini 
444ad6194a2SStefano Zampini #undef __FUNCT__
44569796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
44669796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
44769796d55SStefano Zampini {
44869796d55SStefano Zampini   PetscErrorCode ierr;
44969796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
45069796d55SStefano Zampini   PetscBool      local_sym;
45169796d55SStefano Zampini 
45269796d55SStefano Zampini   PetscFunctionBegin;
45369796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
45469796d55SStefano Zampini   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
45569796d55SStefano Zampini   PetscFunctionReturn(0);
45669796d55SStefano Zampini }
45769796d55SStefano Zampini 
45869796d55SStefano Zampini #undef __FUNCT__
45969796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
46069796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
46169796d55SStefano Zampini {
46269796d55SStefano Zampini   PetscErrorCode ierr;
46369796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
46469796d55SStefano Zampini   PetscBool      local_sym;
46569796d55SStefano Zampini 
46669796d55SStefano Zampini   PetscFunctionBegin;
46769796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
46869796d55SStefano Zampini   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
46969796d55SStefano Zampini   PetscFunctionReturn(0);
47069796d55SStefano Zampini }
47169796d55SStefano Zampini 
47269796d55SStefano Zampini #undef __FUNCT__
473b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
474dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
475b4319ba4SBarry Smith {
476dfbe8321SBarry Smith   PetscErrorCode ierr;
477b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
478b4319ba4SBarry Smith 
479b4319ba4SBarry Smith   PetscFunctionBegin;
4806bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
4816bf464f9SBarry Smith   ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr);
4826bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
4836bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
4846bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr);
48528f4e0baSStefano Zampini   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
48628f4e0baSStefano Zampini   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
487bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
488dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
489bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
490b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
491b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
4922e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
493b4319ba4SBarry Smith   PetscFunctionReturn(0);
494b4319ba4SBarry Smith }
495b4319ba4SBarry Smith 
496b4319ba4SBarry Smith #undef __FUNCT__
497b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
498dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
499b4319ba4SBarry Smith {
500dfbe8321SBarry Smith   PetscErrorCode ierr;
501b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
502b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
503b4319ba4SBarry Smith 
504b4319ba4SBarry Smith   PetscFunctionBegin;
505b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
506ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
507ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
508b4319ba4SBarry Smith 
509b4319ba4SBarry Smith   /* multiply the local matrix */
510b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
511b4319ba4SBarry Smith 
512b4319ba4SBarry Smith   /* scatter product back into global memory */
5132dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
514ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
515ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
516b4319ba4SBarry Smith   PetscFunctionReturn(0);
517b4319ba4SBarry Smith }
518b4319ba4SBarry Smith 
519b4319ba4SBarry Smith #undef __FUNCT__
5202e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
521b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5222e74eeadSLisandro Dalcin {
523650997f4SStefano Zampini   Vec            temp_vec;
5242e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5252e74eeadSLisandro Dalcin 
5262e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
527650997f4SStefano Zampini   if (v3 != v2) {
528650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
529650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
530650997f4SStefano Zampini   } else {
531650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
532650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
533650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
534650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
535650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
536650997f4SStefano Zampini   }
5372e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5382e74eeadSLisandro Dalcin }
5392e74eeadSLisandro Dalcin 
5402e74eeadSLisandro Dalcin #undef __FUNCT__
5412e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
5422e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
5432e74eeadSLisandro Dalcin {
5442e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
5452e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5462e74eeadSLisandro Dalcin 
5472e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  y = A' * x */
5482e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
549ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
550ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5512e74eeadSLisandro Dalcin 
5522e74eeadSLisandro Dalcin   /* multiply the local matrix */
5532e74eeadSLisandro Dalcin   ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr);
5542e74eeadSLisandro Dalcin 
5552e74eeadSLisandro Dalcin   /* scatter product back into global vector */
5562e74eeadSLisandro Dalcin   ierr = VecSet(y,0);CHKERRQ(ierr);
557ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
558ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5592e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5602e74eeadSLisandro Dalcin }
5612e74eeadSLisandro Dalcin 
5622e74eeadSLisandro Dalcin #undef __FUNCT__
5632e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
5642e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5652e74eeadSLisandro Dalcin {
566650997f4SStefano Zampini   Vec            temp_vec;
5672e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5682e74eeadSLisandro Dalcin 
5692e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
570650997f4SStefano Zampini   if (v3 != v2) {
571650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
572650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
573650997f4SStefano Zampini   } else {
574650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
575650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
576650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
577650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
578650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
579650997f4SStefano Zampini   }
5802e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5812e74eeadSLisandro Dalcin }
5822e74eeadSLisandro Dalcin 
5832e74eeadSLisandro Dalcin #undef __FUNCT__
584b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
585dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
586b4319ba4SBarry Smith {
587b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
588dfbe8321SBarry Smith   PetscErrorCode ierr;
589b4319ba4SBarry Smith   PetscViewer    sviewer;
590b4319ba4SBarry Smith 
591b4319ba4SBarry Smith   PetscFunctionBegin;
592dd2fa690SBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
593b4319ba4SBarry Smith   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
594b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
595b4319ba4SBarry Smith   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
596b4319ba4SBarry Smith   PetscFunctionReturn(0);
597b4319ba4SBarry Smith }
598b4319ba4SBarry Smith 
599b4319ba4SBarry Smith #undef __FUNCT__
600b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
601784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
602b4319ba4SBarry Smith {
603dfbe8321SBarry Smith   PetscErrorCode ierr;
6044e4c7dbeSStefano Zampini   PetscInt       n,bs;
605b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
606b4319ba4SBarry Smith   IS             from,to;
607b4319ba4SBarry Smith   Vec            global;
608b4319ba4SBarry Smith 
609b4319ba4SBarry Smith   PetscFunctionBegin;
610784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
611ce94432eSBarry Smith   if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical");
61270cf5478SStefano Zampini   if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */
61370cf5478SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
61470cf5478SStefano Zampini     ierr = VecDestroy(&is->x);CHKERRQ(ierr);
61570cf5478SStefano Zampini     ierr = VecDestroy(&is->y);CHKERRQ(ierr);
61670cf5478SStefano Zampini     ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr);
6171c47cb0fSStefano Zampini     ierr = MatDestroy(&is->A);CHKERRQ(ierr);
61828f4e0baSStefano Zampini     ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
61928f4e0baSStefano Zampini     ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
62070cf5478SStefano Zampini   }
621784ac674SJed Brown   ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr);
6226bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
623784ac674SJed Brown   is->mapping = rmapping;
624fa7f1dd8SStefano Zampini /*
625fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
626fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
627fa7f1dd8SStefano Zampini */
628b4319ba4SBarry Smith 
629b4319ba4SBarry Smith   /* Create the local matrix A */
630784ac674SJed Brown   ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr);
6312e1947a5SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&bs);CHKERRQ(ierr);
632f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
6332e1947a5SStefano Zampini   if (bs > 1) {
6342e1947a5SStefano Zampini     ierr = MatSetType(is->A,MATSEQBAIJ);CHKERRQ(ierr);
6352e1947a5SStefano Zampini   } else {
6362e1947a5SStefano Zampini     ierr = MatSetType(is->A,MATSEQAIJ);CHKERRQ(ierr);
6372e1947a5SStefano Zampini   }
638f69a0ea3SMatthew Knepley   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
6394e4c7dbeSStefano Zampini   ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr);
640ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
641ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
642b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
643b4319ba4SBarry Smith 
644b4319ba4SBarry Smith   /* Create the local work vectors */
6454e4c7dbeSStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr);
6464e4c7dbeSStefano Zampini   ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr);
6474e4c7dbeSStefano Zampini   ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr);
648ff130e51SJed Brown   ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr);
649ff130e51SJed Brown   ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr);
6504e4c7dbeSStefano Zampini   ierr = VecSetFromOptions(is->x);CHKERRQ(ierr);
651b4319ba4SBarry Smith   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
652b4319ba4SBarry Smith 
653b4319ba4SBarry Smith   /* setup the global to local scatter */
654b4319ba4SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
655784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
6562a7a6963SBarry Smith   ierr = MatCreateVecs(A,&global,NULL);CHKERRQ(ierr);
657b4319ba4SBarry Smith   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
6586bf464f9SBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
6596bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
6606bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
661b4319ba4SBarry Smith   PetscFunctionReturn(0);
662b4319ba4SBarry Smith }
663b4319ba4SBarry Smith 
6642e74eeadSLisandro Dalcin #undef __FUNCT__
6652e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
6662e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
6672e74eeadSLisandro Dalcin {
6682e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
6692e74eeadSLisandro Dalcin   PetscInt       rows_l[2048],cols_l[2048];
6702e74eeadSLisandro Dalcin   PetscErrorCode ierr;
6712e74eeadSLisandro Dalcin 
6722e74eeadSLisandro Dalcin   PetscFunctionBegin;
6732e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
674f23aa3ddSBarry 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);
6752e74eeadSLisandro Dalcin #endif
6762e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr);
6772e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr);
6782e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
6792e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
6802e74eeadSLisandro Dalcin }
6812e74eeadSLisandro Dalcin 
6822e74eeadSLisandro Dalcin #undef ISG2LMapSetUp
6832e74eeadSLisandro Dalcin #undef ISG2LMapApply
684b4319ba4SBarry Smith 
685b4319ba4SBarry Smith #undef __FUNCT__
686b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
68713f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
688b4319ba4SBarry Smith {
689dfbe8321SBarry Smith   PetscErrorCode ierr;
690b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
691b4319ba4SBarry Smith 
692b4319ba4SBarry Smith   PetscFunctionBegin;
693b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
694b4319ba4SBarry Smith   PetscFunctionReturn(0);
695b4319ba4SBarry Smith }
696b4319ba4SBarry Smith 
697b4319ba4SBarry Smith #undef __FUNCT__
698f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
699f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
700f0006bf2SLisandro Dalcin {
701f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
702f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
703f0006bf2SLisandro Dalcin 
704f0006bf2SLisandro Dalcin   PetscFunctionBegin;
705f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
706f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
707f0006bf2SLisandro Dalcin }
708f0006bf2SLisandro Dalcin 
709f0006bf2SLisandro Dalcin #undef __FUNCT__
7102e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
7112b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
7122e74eeadSLisandro Dalcin {
7132e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
7140298fd71SBarry Smith   PetscInt       n_l = 0, *rows_l = NULL;
7152e74eeadSLisandro Dalcin   PetscErrorCode ierr;
7162e74eeadSLisandro Dalcin 
7172e74eeadSLisandro Dalcin   PetscFunctionBegin;
718ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
7192e74eeadSLisandro Dalcin   if (n) {
720785e854fSJed Brown     ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr);
7212e74eeadSLisandro Dalcin     ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
7222e74eeadSLisandro Dalcin   }
7232b40b63fSBarry Smith   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
724c31cb41cSBarry Smith   ierr = PetscFree(rows_l);CHKERRQ(ierr);
7252e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7262e74eeadSLisandro Dalcin }
7272e74eeadSLisandro Dalcin 
7282e74eeadSLisandro Dalcin #undef __FUNCT__
729b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
7302b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
731b4319ba4SBarry Smith {
732b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
733dfbe8321SBarry Smith   PetscErrorCode ierr;
734f4df32b1SMatthew Knepley   PetscInt       i;
735b4319ba4SBarry Smith   PetscScalar    *array;
736b4319ba4SBarry Smith 
737b4319ba4SBarry Smith   PetscFunctionBegin;
738ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
739b4319ba4SBarry Smith   {
740b4319ba4SBarry Smith     /*
741b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
742b4319ba4SBarry Smith        work properly in the interface nodes.
743b4319ba4SBarry Smith     */
744b4319ba4SBarry Smith     Vec         counter;
745b4319ba4SBarry Smith     PetscScalar one=1.0, zero=0.0;
7462a7a6963SBarry Smith     ierr = MatCreateVecs(A,&counter,NULL);CHKERRQ(ierr);
7472dcb1b2aSMatthew Knepley     ierr = VecSet(counter,zero);CHKERRQ(ierr);
7482dcb1b2aSMatthew Knepley     ierr = VecSet(is->x,one);CHKERRQ(ierr);
749ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
750ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
751ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
752ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7536bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
754b4319ba4SBarry Smith   }
755958c9bccSBarry Smith   if (!n) {
756b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
757b4319ba4SBarry Smith   } else {
758b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
7592205254eSKarl Rupp 
760b4319ba4SBarry Smith     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
7612b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
762b4319ba4SBarry Smith     for (i=0; i<n; i++) {
763f4df32b1SMatthew Knepley       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
764b4319ba4SBarry Smith     }
765b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
766b4319ba4SBarry Smith     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
767b4319ba4SBarry Smith     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
768b4319ba4SBarry Smith   }
769b4319ba4SBarry Smith   PetscFunctionReturn(0);
770b4319ba4SBarry Smith }
771b4319ba4SBarry Smith 
772b4319ba4SBarry Smith #undef __FUNCT__
773b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
774dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
775b4319ba4SBarry Smith {
776b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
777dfbe8321SBarry Smith   PetscErrorCode ierr;
778b4319ba4SBarry Smith 
779b4319ba4SBarry Smith   PetscFunctionBegin;
780b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
781b4319ba4SBarry Smith   PetscFunctionReturn(0);
782b4319ba4SBarry Smith }
783b4319ba4SBarry Smith 
784b4319ba4SBarry Smith #undef __FUNCT__
785b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
786dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
787b4319ba4SBarry Smith {
788b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
789dfbe8321SBarry Smith   PetscErrorCode ierr;
790b4319ba4SBarry Smith 
791b4319ba4SBarry Smith   PetscFunctionBegin;
792b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
793b4319ba4SBarry Smith   PetscFunctionReturn(0);
794b4319ba4SBarry Smith }
795b4319ba4SBarry Smith 
796b4319ba4SBarry Smith #undef __FUNCT__
797b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
7987087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
799b4319ba4SBarry Smith {
800b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
801b4319ba4SBarry Smith 
802b4319ba4SBarry Smith   PetscFunctionBegin;
803b4319ba4SBarry Smith   *local = is->A;
804b4319ba4SBarry Smith   PetscFunctionReturn(0);
805b4319ba4SBarry Smith }
806b4319ba4SBarry Smith 
807b4319ba4SBarry Smith #undef __FUNCT__
808b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
809b4319ba4SBarry Smith /*@
810b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
811b4319ba4SBarry Smith 
812b4319ba4SBarry Smith   Input Parameter:
813b4319ba4SBarry Smith .  mat - the matrix
814b4319ba4SBarry Smith 
815b4319ba4SBarry Smith   Output Parameter:
816b4319ba4SBarry Smith .  local - the local matrix usually MATSEQAIJ
817b4319ba4SBarry Smith 
818b4319ba4SBarry Smith   Level: advanced
819b4319ba4SBarry Smith 
820b4319ba4SBarry Smith   Notes:
821b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
822b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
823b4319ba4SBarry Smith   of the MatSetValues() operation.
824b4319ba4SBarry Smith 
825b4319ba4SBarry Smith .seealso: MATIS
826b4319ba4SBarry Smith @*/
8277087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
828b4319ba4SBarry Smith {
8294ac538c5SBarry Smith   PetscErrorCode ierr;
830b4319ba4SBarry Smith 
831b4319ba4SBarry Smith   PetscFunctionBegin;
8320700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
833b4319ba4SBarry Smith   PetscValidPointer(local,2);
8344ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
835b4319ba4SBarry Smith   PetscFunctionReturn(0);
836b4319ba4SBarry Smith }
837b4319ba4SBarry Smith 
8383b03a366Sstefano_zampini #undef __FUNCT__
8393b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
8403b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
8413b03a366Sstefano_zampini {
8423b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
8433b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
8443b03a366Sstefano_zampini   PetscErrorCode ierr;
8453b03a366Sstefano_zampini 
8463b03a366Sstefano_zampini   PetscFunctionBegin;
8474e4c7dbeSStefano Zampini   if (is->A) {
8483b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
8493b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
850f23aa3ddSBarry 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);
8514e4c7dbeSStefano Zampini   }
8523b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
8533b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
8543b03a366Sstefano_zampini   is->A = local;
8553b03a366Sstefano_zampini   PetscFunctionReturn(0);
8563b03a366Sstefano_zampini }
8573b03a366Sstefano_zampini 
8583b03a366Sstefano_zampini #undef __FUNCT__
8593b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
8603b03a366Sstefano_zampini /*@
8613b03a366Sstefano_zampini     MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix.
8623b03a366Sstefano_zampini 
8633b03a366Sstefano_zampini   Input Parameter:
8643b03a366Sstefano_zampini .  mat - the matrix
8653b03a366Sstefano_zampini .  local - the local matrix usually MATSEQAIJ
8663b03a366Sstefano_zampini 
8673b03a366Sstefano_zampini   Output Parameter:
8683b03a366Sstefano_zampini 
8693b03a366Sstefano_zampini   Level: advanced
8703b03a366Sstefano_zampini 
8713b03a366Sstefano_zampini   Notes:
8723b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
8733b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
8743b03a366Sstefano_zampini 
8753b03a366Sstefano_zampini .seealso: MATIS
8763b03a366Sstefano_zampini @*/
8773b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
8783b03a366Sstefano_zampini {
8793b03a366Sstefano_zampini   PetscErrorCode ierr;
8803b03a366Sstefano_zampini 
8813b03a366Sstefano_zampini   PetscFunctionBegin;
8823b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
883b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
8843b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
8853b03a366Sstefano_zampini   PetscFunctionReturn(0);
8863b03a366Sstefano_zampini }
8873b03a366Sstefano_zampini 
8886726f965SBarry Smith #undef __FUNCT__
8896726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
8906726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
8916726f965SBarry Smith {
8926726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
8936726f965SBarry Smith   PetscErrorCode ierr;
8946726f965SBarry Smith 
8956726f965SBarry Smith   PetscFunctionBegin;
8966726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
8976726f965SBarry Smith   PetscFunctionReturn(0);
8986726f965SBarry Smith }
8996726f965SBarry Smith 
9006726f965SBarry Smith #undef __FUNCT__
9012e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
9022e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
9032e74eeadSLisandro Dalcin {
9042e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9052e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9062e74eeadSLisandro Dalcin 
9072e74eeadSLisandro Dalcin   PetscFunctionBegin;
9082e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
9092e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9102e74eeadSLisandro Dalcin }
9112e74eeadSLisandro Dalcin 
9122e74eeadSLisandro Dalcin #undef __FUNCT__
9132e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
9142e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
9152e74eeadSLisandro Dalcin {
9162e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9172e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9182e74eeadSLisandro Dalcin 
9192e74eeadSLisandro Dalcin   PetscFunctionBegin;
9202e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
9212e74eeadSLisandro Dalcin   ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr);
9222e74eeadSLisandro Dalcin 
9232e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
9242e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
925ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
926ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9272e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9282e74eeadSLisandro Dalcin }
9292e74eeadSLisandro Dalcin 
9302e74eeadSLisandro Dalcin #undef __FUNCT__
9316726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
932ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
9336726f965SBarry Smith {
9346726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
9356726f965SBarry Smith   PetscErrorCode ierr;
9366726f965SBarry Smith 
9376726f965SBarry Smith   PetscFunctionBegin;
9384e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
9396726f965SBarry Smith   PetscFunctionReturn(0);
9406726f965SBarry Smith }
9416726f965SBarry Smith 
942284134d9SBarry Smith #undef __FUNCT__
943284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
944284134d9SBarry Smith /*@
945284134d9SBarry Smith     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
946284134d9SBarry Smith        process but not across processes.
947284134d9SBarry Smith 
948284134d9SBarry Smith    Input Parameters:
949284134d9SBarry Smith +     comm - MPI communicator that will share the matrix
9504e4c7dbeSStefano Zampini .     bs - local and global block size of the matrix
951df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
952284134d9SBarry Smith -     map - mapping that defines the global number for each local number
953284134d9SBarry Smith 
954284134d9SBarry Smith    Output Parameter:
955284134d9SBarry Smith .    A - the resulting matrix
956284134d9SBarry Smith 
9578e6c10adSSatish Balay    Level: advanced
9588e6c10adSSatish Balay 
959284134d9SBarry Smith    Notes: See MATIS for more details
9608cda6cd7SBarry Smith           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
9618cda6cd7SBarry Smith           by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
9628cda6cd7SBarry Smith           plus the ghost points to global indices.
963284134d9SBarry Smith 
964284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
965284134d9SBarry Smith @*/
9664e4c7dbeSStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
967284134d9SBarry Smith {
968284134d9SBarry Smith   PetscErrorCode ierr;
969284134d9SBarry Smith 
970284134d9SBarry Smith   PetscFunctionBegin;
971284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
972d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
973284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
974284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
9753b03a366Sstefano_zampini   ierr = MatSetUp(*A);CHKERRQ(ierr);
976784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr);
977284134d9SBarry Smith   PetscFunctionReturn(0);
978284134d9SBarry Smith }
979284134d9SBarry Smith 
980b4319ba4SBarry Smith /*MC
9816a818285SBarry Smith    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners, see PCBDDC.
982b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
983b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
984b4319ba4SBarry Smith    product is handled "implicitly".
985b4319ba4SBarry Smith 
986b4319ba4SBarry Smith    Operations Provided:
9876726f965SBarry Smith +  MatMult()
9882e74eeadSLisandro Dalcin .  MatMultAdd()
9892e74eeadSLisandro Dalcin .  MatMultTranspose()
9902e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
9916726f965SBarry Smith .  MatZeroEntries()
9926726f965SBarry Smith .  MatSetOption()
9932e74eeadSLisandro Dalcin .  MatZeroRows()
9946726f965SBarry Smith .  MatZeroRowsLocal()
9952e74eeadSLisandro Dalcin .  MatSetValues()
9966726f965SBarry Smith .  MatSetValuesLocal()
9972e74eeadSLisandro Dalcin .  MatScale()
9982e74eeadSLisandro Dalcin .  MatGetDiagonal()
9996726f965SBarry Smith -  MatSetLocalToGlobalMapping()
1000b4319ba4SBarry Smith 
1001b4319ba4SBarry Smith    Options Database Keys:
1002b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
1003b4319ba4SBarry Smith 
1004b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
1005b4319ba4SBarry Smith 
1006b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
1007b4319ba4SBarry Smith 
1008b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
1009b4319ba4SBarry Smith           MatISGetLocalMat()
1010b4319ba4SBarry Smith 
1011b4319ba4SBarry Smith   Level: advanced
1012b4319ba4SBarry Smith 
10136a818285SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), PCBDDC
1014b4319ba4SBarry Smith 
1015b4319ba4SBarry Smith M*/
1016b4319ba4SBarry Smith 
1017b4319ba4SBarry Smith #undef __FUNCT__
1018b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
10198cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
1020b4319ba4SBarry Smith {
1021dfbe8321SBarry Smith   PetscErrorCode ierr;
1022b4319ba4SBarry Smith   Mat_IS         *b;
1023b4319ba4SBarry Smith 
1024b4319ba4SBarry Smith   PetscFunctionBegin;
1025b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1026b4319ba4SBarry Smith   A->data = (void*)b;
1027b4319ba4SBarry Smith   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1028b4319ba4SBarry Smith 
1029b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
10302e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
10312e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
10322e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1033b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1034b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
10352e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
1036b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1037f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
10382e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1039b4319ba4SBarry Smith   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
1040b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1041b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1042b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
10436726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
10442e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
10452e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
10466726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
104769796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
104869796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1049ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
1050b4319ba4SBarry Smith 
105126283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
105226283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1053b4319ba4SBarry Smith 
1054b7ce53b6SStefano Zampini   /* zeros pointer */
1055b4319ba4SBarry Smith   b->A           = 0;
1056b4319ba4SBarry Smith   b->ctx         = 0;
1057b4319ba4SBarry Smith   b->x           = 0;
1058b4319ba4SBarry Smith   b->y           = 0;
1059b09366fdSStefano Zampini   b->mapping     = 0;
106028f4e0baSStefano Zampini   b->sf          = 0;
106128f4e0baSStefano Zampini   b->sf_rootdata = 0;
106228f4e0baSStefano Zampini   b->sf_leafdata = 0;
1063b7ce53b6SStefano Zampini 
1064b7ce53b6SStefano Zampini   /* special MATIS functions */
1065bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1066bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1067b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
10682e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
106917667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1070b4319ba4SBarry Smith   PetscFunctionReturn(0);
1071b4319ba4SBarry Smith }
1072