xref: /petsc/src/mat/impls/is/matis.c (revision 3f08860eea72c60ab5998520c0bd5a036e7e882e)
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);
313bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(B->rmap->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);
343bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->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"
41eb82efa4SStefano 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__
1323927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private"
1333927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
1343927de2eSStefano Zampini {
1353927de2eSStefano Zampini   Mat_IS          *matis = (Mat_IS*)(A->data);
1363927de2eSStefano Zampini   PetscInt        *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership;
137ecf5a873SStefano Zampini   const PetscInt  *global_indices_r,*global_indices_c;
1383927de2eSStefano Zampini   PetscInt        i,j,bs,rows,cols;
1393927de2eSStefano Zampini   PetscInt        lrows,lcols;
1403927de2eSStefano Zampini   PetscInt        local_rows,local_cols;
1413927de2eSStefano Zampini   PetscMPIInt     nsubdomains;
1423927de2eSStefano Zampini   PetscBool       isdense,issbaij;
1433927de2eSStefano Zampini   PetscErrorCode  ierr;
1443927de2eSStefano Zampini 
1453927de2eSStefano Zampini   PetscFunctionBegin;
1463927de2eSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr);
1473927de2eSStefano Zampini   ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr);
1483927de2eSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1493927de2eSStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
1503927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
1513927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
152ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
153ecf5a873SStefano Zampini   if (A->rmap->mapping != A->cmap->mapping) {
154ecf5a873SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_c);CHKERRQ(ierr);
155ecf5a873SStefano Zampini   } else {
156ecf5a873SStefano Zampini     global_indices_c = global_indices_r;
157ecf5a873SStefano Zampini   }
158ecf5a873SStefano Zampini 
1593927de2eSStefano Zampini   if (issbaij) {
1603927de2eSStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
1613927de2eSStefano Zampini   }
1623927de2eSStefano Zampini   /*
163ecf5a873SStefano Zampini      An SF reduce is needed to sum up properly on shared rows.
1643927de2eSStefano Zampini      Note that generally preallocation is not exact, since it overestimates nonzeros
1653927de2eSStefano Zampini   */
1663927de2eSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
1673927de2eSStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
1683927de2eSStefano Zampini   }
1693927de2eSStefano Zampini   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
1703927de2eSStefano Zampini   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
1713927de2eSStefano Zampini   /* All processes need to compute entire row ownership */
1723927de2eSStefano Zampini   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
1733927de2eSStefano Zampini   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
1743927de2eSStefano Zampini   for (i=0;i<nsubdomains;i++) {
1753927de2eSStefano Zampini     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
1763927de2eSStefano Zampini       row_ownership[j] = i;
1773927de2eSStefano Zampini     }
1783927de2eSStefano Zampini   }
1793927de2eSStefano Zampini 
1803927de2eSStefano Zampini   /*
1813927de2eSStefano Zampini      my_dnz and my_onz contains exact contribution to preallocation from each local mat
1823927de2eSStefano Zampini      then, they will be summed up properly. This way, preallocation is always sufficient
1833927de2eSStefano Zampini   */
1843927de2eSStefano Zampini   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
1853927de2eSStefano Zampini   /* preallocation as a MATAIJ */
1863927de2eSStefano Zampini   if (isdense) { /* special case for dense local matrices */
1873927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
188ecf5a873SStefano Zampini       PetscInt index_row = global_indices_r[i];
1893927de2eSStefano Zampini       for (j=i;j<local_rows;j++) {
1903927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
191ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[j];
1923927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
1933927de2eSStefano Zampini           my_dnz[i] += 1;
1943927de2eSStefano Zampini         } else { /* offdiag block */
1953927de2eSStefano Zampini           my_onz[i] += 1;
1963927de2eSStefano Zampini         }
1973927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
1983927de2eSStefano Zampini         if (i != j) {
1993927de2eSStefano Zampini           owner = row_ownership[index_col];
2003927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
2013927de2eSStefano Zampini             my_dnz[j] += 1;
2023927de2eSStefano Zampini           } else {
2033927de2eSStefano Zampini             my_onz[j] += 1;
2043927de2eSStefano Zampini           }
2053927de2eSStefano Zampini         }
2063927de2eSStefano Zampini       }
2073927de2eSStefano Zampini     }
208ecf5a873SStefano Zampini   } else { /* TODO: this could be optimized using MatGetRowIJ */
2093927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
2103927de2eSStefano Zampini       const PetscInt *cols;
211ecf5a873SStefano Zampini       PetscInt       ncols,index_row = global_indices_r[i];
2123927de2eSStefano Zampini       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2133927de2eSStefano Zampini       for (j=0;j<ncols;j++) {
2143927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
215ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[cols[j]];
2163927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
2173927de2eSStefano Zampini           my_dnz[i] += 1;
2183927de2eSStefano Zampini         } else { /* offdiag block */
2193927de2eSStefano Zampini           my_onz[i] += 1;
2203927de2eSStefano Zampini         }
2213927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
222d9a9e74cSStefano Zampini         if (issbaij && index_col != index_row) {
2233927de2eSStefano Zampini           owner = row_ownership[index_col];
2243927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
225d9a9e74cSStefano Zampini             my_dnz[cols[j]] += 1;
2263927de2eSStefano Zampini           } else {
227d9a9e74cSStefano Zampini             my_onz[cols[j]] += 1;
2283927de2eSStefano Zampini           }
2293927de2eSStefano Zampini         }
2303927de2eSStefano Zampini       }
2313927de2eSStefano Zampini       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2323927de2eSStefano Zampini     }
2333927de2eSStefano Zampini   }
234ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
235ecf5a873SStefano Zampini   if (global_indices_c != global_indices_r) {
236ecf5a873SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_c);CHKERRQ(ierr);
237ecf5a873SStefano Zampini   }
2383927de2eSStefano Zampini   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
239ecf5a873SStefano Zampini 
240ecf5a873SStefano Zampini   /* Reduce my_dnz and my_onz */
2413927de2eSStefano Zampini   if (maxreduce) {
2423927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2433927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2443927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2453927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2463927de2eSStefano Zampini   } else {
2473927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2483927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2493927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2503927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2513927de2eSStefano Zampini   }
2523927de2eSStefano Zampini   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
2533927de2eSStefano Zampini 
2543927de2eSStefano Zampini   /* Resize preallocation if overestimated */
2553927de2eSStefano Zampini   for (i=0;i<lrows;i++) {
2563927de2eSStefano Zampini     dnz[i] = PetscMin(dnz[i],lcols);
2573927de2eSStefano Zampini     onz[i] = PetscMin(onz[i],cols-lcols);
2583927de2eSStefano Zampini   }
2593927de2eSStefano Zampini   /* set preallocation */
2603927de2eSStefano Zampini   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
2613927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
2623927de2eSStefano Zampini     dnz[i] = dnz[i*bs]/bs;
2633927de2eSStefano Zampini     onz[i] = onz[i*bs]/bs;
2643927de2eSStefano Zampini   }
2653927de2eSStefano Zampini   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2663927de2eSStefano Zampini   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2673927de2eSStefano Zampini   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
2683927de2eSStefano Zampini   if (issbaij) {
2693927de2eSStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
2703927de2eSStefano Zampini   }
2713927de2eSStefano Zampini   PetscFunctionReturn(0);
2723927de2eSStefano Zampini }
2733927de2eSStefano Zampini 
2743927de2eSStefano Zampini #undef __FUNCT__
275b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
276b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
277b7ce53b6SStefano Zampini {
278b7ce53b6SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
279d9a9e74cSStefano Zampini   Mat            local_mat;
280b7ce53b6SStefano Zampini   /* info on mat */
2813cfa4ea4SStefano Zampini   PetscInt       bs,rows,cols,lrows,lcols;
282b7ce53b6SStefano Zampini   PetscInt       local_rows,local_cols;
283686e3a49SStefano Zampini   PetscBool      isdense,issbaij,isseqaij;
2847c03b4e8SStefano Zampini   PetscMPIInt    nsubdomains;
285b7ce53b6SStefano Zampini   /* values insertion */
286b7ce53b6SStefano Zampini   PetscScalar    *array;
287b7ce53b6SStefano Zampini   /* work */
288b7ce53b6SStefano Zampini   PetscErrorCode ierr;
289b7ce53b6SStefano Zampini 
290b7ce53b6SStefano Zampini   PetscFunctionBegin;
291b7ce53b6SStefano Zampini   /* get info from mat */
2927c03b4e8SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
2937c03b4e8SStefano Zampini   if (nsubdomains == 1) {
2947c03b4e8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
2957c03b4e8SStefano Zampini       ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr);
2967c03b4e8SStefano Zampini     } else {
2977c03b4e8SStefano Zampini       ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
2987c03b4e8SStefano Zampini     }
2997c03b4e8SStefano Zampini     PetscFunctionReturn(0);
3007c03b4e8SStefano Zampini   }
301b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
302b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
3033cfa4ea4SStefano Zampini   ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr);
304b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
305b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
306686e3a49SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
307b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
308b7ce53b6SStefano Zampini 
309b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
310b7ce53b6SStefano Zampini     MatType     new_mat_type;
3113927de2eSStefano Zampini     PetscBool   issbaij_red;
312b7ce53b6SStefano Zampini 
313b7ce53b6SStefano Zampini     /* determining new matrix type */
314b7ce53b6SStefano Zampini     ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
315b7ce53b6SStefano Zampini     if (issbaij_red) {
316b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
317b7ce53b6SStefano Zampini     } else {
318b7ce53b6SStefano Zampini       if (bs>1) {
319b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
320b7ce53b6SStefano Zampini       } else {
321b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
322b7ce53b6SStefano Zampini       }
323b7ce53b6SStefano Zampini     }
324b7ce53b6SStefano Zampini 
3253927de2eSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
3263cfa4ea4SStefano Zampini     ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr);
3273927de2eSStefano Zampini     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
3283927de2eSStefano Zampini     ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr);
3293927de2eSStefano Zampini     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
330b7ce53b6SStefano Zampini   } else {
3313cfa4ea4SStefano Zampini     PetscInt mbs,mrows,mcols,mlrows,mlcols;
332b7ce53b6SStefano Zampini     /* some checks */
333b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
334b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
3353cfa4ea4SStefano Zampini     ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr);
336b7ce53b6SStefano Zampini     if (mrows != rows) {
337b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
338b7ce53b6SStefano Zampini     }
3393cfa4ea4SStefano Zampini     if (mcols != cols) {
340b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
341b7ce53b6SStefano Zampini     }
3423cfa4ea4SStefano Zampini     if (mlrows != lrows) {
3433cfa4ea4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows);
3443cfa4ea4SStefano Zampini     }
3453cfa4ea4SStefano Zampini     if (mlcols != lcols) {
3463cfa4ea4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols);
3473cfa4ea4SStefano Zampini     }
348b7ce53b6SStefano Zampini     if (mbs != bs) {
349b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
350b7ce53b6SStefano Zampini     }
351b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
352b7ce53b6SStefano Zampini   }
353d9a9e74cSStefano Zampini 
354d9a9e74cSStefano Zampini   if (issbaij) {
355d9a9e74cSStefano Zampini     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr);
356d9a9e74cSStefano Zampini   } else {
357d9a9e74cSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
358d9a9e74cSStefano Zampini     local_mat = matis->A;
359d9a9e74cSStefano Zampini   }
360686e3a49SStefano Zampini 
361b7ce53b6SStefano Zampini   /* Set values */
362ecf5a873SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr);
363b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
364ecf5a873SStefano Zampini     PetscInt i,*dummy_rows,*dummy_cols;
365ecf5a873SStefano Zampini 
366ecf5a873SStefano Zampini     if (local_rows != local_cols) {
367ecf5a873SStefano Zampini       ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr);
368ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
369ecf5a873SStefano Zampini       for (i=0;i<local_cols;i++) dummy_cols[i] = i;
370ecf5a873SStefano Zampini     } else {
371ecf5a873SStefano Zampini       ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr);
372ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
373ecf5a873SStefano Zampini       dummy_cols = dummy_rows;
374ecf5a873SStefano Zampini     }
375b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
376d9a9e74cSStefano Zampini     ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr);
377ecf5a873SStefano Zampini     ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr);
378d9a9e74cSStefano Zampini     ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr);
379ecf5a873SStefano Zampini     if (dummy_rows != dummy_cols) {
380ecf5a873SStefano Zampini       ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr);
381ecf5a873SStefano Zampini     } else {
382ecf5a873SStefano Zampini       ierr = PetscFree(dummy_rows);CHKERRQ(ierr);
383ecf5a873SStefano Zampini     }
384686e3a49SStefano Zampini   } else if (isseqaij) {
385ecf5a873SStefano Zampini     PetscInt  i,nvtxs,*xadj,*adjncy;
386686e3a49SStefano Zampini     PetscBool done;
387686e3a49SStefano Zampini 
388d9a9e74cSStefano Zampini     ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
389686e3a49SStefano Zampini     if (!done) {
390d9a9e74cSStefano Zampini       SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
391b7ce53b6SStefano Zampini     }
392d9a9e74cSStefano Zampini     ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr);
393686e3a49SStefano Zampini     for (i=0;i<nvtxs;i++) {
394ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr);
395686e3a49SStefano Zampini     }
396d9a9e74cSStefano Zampini     ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
397686e3a49SStefano Zampini     if (!done) {
398d9a9e74cSStefano Zampini       SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
399686e3a49SStefano Zampini     }
400d9a9e74cSStefano Zampini     ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr);
401686e3a49SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
402ecf5a873SStefano Zampini     PetscInt i;
403c0962df8SStefano Zampini 
404686e3a49SStefano Zampini     for (i=0;i<local_rows;i++) {
405686e3a49SStefano Zampini       PetscInt       j;
406ecf5a873SStefano Zampini       const PetscInt *local_indices_cols;
407686e3a49SStefano Zampini 
408ecf5a873SStefano Zampini       ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
409ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr);
410ecf5a873SStefano Zampini       ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
411686e3a49SStefano Zampini     }
412b7ce53b6SStefano Zampini   }
413b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
414d9a9e74cSStefano Zampini   ierr = MatDestroy(&local_mat);CHKERRQ(ierr);
415b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
416b7ce53b6SStefano Zampini   if (isdense) {
417b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
418b7ce53b6SStefano Zampini   }
419b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
420b7ce53b6SStefano Zampini }
421b7ce53b6SStefano Zampini 
422b7ce53b6SStefano Zampini #undef __FUNCT__
423b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
424b7ce53b6SStefano Zampini /*@
425b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
426b7ce53b6SStefano Zampini 
427b7ce53b6SStefano Zampini   Input Parameter:
428b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
429b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
430b7ce53b6SStefano Zampini 
431b7ce53b6SStefano Zampini   Output Parameter:
432b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
433b7ce53b6SStefano Zampini 
434b7ce53b6SStefano Zampini   Level: developer
435b7ce53b6SStefano Zampini 
436eb82efa4SStefano Zampini   Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested.
437b7ce53b6SStefano Zampini 
438b7ce53b6SStefano Zampini .seealso: MATIS
439b7ce53b6SStefano Zampini @*/
440b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
441b7ce53b6SStefano Zampini {
442b7ce53b6SStefano Zampini   PetscErrorCode ierr;
443b7ce53b6SStefano Zampini 
444b7ce53b6SStefano Zampini   PetscFunctionBegin;
445b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
446b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
447b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
448b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
449b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
450b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
451eb82efa4SStefano Zampini     if (mat == *newmat) {
452eb82efa4SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix");
453eb82efa4SStefano Zampini     }
454b7ce53b6SStefano Zampini   }
455b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
456b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
457b7ce53b6SStefano Zampini }
458b7ce53b6SStefano Zampini 
459b7ce53b6SStefano Zampini #undef __FUNCT__
460ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
461ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
462ad6194a2SStefano Zampini {
463ad6194a2SStefano Zampini   PetscErrorCode ierr;
464ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
465ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
466ad6194a2SStefano Zampini   Mat            B,localmat;
467ad6194a2SStefano Zampini 
468ad6194a2SStefano Zampini   PetscFunctionBegin;
469ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
470ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
471ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
472e176bc59SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr);
473ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
474ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
475b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
476ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
477ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
478ad6194a2SStefano Zampini   *newmat = B;
479ad6194a2SStefano Zampini   PetscFunctionReturn(0);
480ad6194a2SStefano Zampini }
481ad6194a2SStefano Zampini 
482ad6194a2SStefano Zampini #undef __FUNCT__
48369796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
48469796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
48569796d55SStefano Zampini {
48669796d55SStefano Zampini   PetscErrorCode ierr;
48769796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
48869796d55SStefano Zampini   PetscBool      local_sym;
48969796d55SStefano Zampini 
49069796d55SStefano Zampini   PetscFunctionBegin;
49169796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
49269796d55SStefano Zampini   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
49369796d55SStefano Zampini   PetscFunctionReturn(0);
49469796d55SStefano Zampini }
49569796d55SStefano Zampini 
49669796d55SStefano Zampini #undef __FUNCT__
49769796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
49869796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
49969796d55SStefano Zampini {
50069796d55SStefano Zampini   PetscErrorCode ierr;
50169796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
50269796d55SStefano Zampini   PetscBool      local_sym;
50369796d55SStefano Zampini 
50469796d55SStefano Zampini   PetscFunctionBegin;
50569796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
50669796d55SStefano Zampini   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
50769796d55SStefano Zampini   PetscFunctionReturn(0);
50869796d55SStefano Zampini }
50969796d55SStefano Zampini 
51069796d55SStefano Zampini #undef __FUNCT__
511b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
512dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
513b4319ba4SBarry Smith {
514dfbe8321SBarry Smith   PetscErrorCode ierr;
515b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
516b4319ba4SBarry Smith 
517b4319ba4SBarry Smith   PetscFunctionBegin;
5186bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
519e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr);
520e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr);
5216bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
5226bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
52328f4e0baSStefano Zampini   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
52428f4e0baSStefano Zampini   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
525bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
526dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
527bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
528b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
529b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
5302e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
531b4319ba4SBarry Smith   PetscFunctionReturn(0);
532b4319ba4SBarry Smith }
533b4319ba4SBarry Smith 
534b4319ba4SBarry Smith #undef __FUNCT__
535b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
536dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
537b4319ba4SBarry Smith {
538dfbe8321SBarry Smith   PetscErrorCode ierr;
539b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
540b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
541b4319ba4SBarry Smith 
542b4319ba4SBarry Smith   PetscFunctionBegin;
543b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
544e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
545e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
546b4319ba4SBarry Smith 
547b4319ba4SBarry Smith   /* multiply the local matrix */
548b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
549b4319ba4SBarry Smith 
550b4319ba4SBarry Smith   /* scatter product back into global memory */
5512dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
552e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
553e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
554b4319ba4SBarry Smith   PetscFunctionReturn(0);
555b4319ba4SBarry Smith }
556b4319ba4SBarry Smith 
557b4319ba4SBarry Smith #undef __FUNCT__
5582e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
559b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5602e74eeadSLisandro Dalcin {
561650997f4SStefano Zampini   Vec            temp_vec;
5622e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5632e74eeadSLisandro Dalcin 
5642e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
565650997f4SStefano Zampini   if (v3 != v2) {
566650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
567650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
568650997f4SStefano Zampini   } else {
569650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
570650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
571650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
572650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
573650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
574650997f4SStefano Zampini   }
5752e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5762e74eeadSLisandro Dalcin }
5772e74eeadSLisandro Dalcin 
5782e74eeadSLisandro Dalcin #undef __FUNCT__
5792e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
580e176bc59SStefano Zampini PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x)
5812e74eeadSLisandro Dalcin {
5822e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
5832e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5842e74eeadSLisandro Dalcin 
585e176bc59SStefano Zampini   PetscFunctionBegin;
5862e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
587e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
588e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5892e74eeadSLisandro Dalcin 
5902e74eeadSLisandro Dalcin   /* multiply the local matrix */
591e176bc59SStefano Zampini   ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr);
5922e74eeadSLisandro Dalcin 
5932e74eeadSLisandro Dalcin   /* scatter product back into global vector */
594e176bc59SStefano Zampini   ierr = VecSet(x,0);CHKERRQ(ierr);
595e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
596e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5972e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5982e74eeadSLisandro Dalcin }
5992e74eeadSLisandro Dalcin 
6002e74eeadSLisandro Dalcin #undef __FUNCT__
6012e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
6022e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
6032e74eeadSLisandro Dalcin {
604650997f4SStefano Zampini   Vec            temp_vec;
6052e74eeadSLisandro Dalcin   PetscErrorCode ierr;
6062e74eeadSLisandro Dalcin 
6072e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
608650997f4SStefano Zampini   if (v3 != v2) {
609650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
610650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
611650997f4SStefano Zampini   } else {
612650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
613650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
614650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
615650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
616650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
617650997f4SStefano Zampini   }
6182e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
6192e74eeadSLisandro Dalcin }
6202e74eeadSLisandro Dalcin 
6212e74eeadSLisandro Dalcin #undef __FUNCT__
622b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
623dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
624b4319ba4SBarry Smith {
625b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
626dfbe8321SBarry Smith   PetscErrorCode ierr;
627b4319ba4SBarry Smith   PetscViewer    sviewer;
628b4319ba4SBarry Smith 
629b4319ba4SBarry Smith   PetscFunctionBegin;
630dd2fa690SBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
631*3f08860eSBarry Smith   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
632b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
633*3f08860eSBarry Smith   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
634b4319ba4SBarry Smith   PetscFunctionReturn(0);
635b4319ba4SBarry Smith }
636b4319ba4SBarry Smith 
637b4319ba4SBarry Smith #undef __FUNCT__
638b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
639784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
640b4319ba4SBarry Smith {
641dfbe8321SBarry Smith   PetscErrorCode ierr;
642e176bc59SStefano Zampini   PetscInt       nr,rbs,nc,cbs;
643b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
644b4319ba4SBarry Smith   IS             from,to;
645e176bc59SStefano Zampini   Vec            cglobal,rglobal;
646b4319ba4SBarry Smith 
647b4319ba4SBarry Smith   PetscFunctionBegin;
648784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
649e176bc59SStefano Zampini   PetscCheckSameComm(A,1,cmapping,3);
6503bbff08aSStefano Zampini   /* Destroy any previous data */
65170cf5478SStefano Zampini   ierr = VecDestroy(&is->x);CHKERRQ(ierr);
65270cf5478SStefano Zampini   ierr = VecDestroy(&is->y);CHKERRQ(ierr);
653e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr);
654e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr);
6551c47cb0fSStefano Zampini   ierr = MatDestroy(&is->A);CHKERRQ(ierr);
65628f4e0baSStefano Zampini   ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
65728f4e0baSStefano Zampini   ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
6583bbff08aSStefano Zampini 
6593bbff08aSStefano Zampini   /* Setup Layout and set local to global maps */
660fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
661fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
662fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
663fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
664b4319ba4SBarry Smith 
665b4319ba4SBarry Smith   /* Create the local matrix A */
666e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr);
667e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr);
668e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr);
669e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr);
670f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
671e176bc59SStefano Zampini   ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr);
672e176bc59SStefano Zampini   ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr);
673e176bc59SStefano Zampini   ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr);
674ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
675ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
676b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
677b4319ba4SBarry Smith 
678b4319ba4SBarry Smith   /* Create the local work vectors */
6793bbff08aSStefano Zampini   ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr);
680b4319ba4SBarry Smith 
681e176bc59SStefano Zampini   /* setup the global to local scatters */
682e176bc59SStefano Zampini   ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr);
683e176bc59SStefano Zampini   ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr);
684784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
685e176bc59SStefano Zampini   ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr);
686e176bc59SStefano Zampini   if (rmapping != cmapping) {
687e176bc59SStefano Zampini     ierr = ISDestroy(&to);CHKERRQ(ierr);
688e176bc59SStefano Zampini     ierr = ISDestroy(&from);CHKERRQ(ierr);
689e176bc59SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr);
690e176bc59SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr);
691e176bc59SStefano Zampini     ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr);
692e176bc59SStefano Zampini   } else {
693e176bc59SStefano Zampini     ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr);
694e176bc59SStefano Zampini     is->cctx = is->rctx;
695e176bc59SStefano Zampini   }
696e176bc59SStefano Zampini   ierr = VecDestroy(&rglobal);CHKERRQ(ierr);
697e176bc59SStefano Zampini   ierr = VecDestroy(&cglobal);CHKERRQ(ierr);
6986bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
6996bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
700b4319ba4SBarry Smith   PetscFunctionReturn(0);
701b4319ba4SBarry Smith }
702b4319ba4SBarry Smith 
7032e74eeadSLisandro Dalcin #undef __FUNCT__
7042e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
7052e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
7062e74eeadSLisandro Dalcin {
7072e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
7082e74eeadSLisandro Dalcin   PetscInt       rows_l[2048],cols_l[2048];
7092e74eeadSLisandro Dalcin   PetscErrorCode ierr;
7102e74eeadSLisandro Dalcin 
7112e74eeadSLisandro Dalcin   PetscFunctionBegin;
7122e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
713f23aa3ddSBarry 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);
7142e74eeadSLisandro Dalcin #endif
7153bbff08aSStefano Zampini   ierr = ISG2LMapApply(mat->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr);
7163bbff08aSStefano Zampini   ierr = ISG2LMapApply(mat->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr);
7172e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
7182e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7192e74eeadSLisandro Dalcin }
7202e74eeadSLisandro Dalcin 
721b4319ba4SBarry Smith #undef __FUNCT__
722b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
72313f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
724b4319ba4SBarry Smith {
725dfbe8321SBarry Smith   PetscErrorCode ierr;
726b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
727b4319ba4SBarry Smith 
728b4319ba4SBarry Smith   PetscFunctionBegin;
729b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
730b4319ba4SBarry Smith   PetscFunctionReturn(0);
731b4319ba4SBarry Smith }
732b4319ba4SBarry Smith 
733b4319ba4SBarry Smith #undef __FUNCT__
734f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
735f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
736f0006bf2SLisandro Dalcin {
737f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
738f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
739f0006bf2SLisandro Dalcin 
740f0006bf2SLisandro Dalcin   PetscFunctionBegin;
741f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
742f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
743f0006bf2SLisandro Dalcin }
744f0006bf2SLisandro Dalcin 
745f0006bf2SLisandro Dalcin #undef __FUNCT__
7462e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
7472b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
7482e74eeadSLisandro Dalcin {
7490298fd71SBarry Smith   PetscInt       n_l = 0, *rows_l = NULL;
7502e74eeadSLisandro Dalcin   PetscErrorCode ierr;
7512e74eeadSLisandro Dalcin 
7522e74eeadSLisandro Dalcin   PetscFunctionBegin;
753ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
7542e74eeadSLisandro Dalcin   if (n) {
755785e854fSJed Brown     ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr);
7563bbff08aSStefano Zampini     ierr = ISGlobalToLocalMappingApply(A->rmap->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
7572e74eeadSLisandro Dalcin   }
7582b40b63fSBarry Smith   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
759c31cb41cSBarry Smith   ierr = PetscFree(rows_l);CHKERRQ(ierr);
7602e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7612e74eeadSLisandro Dalcin }
7622e74eeadSLisandro Dalcin 
7632e74eeadSLisandro Dalcin #undef __FUNCT__
764b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
7652b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
766b4319ba4SBarry Smith {
767b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
768dfbe8321SBarry Smith   PetscErrorCode ierr;
769f4df32b1SMatthew Knepley   PetscInt       i;
770b4319ba4SBarry Smith   PetscScalar    *array;
771b4319ba4SBarry Smith 
772b4319ba4SBarry Smith   PetscFunctionBegin;
773ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
774b4319ba4SBarry Smith   {
775b4319ba4SBarry Smith     /*
776b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
777b4319ba4SBarry Smith        work properly in the interface nodes.
778b4319ba4SBarry Smith     */
779b4319ba4SBarry Smith     Vec counter;
780e176bc59SStefano Zampini     ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr);
781e176bc59SStefano Zampini     ierr = VecSet(counter,0.);CHKERRQ(ierr);
782e176bc59SStefano Zampini     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
783e176bc59SStefano Zampini     ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
784e176bc59SStefano Zampini     ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
785e176bc59SStefano Zampini     ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
786e176bc59SStefano Zampini     ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7876bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
788b4319ba4SBarry Smith   }
789958c9bccSBarry Smith   if (!n) {
790b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
791b4319ba4SBarry Smith   } else {
792b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
7932205254eSKarl Rupp 
794e176bc59SStefano Zampini     ierr = VecGetArray(is->y,&array);CHKERRQ(ierr);
7952b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
796b4319ba4SBarry Smith     for (i=0; i<n; i++) {
797f4df32b1SMatthew Knepley       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
798b4319ba4SBarry Smith     }
799b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
800b4319ba4SBarry Smith     ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
801e176bc59SStefano Zampini     ierr = VecRestoreArray(is->y,&array);CHKERRQ(ierr);
802b4319ba4SBarry Smith   }
803b4319ba4SBarry Smith   PetscFunctionReturn(0);
804b4319ba4SBarry Smith }
805b4319ba4SBarry Smith 
806b4319ba4SBarry Smith #undef __FUNCT__
807b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
808dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
809b4319ba4SBarry Smith {
810b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
811dfbe8321SBarry Smith   PetscErrorCode ierr;
812b4319ba4SBarry Smith 
813b4319ba4SBarry Smith   PetscFunctionBegin;
814b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
815b4319ba4SBarry Smith   PetscFunctionReturn(0);
816b4319ba4SBarry Smith }
817b4319ba4SBarry Smith 
818b4319ba4SBarry Smith #undef __FUNCT__
819b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
820dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
821b4319ba4SBarry Smith {
822b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
823dfbe8321SBarry Smith   PetscErrorCode ierr;
824b4319ba4SBarry Smith 
825b4319ba4SBarry Smith   PetscFunctionBegin;
826b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
827b4319ba4SBarry Smith   PetscFunctionReturn(0);
828b4319ba4SBarry Smith }
829b4319ba4SBarry Smith 
830b4319ba4SBarry Smith #undef __FUNCT__
831b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
8327087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
833b4319ba4SBarry Smith {
834b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
835b4319ba4SBarry Smith 
836b4319ba4SBarry Smith   PetscFunctionBegin;
837b4319ba4SBarry Smith   *local = is->A;
838b4319ba4SBarry Smith   PetscFunctionReturn(0);
839b4319ba4SBarry Smith }
840b4319ba4SBarry Smith 
841b4319ba4SBarry Smith #undef __FUNCT__
842b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
843b4319ba4SBarry Smith /*@
844b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
845b4319ba4SBarry Smith 
846b4319ba4SBarry Smith   Input Parameter:
847b4319ba4SBarry Smith .  mat - the matrix
848b4319ba4SBarry Smith 
849b4319ba4SBarry Smith   Output Parameter:
850eb82efa4SStefano Zampini .  local - the local matrix
851b4319ba4SBarry Smith 
852b4319ba4SBarry Smith   Level: advanced
853b4319ba4SBarry Smith 
854b4319ba4SBarry Smith   Notes:
855b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
856b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
857b4319ba4SBarry Smith   of the MatSetValues() operation.
858b4319ba4SBarry Smith 
859b4319ba4SBarry Smith .seealso: MATIS
860b4319ba4SBarry Smith @*/
8617087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
862b4319ba4SBarry Smith {
8634ac538c5SBarry Smith   PetscErrorCode ierr;
864b4319ba4SBarry Smith 
865b4319ba4SBarry Smith   PetscFunctionBegin;
8660700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
867b4319ba4SBarry Smith   PetscValidPointer(local,2);
8684ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
869b4319ba4SBarry Smith   PetscFunctionReturn(0);
870b4319ba4SBarry Smith }
871b4319ba4SBarry Smith 
8723b03a366Sstefano_zampini #undef __FUNCT__
8733b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
8743b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
8753b03a366Sstefano_zampini {
8763b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
8773b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
8783b03a366Sstefano_zampini   PetscErrorCode ierr;
8793b03a366Sstefano_zampini 
8803b03a366Sstefano_zampini   PetscFunctionBegin;
8814e4c7dbeSStefano Zampini   if (is->A) {
8823b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
8833b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
884f23aa3ddSBarry 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);
8854e4c7dbeSStefano Zampini   }
8863b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
8873b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
8883b03a366Sstefano_zampini   is->A = local;
8893b03a366Sstefano_zampini   PetscFunctionReturn(0);
8903b03a366Sstefano_zampini }
8913b03a366Sstefano_zampini 
8923b03a366Sstefano_zampini #undef __FUNCT__
8933b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
8943b03a366Sstefano_zampini /*@
895eb82efa4SStefano Zampini     MatISSetLocalMat - Replace the local matrix stored inside a MATIS object.
8963b03a366Sstefano_zampini 
8973b03a366Sstefano_zampini   Input Parameter:
8983b03a366Sstefano_zampini .  mat - the matrix
899eb82efa4SStefano Zampini .  local - the local matrix
9003b03a366Sstefano_zampini 
9013b03a366Sstefano_zampini   Output Parameter:
9023b03a366Sstefano_zampini 
9033b03a366Sstefano_zampini   Level: advanced
9043b03a366Sstefano_zampini 
9053b03a366Sstefano_zampini   Notes:
9063b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
9073b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
9083b03a366Sstefano_zampini 
9093b03a366Sstefano_zampini .seealso: MATIS
9103b03a366Sstefano_zampini @*/
9113b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
9123b03a366Sstefano_zampini {
9133b03a366Sstefano_zampini   PetscErrorCode ierr;
9143b03a366Sstefano_zampini 
9153b03a366Sstefano_zampini   PetscFunctionBegin;
9163b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
917b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
9183b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
9193b03a366Sstefano_zampini   PetscFunctionReturn(0);
9203b03a366Sstefano_zampini }
9213b03a366Sstefano_zampini 
9226726f965SBarry Smith #undef __FUNCT__
9236726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
9246726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
9256726f965SBarry Smith {
9266726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
9276726f965SBarry Smith   PetscErrorCode ierr;
9286726f965SBarry Smith 
9296726f965SBarry Smith   PetscFunctionBegin;
9306726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
9316726f965SBarry Smith   PetscFunctionReturn(0);
9326726f965SBarry Smith }
9336726f965SBarry Smith 
9346726f965SBarry Smith #undef __FUNCT__
9352e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
9362e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
9372e74eeadSLisandro Dalcin {
9382e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9392e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9402e74eeadSLisandro Dalcin 
9412e74eeadSLisandro Dalcin   PetscFunctionBegin;
9422e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
9432e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9442e74eeadSLisandro Dalcin }
9452e74eeadSLisandro Dalcin 
9462e74eeadSLisandro Dalcin #undef __FUNCT__
9472e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
9482e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
9492e74eeadSLisandro Dalcin {
9502e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9512e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9522e74eeadSLisandro Dalcin 
9532e74eeadSLisandro Dalcin   PetscFunctionBegin;
9542e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
955e176bc59SStefano Zampini   ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr);
9562e74eeadSLisandro Dalcin 
9572e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
9582e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
959e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
960e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9612e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9622e74eeadSLisandro Dalcin }
9632e74eeadSLisandro Dalcin 
9642e74eeadSLisandro Dalcin #undef __FUNCT__
9656726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
966ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
9676726f965SBarry Smith {
9686726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
9696726f965SBarry Smith   PetscErrorCode ierr;
9706726f965SBarry Smith 
9716726f965SBarry Smith   PetscFunctionBegin;
9724e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
9736726f965SBarry Smith   PetscFunctionReturn(0);
9746726f965SBarry Smith }
9756726f965SBarry Smith 
976284134d9SBarry Smith #undef __FUNCT__
977284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
978284134d9SBarry Smith /*@
979284134d9SBarry Smith     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
980284134d9SBarry Smith        process but not across processes.
981284134d9SBarry Smith 
982284134d9SBarry Smith    Input Parameters:
983284134d9SBarry Smith +     comm    - MPI communicator that will share the matrix
984e176bc59SStefano Zampini .     bs      - block size of the matrix
985df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
986e176bc59SStefano Zampini .     rmap    - local to global map for rows
987e176bc59SStefano Zampini -     cmap    - local to global map for cols
988284134d9SBarry Smith 
989284134d9SBarry Smith    Output Parameter:
990284134d9SBarry Smith .    A - the resulting matrix
991284134d9SBarry Smith 
9928e6c10adSSatish Balay    Level: advanced
9938e6c10adSSatish Balay 
994284134d9SBarry Smith    Notes: See MATIS for more details
9958cda6cd7SBarry Smith           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
996e176bc59SStefano Zampini           by that process. The sizes of rmap and cmap define the size of the local matrices.
997e176bc59SStefano Zampini           If either rmap or cmap are NULL, than the matrix is assumed to be square
998284134d9SBarry Smith 
999284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
1000284134d9SBarry Smith @*/
1001e176bc59SStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A)
1002284134d9SBarry Smith {
1003284134d9SBarry Smith   PetscErrorCode ierr;
1004284134d9SBarry Smith 
1005284134d9SBarry Smith   PetscFunctionBegin;
1006e176bc59SStefano Zampini   if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping");
1007284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1008d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
1009284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1010284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
10113b03a366Sstefano_zampini   ierr = MatSetUp(*A);CHKERRQ(ierr);
1012e176bc59SStefano Zampini   if (rmap && cmap) {
1013e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr);
1014e176bc59SStefano Zampini   } else if (!rmap) {
1015e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr);
1016e176bc59SStefano Zampini   } else {
1017e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr);
1018e176bc59SStefano Zampini   }
1019284134d9SBarry Smith   PetscFunctionReturn(0);
1020284134d9SBarry Smith }
1021284134d9SBarry Smith 
1022b4319ba4SBarry Smith /*MC
1023eb82efa4SStefano Zampini    MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC).
1024b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
1025b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
1026b4319ba4SBarry Smith    product is handled "implicitly".
1027b4319ba4SBarry Smith 
1028b4319ba4SBarry Smith    Operations Provided:
10296726f965SBarry Smith +  MatMult()
10302e74eeadSLisandro Dalcin .  MatMultAdd()
10312e74eeadSLisandro Dalcin .  MatMultTranspose()
10322e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
10336726f965SBarry Smith .  MatZeroEntries()
10346726f965SBarry Smith .  MatSetOption()
10352e74eeadSLisandro Dalcin .  MatZeroRows()
10366726f965SBarry Smith .  MatZeroRowsLocal()
10372e74eeadSLisandro Dalcin .  MatSetValues()
10386726f965SBarry Smith .  MatSetValuesLocal()
10392e74eeadSLisandro Dalcin .  MatScale()
10402e74eeadSLisandro Dalcin .  MatGetDiagonal()
10416726f965SBarry Smith -  MatSetLocalToGlobalMapping()
1042b4319ba4SBarry Smith 
1043b4319ba4SBarry Smith    Options Database Keys:
1044b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
1045b4319ba4SBarry Smith 
1046b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
1047b4319ba4SBarry Smith 
1048b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
1049b4319ba4SBarry Smith 
1050b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
1051eb82efa4SStefano Zampini           MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation()
1052b4319ba4SBarry Smith 
1053b4319ba4SBarry Smith   Level: advanced
1054b4319ba4SBarry Smith 
1055eb82efa4SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC
1056b4319ba4SBarry Smith 
1057b4319ba4SBarry Smith M*/
1058b4319ba4SBarry Smith 
1059b4319ba4SBarry Smith #undef __FUNCT__
1060b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
10618cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
1062b4319ba4SBarry Smith {
1063dfbe8321SBarry Smith   PetscErrorCode ierr;
1064b4319ba4SBarry Smith   Mat_IS         *b;
1065b4319ba4SBarry Smith 
1066b4319ba4SBarry Smith   PetscFunctionBegin;
1067b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1068b4319ba4SBarry Smith   A->data = (void*)b;
1069b4319ba4SBarry Smith 
1070e176bc59SStefano Zampini   /* matrix ops */
1071e176bc59SStefano Zampini   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1072b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
10732e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
10742e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
10752e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1076b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1077b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
10782e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
1079b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1080f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
10812e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1082b4319ba4SBarry Smith   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
1083b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1084b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1085b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
10866726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
10872e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
10882e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
10896726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
109069796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
109169796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1092ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
1093b4319ba4SBarry Smith 
1094b7ce53b6SStefano Zampini   /* special MATIS functions */
1095bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1096bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1097b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
10982e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
109917667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1100b4319ba4SBarry Smith   PetscFunctionReturn(0);
1101b4319ba4SBarry Smith }
1102