xref: /petsc/src/mat/impls/is/matis.c (revision 3cfa4ea460eed5c3d3ab4210924c568cabc7705a)
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__
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;
1373927de2eSStefano Zampini   const PetscInt* global_indices;
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);
1523927de2eSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&global_indices);CHKERRQ(ierr);
1533927de2eSStefano Zampini   if (issbaij) {
1543927de2eSStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
1553927de2eSStefano Zampini   }
1563927de2eSStefano Zampini   /*
1573927de2eSStefano Zampini      An SF reduce is needed to sum up properly on shared interface dofs.
1583927de2eSStefano Zampini      Note that generally preallocation is not exact, since it overestimates nonzeros
1593927de2eSStefano Zampini   */
1603927de2eSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
1613927de2eSStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
1623927de2eSStefano Zampini   }
1633927de2eSStefano Zampini   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
1643927de2eSStefano Zampini   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
1653927de2eSStefano Zampini   /* All processes need to compute entire row ownership */
1663927de2eSStefano Zampini   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
1673927de2eSStefano Zampini   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
1683927de2eSStefano Zampini   for (i=0;i<nsubdomains;i++) {
1693927de2eSStefano Zampini     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
1703927de2eSStefano Zampini       row_ownership[j]=i;
1713927de2eSStefano Zampini     }
1723927de2eSStefano Zampini   }
1733927de2eSStefano Zampini 
1743927de2eSStefano Zampini   /*
1753927de2eSStefano Zampini      my_dnz and my_onz contains exact contribution to preallocation from each local mat
1763927de2eSStefano Zampini      then, they will be summed up properly. This way, preallocation is always sufficient
1773927de2eSStefano Zampini   */
1783927de2eSStefano Zampini   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
1793927de2eSStefano Zampini   /* preallocation as a MATAIJ */
1803927de2eSStefano Zampini   if (isdense) { /* special case for dense local matrices */
1813927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
1823927de2eSStefano Zampini       PetscInt index_row = global_indices[i];
1833927de2eSStefano Zampini       for (j=i;j<local_rows;j++) {
1843927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
1853927de2eSStefano Zampini         PetscInt index_col = global_indices[j];
1863927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
1873927de2eSStefano Zampini           my_dnz[i] += 1;
1883927de2eSStefano Zampini         } else { /* offdiag block */
1893927de2eSStefano Zampini           my_onz[i] += 1;
1903927de2eSStefano Zampini         }
1913927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
1923927de2eSStefano Zampini         if (i != j) {
1933927de2eSStefano Zampini           owner = row_ownership[index_col];
1943927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
1953927de2eSStefano Zampini             my_dnz[j] += 1;
1963927de2eSStefano Zampini           } else {
1973927de2eSStefano Zampini             my_onz[j] += 1;
1983927de2eSStefano Zampini           }
1993927de2eSStefano Zampini         }
2003927de2eSStefano Zampini       }
2013927de2eSStefano Zampini     }
2023927de2eSStefano Zampini   } else {
2033927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
2043927de2eSStefano Zampini       const PetscInt *cols;
2053927de2eSStefano Zampini       PetscInt       ncols,index_row = global_indices[i];
2063927de2eSStefano Zampini       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2073927de2eSStefano Zampini       for (j=0;j<ncols;j++) {
2083927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
2093927de2eSStefano Zampini         PetscInt index_col = global_indices[cols[j]];
2103927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
2113927de2eSStefano Zampini           my_dnz[i] += 1;
2123927de2eSStefano Zampini         } else { /* offdiag block */
2133927de2eSStefano Zampini           my_onz[i] += 1;
2143927de2eSStefano Zampini         }
2153927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
2163927de2eSStefano Zampini         if (issbaij) {
2173927de2eSStefano Zampini           owner = row_ownership[index_col];
2183927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
2193927de2eSStefano Zampini             my_dnz[j] += 1;
2203927de2eSStefano Zampini           } else {
2213927de2eSStefano Zampini             my_onz[j] += 1;
2223927de2eSStefano Zampini           }
2233927de2eSStefano Zampini         }
2243927de2eSStefano Zampini       }
2253927de2eSStefano Zampini       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2263927de2eSStefano Zampini     }
2273927de2eSStefano Zampini   }
2283927de2eSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&global_indices);CHKERRQ(ierr);
2293927de2eSStefano Zampini   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
2303927de2eSStefano Zampini   if (maxreduce) {
2313927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2323927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2333927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2343927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2353927de2eSStefano Zampini   } else {
2363927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2373927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2383927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2393927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2403927de2eSStefano Zampini   }
2413927de2eSStefano Zampini   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
2423927de2eSStefano Zampini 
2433927de2eSStefano Zampini   /* Resize preallocation if overestimated */
2443927de2eSStefano Zampini   for (i=0;i<lrows;i++) {
2453927de2eSStefano Zampini     dnz[i] = PetscMin(dnz[i],lcols);
2463927de2eSStefano Zampini     onz[i] = PetscMin(onz[i],cols-lcols);
2473927de2eSStefano Zampini   }
2483927de2eSStefano Zampini   /* set preallocation */
2493927de2eSStefano Zampini   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
2503927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
2513927de2eSStefano Zampini     dnz[i] = dnz[i*bs]/bs;
2523927de2eSStefano Zampini     onz[i] = onz[i*bs]/bs;
2533927de2eSStefano Zampini   }
2543927de2eSStefano Zampini   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2553927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
2563927de2eSStefano Zampini     dnz[i] = dnz[i]-i;
2573927de2eSStefano Zampini   }
2583927de2eSStefano Zampini   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2593927de2eSStefano Zampini   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
2603927de2eSStefano Zampini   if (issbaij) {
2613927de2eSStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
2623927de2eSStefano Zampini   }
2633927de2eSStefano Zampini   PetscFunctionReturn(0);
2643927de2eSStefano Zampini }
2653927de2eSStefano Zampini 
2663927de2eSStefano 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; */
273*3cfa4ea4SStefano Zampini   PetscInt               bs,rows,cols,lrows,lcols;
274b7ce53b6SStefano Zampini   PetscInt               local_rows,local_cols;
275686e3a49SStefano Zampini   PetscBool              isdense,issbaij,isseqaij;
276686e3a49SStefano Zampini   const PetscInt         *global_indices_rows;
2777c03b4e8SStefano Zampini   PetscMPIInt            nsubdomains;
278b7ce53b6SStefano Zampini   /* values insertion */
279b7ce53b6SStefano Zampini   PetscScalar            *array;
280b7ce53b6SStefano Zampini   /* work */
281b7ce53b6SStefano Zampini   PetscErrorCode         ierr;
282b7ce53b6SStefano Zampini 
283b7ce53b6SStefano Zampini   PetscFunctionBegin;
284b7ce53b6SStefano Zampini   /* MISSING CHECKS
285b7ce53b6SStefano Zampini     - rectangular case not covered (it is not allowed by MATIS)
286b7ce53b6SStefano Zampini   */
287b7ce53b6SStefano Zampini   /* get info from mat */
2887c03b4e8SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
2897c03b4e8SStefano Zampini   if (nsubdomains == 1) {
2907c03b4e8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
2917c03b4e8SStefano Zampini       ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr);
2927c03b4e8SStefano Zampini     } else {
2937c03b4e8SStefano Zampini       ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
2947c03b4e8SStefano Zampini     }
2957c03b4e8SStefano Zampini     PetscFunctionReturn(0);
2967c03b4e8SStefano Zampini   }
297b7ce53b6SStefano Zampini   /* ierr = MatGetLocalToGlobalMapping(mat,&rmapping,&cmapping);CHKERRQ(ierr); */
298b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
299b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
300*3cfa4ea4SStefano Zampini   ierr = MatGetLocalSize(mat,&lrows,&lcols);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);
303686e3a49SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
304b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
305b7ce53b6SStefano Zampini 
306686e3a49SStefano Zampini   /* map indices of local mat rows to global values */
307686e3a49SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&global_indices_rows);CHKERRQ(ierr);
308b7ce53b6SStefano Zampini 
309b7ce53b6SStefano Zampini   if (issbaij) {
310b7ce53b6SStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
311b7ce53b6SStefano Zampini   }
312b7ce53b6SStefano Zampini 
313b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
314b7ce53b6SStefano Zampini     MatType     new_mat_type;
3153927de2eSStefano Zampini     PetscBool   issbaij_red;
316b7ce53b6SStefano Zampini 
317b7ce53b6SStefano Zampini     /* determining new matrix type */
318b7ce53b6SStefano Zampini     ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
319b7ce53b6SStefano Zampini     if (issbaij_red) {
320b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
321b7ce53b6SStefano Zampini     } else {
322b7ce53b6SStefano Zampini       if (bs>1) {
323b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
324b7ce53b6SStefano Zampini       } else {
325b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
326b7ce53b6SStefano Zampini       }
327b7ce53b6SStefano Zampini     }
328b7ce53b6SStefano Zampini 
3293927de2eSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
330*3cfa4ea4SStefano Zampini     ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr);
3313927de2eSStefano Zampini     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
3323927de2eSStefano Zampini     ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr);
3333927de2eSStefano Zampini     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
334b7ce53b6SStefano Zampini   } else {
335*3cfa4ea4SStefano Zampini     PetscInt mbs,mrows,mcols,mlrows,mlcols;
336b7ce53b6SStefano Zampini     /* some checks */
337b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
338b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
339*3cfa4ea4SStefano Zampini     ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr);
340b7ce53b6SStefano Zampini     if (mrows != rows) {
341b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
342b7ce53b6SStefano Zampini     }
343*3cfa4ea4SStefano Zampini     if (mcols != cols) {
344b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
345b7ce53b6SStefano Zampini     }
346*3cfa4ea4SStefano Zampini     if (mlrows != lrows) {
347*3cfa4ea4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows);
348*3cfa4ea4SStefano Zampini     }
349*3cfa4ea4SStefano Zampini     if (mlcols != lcols) {
350*3cfa4ea4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols);
351*3cfa4ea4SStefano Zampini     }
352b7ce53b6SStefano Zampini     if (mbs != bs) {
353b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
354b7ce53b6SStefano Zampini     }
355b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
356b7ce53b6SStefano Zampini   }
357b7ce53b6SStefano Zampini   /* set local to global mappings */
358c0962df8SStefano Zampini   /* ierr = MatSetLocalToGlobalMapping(*M,matis->mapping,matis->mapping);CHKERRQ(ierr); */
359686e3a49SStefano Zampini 
360b7ce53b6SStefano Zampini   /* Set values */
361b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
362b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
363b7ce53b6SStefano Zampini     ierr = MatDenseGetArray(matis->A,&array);CHKERRQ(ierr);
364686e3a49SStefano Zampini     ierr = MatSetValues(*M,local_rows,global_indices_rows,local_cols,global_indices_rows,array,ADD_VALUES);CHKERRQ(ierr);
365b7ce53b6SStefano Zampini     ierr = MatDenseRestoreArray(matis->A,&array);CHKERRQ(ierr);
366686e3a49SStefano Zampini   } else if (isseqaij) {
367686e3a49SStefano Zampini     PetscInt  i,nvtxs,*xadj,*adjncy,*global_indices_cols;
368686e3a49SStefano Zampini     PetscBool done;
369686e3a49SStefano Zampini 
370686e3a49SStefano Zampini     ierr = MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
371686e3a49SStefano Zampini     if (!done) {
372686e3a49SStefano Zampini       SETERRQ1(PetscObjectComm((PetscObject)matis->A),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
373b7ce53b6SStefano Zampini     }
374686e3a49SStefano Zampini     ierr = MatSeqAIJGetArray(matis->A,&array);CHKERRQ(ierr);
375686e3a49SStefano Zampini     ierr = PetscMalloc1(xadj[nvtxs],&global_indices_cols);CHKERRQ(ierr);
376686e3a49SStefano Zampini     ierr = ISLocalToGlobalMappingApply(matis->mapping,xadj[nvtxs],adjncy,global_indices_cols);CHKERRQ(ierr);
377686e3a49SStefano Zampini     for (i=0;i<nvtxs;i++) {
378686e3a49SStefano Zampini       PetscInt    row,*cols,ncols;
379686e3a49SStefano Zampini       PetscScalar *mat_vals;
380686e3a49SStefano Zampini 
381686e3a49SStefano Zampini       row = global_indices_rows[i];
382686e3a49SStefano Zampini       ncols = xadj[i+1]-xadj[i];
383686e3a49SStefano Zampini       cols = global_indices_cols+xadj[i];
384686e3a49SStefano Zampini       mat_vals = array+xadj[i];
385686e3a49SStefano Zampini       ierr = MatSetValues(*M,1,&row,ncols,cols,mat_vals,ADD_VALUES);CHKERRQ(ierr);
386686e3a49SStefano Zampini     }
387686e3a49SStefano Zampini     ierr = MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
388686e3a49SStefano Zampini     if (!done) {
389686e3a49SStefano Zampini       SETERRQ1(PetscObjectComm((PetscObject)matis->A),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
390686e3a49SStefano Zampini     }
391686e3a49SStefano Zampini     ierr = MatSeqAIJRestoreArray(matis->A,&array);CHKERRQ(ierr);
392686e3a49SStefano Zampini     ierr = PetscFree(global_indices_cols);CHKERRQ(ierr);
393686e3a49SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
394c0962df8SStefano Zampini     PetscInt i,*global_indices_cols;
395c0962df8SStefano Zampini 
396c0962df8SStefano Zampini     ierr = PetscMalloc1(local_cols,&global_indices_cols);CHKERRQ(ierr);
397686e3a49SStefano Zampini     for (i=0;i<local_rows;i++) {
398686e3a49SStefano Zampini       PetscInt       j;
399686e3a49SStefano Zampini       const PetscInt *local_indices;
400686e3a49SStefano Zampini 
401686e3a49SStefano Zampini       ierr = MatGetRow(matis->A,i,&j,&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
402c0962df8SStefano Zampini       ierr = ISLocalToGlobalMappingApply(matis->mapping,j,local_indices,global_indices_cols);CHKERRQ(ierr);
403c0962df8SStefano Zampini       ierr = MatSetValues(*M,1,&global_indices_rows[i],j,global_indices_cols,array,ADD_VALUES);CHKERRQ(ierr);
404686e3a49SStefano Zampini       ierr = MatRestoreRow(matis->A,i,&j,&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
405686e3a49SStefano Zampini     }
406c0962df8SStefano Zampini     ierr = PetscFree(global_indices_cols);CHKERRQ(ierr);
407b7ce53b6SStefano Zampini   }
408b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
409b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
410b7ce53b6SStefano Zampini   if (isdense) {
411b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
412b7ce53b6SStefano Zampini   }
413b7ce53b6SStefano Zampini   if (issbaij) {
414b7ce53b6SStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
415b7ce53b6SStefano Zampini   }
416686e3a49SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&global_indices_rows);CHKERRQ(ierr);
417b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
418b7ce53b6SStefano Zampini }
419b7ce53b6SStefano Zampini 
420b7ce53b6SStefano Zampini #undef __FUNCT__
421b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
422b7ce53b6SStefano Zampini /*@
423b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
424b7ce53b6SStefano Zampini 
425b7ce53b6SStefano Zampini   Input Parameter:
426b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
427b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
428b7ce53b6SStefano Zampini 
429b7ce53b6SStefano Zampini   Output Parameter:
430b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
431b7ce53b6SStefano Zampini 
432b7ce53b6SStefano Zampini   Level: developer
433b7ce53b6SStefano Zampini 
434b7ce53b6SStefano Zampini   Notes:
435b7ce53b6SStefano Zampini 
436b7ce53b6SStefano Zampini .seealso: MATIS
437b7ce53b6SStefano Zampini @*/
438b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
439b7ce53b6SStefano Zampini {
440b7ce53b6SStefano Zampini   PetscErrorCode ierr;
441b7ce53b6SStefano Zampini 
442b7ce53b6SStefano Zampini   PetscFunctionBegin;
443b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
444b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
445b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
446b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
447b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
448b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
449b7ce53b6SStefano Zampini   }
450b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
451b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
452b7ce53b6SStefano Zampini }
453b7ce53b6SStefano Zampini 
454b7ce53b6SStefano Zampini #undef __FUNCT__
455ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
456ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
457ad6194a2SStefano Zampini {
458ad6194a2SStefano Zampini   PetscErrorCode ierr;
459ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
460ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
461ad6194a2SStefano Zampini   Mat            B,localmat;
462ad6194a2SStefano Zampini 
463ad6194a2SStefano Zampini   PetscFunctionBegin;
464ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
465ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
466ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
467ad6194a2SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr);
468ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
469ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
470b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
471ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
472ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
473ad6194a2SStefano Zampini   *newmat = B;
474ad6194a2SStefano Zampini   PetscFunctionReturn(0);
475ad6194a2SStefano Zampini }
476ad6194a2SStefano Zampini 
477ad6194a2SStefano Zampini #undef __FUNCT__
47869796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
47969796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
48069796d55SStefano Zampini {
48169796d55SStefano Zampini   PetscErrorCode ierr;
48269796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
48369796d55SStefano Zampini   PetscBool      local_sym;
48469796d55SStefano Zampini 
48569796d55SStefano Zampini   PetscFunctionBegin;
48669796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
48769796d55SStefano Zampini   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
48869796d55SStefano Zampini   PetscFunctionReturn(0);
48969796d55SStefano Zampini }
49069796d55SStefano Zampini 
49169796d55SStefano Zampini #undef __FUNCT__
49269796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
49369796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
49469796d55SStefano Zampini {
49569796d55SStefano Zampini   PetscErrorCode ierr;
49669796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
49769796d55SStefano Zampini   PetscBool      local_sym;
49869796d55SStefano Zampini 
49969796d55SStefano Zampini   PetscFunctionBegin;
50069796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
50169796d55SStefano Zampini   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
50269796d55SStefano Zampini   PetscFunctionReturn(0);
50369796d55SStefano Zampini }
50469796d55SStefano Zampini 
50569796d55SStefano Zampini #undef __FUNCT__
506b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
507dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
508b4319ba4SBarry Smith {
509dfbe8321SBarry Smith   PetscErrorCode ierr;
510b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
511b4319ba4SBarry Smith 
512b4319ba4SBarry Smith   PetscFunctionBegin;
5136bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
5146bf464f9SBarry Smith   ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr);
5156bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
5166bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
5176bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr);
51828f4e0baSStefano Zampini   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
51928f4e0baSStefano Zampini   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
520bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
521dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
522bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
523b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
524b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
5252e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
526b4319ba4SBarry Smith   PetscFunctionReturn(0);
527b4319ba4SBarry Smith }
528b4319ba4SBarry Smith 
529b4319ba4SBarry Smith #undef __FUNCT__
530b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
531dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
532b4319ba4SBarry Smith {
533dfbe8321SBarry Smith   PetscErrorCode ierr;
534b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
535b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
536b4319ba4SBarry Smith 
537b4319ba4SBarry Smith   PetscFunctionBegin;
538b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
539ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
540ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
541b4319ba4SBarry Smith 
542b4319ba4SBarry Smith   /* multiply the local matrix */
543b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
544b4319ba4SBarry Smith 
545b4319ba4SBarry Smith   /* scatter product back into global memory */
5462dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
547ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
548ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
549b4319ba4SBarry Smith   PetscFunctionReturn(0);
550b4319ba4SBarry Smith }
551b4319ba4SBarry Smith 
552b4319ba4SBarry Smith #undef __FUNCT__
5532e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
554b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5552e74eeadSLisandro Dalcin {
556650997f4SStefano Zampini   Vec            temp_vec;
5572e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5582e74eeadSLisandro Dalcin 
5592e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
560650997f4SStefano Zampini   if (v3 != v2) {
561650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
562650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
563650997f4SStefano Zampini   } else {
564650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
565650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
566650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
567650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
568650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
569650997f4SStefano Zampini   }
5702e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5712e74eeadSLisandro Dalcin }
5722e74eeadSLisandro Dalcin 
5732e74eeadSLisandro Dalcin #undef __FUNCT__
5742e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
5752e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
5762e74eeadSLisandro Dalcin {
5772e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
5782e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5792e74eeadSLisandro Dalcin 
5802e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  y = A' * x */
5812e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
582ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
583ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5842e74eeadSLisandro Dalcin 
5852e74eeadSLisandro Dalcin   /* multiply the local matrix */
5862e74eeadSLisandro Dalcin   ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr);
5872e74eeadSLisandro Dalcin 
5882e74eeadSLisandro Dalcin   /* scatter product back into global vector */
5892e74eeadSLisandro Dalcin   ierr = VecSet(y,0);CHKERRQ(ierr);
590ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
591ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5922e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5932e74eeadSLisandro Dalcin }
5942e74eeadSLisandro Dalcin 
5952e74eeadSLisandro Dalcin #undef __FUNCT__
5962e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
5972e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5982e74eeadSLisandro Dalcin {
599650997f4SStefano Zampini   Vec            temp_vec;
6002e74eeadSLisandro Dalcin   PetscErrorCode ierr;
6012e74eeadSLisandro Dalcin 
6022e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
603650997f4SStefano Zampini   if (v3 != v2) {
604650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
605650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
606650997f4SStefano Zampini   } else {
607650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
608650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
609650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
610650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
611650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
612650997f4SStefano Zampini   }
6132e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
6142e74eeadSLisandro Dalcin }
6152e74eeadSLisandro Dalcin 
6162e74eeadSLisandro Dalcin #undef __FUNCT__
617b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
618dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
619b4319ba4SBarry Smith {
620b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
621dfbe8321SBarry Smith   PetscErrorCode ierr;
622b4319ba4SBarry Smith   PetscViewer    sviewer;
623b4319ba4SBarry Smith 
624b4319ba4SBarry Smith   PetscFunctionBegin;
625dd2fa690SBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
626b4319ba4SBarry Smith   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
627b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
628b4319ba4SBarry Smith   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
629b4319ba4SBarry Smith   PetscFunctionReturn(0);
630b4319ba4SBarry Smith }
631b4319ba4SBarry Smith 
632b4319ba4SBarry Smith #undef __FUNCT__
633b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
634784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
635b4319ba4SBarry Smith {
636dfbe8321SBarry Smith   PetscErrorCode ierr;
6374e4c7dbeSStefano Zampini   PetscInt       n,bs;
638b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
639b4319ba4SBarry Smith   IS             from,to;
640b4319ba4SBarry Smith   Vec            global;
641b4319ba4SBarry Smith 
642b4319ba4SBarry Smith   PetscFunctionBegin;
643784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
644ce94432eSBarry Smith   if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical");
64570cf5478SStefano Zampini   if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */
64670cf5478SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
64770cf5478SStefano Zampini     ierr = VecDestroy(&is->x);CHKERRQ(ierr);
64870cf5478SStefano Zampini     ierr = VecDestroy(&is->y);CHKERRQ(ierr);
64970cf5478SStefano Zampini     ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr);
6501c47cb0fSStefano Zampini     ierr = MatDestroy(&is->A);CHKERRQ(ierr);
65128f4e0baSStefano Zampini     ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
65228f4e0baSStefano Zampini     ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
65370cf5478SStefano Zampini   }
654784ac674SJed Brown   ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr);
6556bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
656784ac674SJed Brown   is->mapping = rmapping;
657fa7f1dd8SStefano Zampini /*
658fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
659fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
660fa7f1dd8SStefano Zampini */
661b4319ba4SBarry Smith 
662b4319ba4SBarry Smith   /* Create the local matrix A */
663784ac674SJed Brown   ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr);
6642e1947a5SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&bs);CHKERRQ(ierr);
665f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
6662e1947a5SStefano Zampini   if (bs > 1) {
6672e1947a5SStefano Zampini     ierr = MatSetType(is->A,MATSEQBAIJ);CHKERRQ(ierr);
6682e1947a5SStefano Zampini   } else {
6692e1947a5SStefano Zampini     ierr = MatSetType(is->A,MATSEQAIJ);CHKERRQ(ierr);
6702e1947a5SStefano Zampini   }
671f69a0ea3SMatthew Knepley   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
6724e4c7dbeSStefano Zampini   ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr);
673ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
674ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
675b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
676b4319ba4SBarry Smith 
677b4319ba4SBarry Smith   /* Create the local work vectors */
6784e4c7dbeSStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr);
6794e4c7dbeSStefano Zampini   ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr);
6804e4c7dbeSStefano Zampini   ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr);
681ff130e51SJed Brown   ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr);
682ff130e51SJed Brown   ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr);
6834e4c7dbeSStefano Zampini   ierr = VecSetFromOptions(is->x);CHKERRQ(ierr);
684b4319ba4SBarry Smith   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
685b4319ba4SBarry Smith 
686b4319ba4SBarry Smith   /* setup the global to local scatter */
687b4319ba4SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
688784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
6892a7a6963SBarry Smith   ierr = MatCreateVecs(A,&global,NULL);CHKERRQ(ierr);
690b4319ba4SBarry Smith   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
6916bf464f9SBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
6926bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
6936bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
694b4319ba4SBarry Smith   PetscFunctionReturn(0);
695b4319ba4SBarry Smith }
696b4319ba4SBarry Smith 
6972e74eeadSLisandro Dalcin #undef __FUNCT__
6982e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
6992e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
7002e74eeadSLisandro Dalcin {
7012e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
7022e74eeadSLisandro Dalcin   PetscInt       rows_l[2048],cols_l[2048];
7032e74eeadSLisandro Dalcin   PetscErrorCode ierr;
7042e74eeadSLisandro Dalcin 
7052e74eeadSLisandro Dalcin   PetscFunctionBegin;
7062e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
707f23aa3ddSBarry 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);
7082e74eeadSLisandro Dalcin #endif
7092e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr);
7102e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr);
7112e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
7122e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7132e74eeadSLisandro Dalcin }
7142e74eeadSLisandro Dalcin 
7152e74eeadSLisandro Dalcin #undef ISG2LMapSetUp
7162e74eeadSLisandro Dalcin #undef ISG2LMapApply
717b4319ba4SBarry Smith 
718b4319ba4SBarry Smith #undef __FUNCT__
719b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
72013f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
721b4319ba4SBarry Smith {
722dfbe8321SBarry Smith   PetscErrorCode ierr;
723b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
724b4319ba4SBarry Smith 
725b4319ba4SBarry Smith   PetscFunctionBegin;
726b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
727b4319ba4SBarry Smith   PetscFunctionReturn(0);
728b4319ba4SBarry Smith }
729b4319ba4SBarry Smith 
730b4319ba4SBarry Smith #undef __FUNCT__
731f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
732f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
733f0006bf2SLisandro Dalcin {
734f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
735f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
736f0006bf2SLisandro Dalcin 
737f0006bf2SLisandro Dalcin   PetscFunctionBegin;
738f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
739f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
740f0006bf2SLisandro Dalcin }
741f0006bf2SLisandro Dalcin 
742f0006bf2SLisandro Dalcin #undef __FUNCT__
7432e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
7442b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
7452e74eeadSLisandro Dalcin {
7462e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
7470298fd71SBarry Smith   PetscInt       n_l = 0, *rows_l = NULL;
7482e74eeadSLisandro Dalcin   PetscErrorCode ierr;
7492e74eeadSLisandro Dalcin 
7502e74eeadSLisandro Dalcin   PetscFunctionBegin;
751ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
7522e74eeadSLisandro Dalcin   if (n) {
753785e854fSJed Brown     ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr);
7542e74eeadSLisandro Dalcin     ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
7552e74eeadSLisandro Dalcin   }
7562b40b63fSBarry Smith   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
757c31cb41cSBarry Smith   ierr = PetscFree(rows_l);CHKERRQ(ierr);
7582e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7592e74eeadSLisandro Dalcin }
7602e74eeadSLisandro Dalcin 
7612e74eeadSLisandro Dalcin #undef __FUNCT__
762b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
7632b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
764b4319ba4SBarry Smith {
765b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
766dfbe8321SBarry Smith   PetscErrorCode ierr;
767f4df32b1SMatthew Knepley   PetscInt       i;
768b4319ba4SBarry Smith   PetscScalar    *array;
769b4319ba4SBarry Smith 
770b4319ba4SBarry Smith   PetscFunctionBegin;
771ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
772b4319ba4SBarry Smith   {
773b4319ba4SBarry Smith     /*
774b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
775b4319ba4SBarry Smith        work properly in the interface nodes.
776b4319ba4SBarry Smith     */
777b4319ba4SBarry Smith     Vec         counter;
778b4319ba4SBarry Smith     PetscScalar one=1.0, zero=0.0;
7792a7a6963SBarry Smith     ierr = MatCreateVecs(A,&counter,NULL);CHKERRQ(ierr);
7802dcb1b2aSMatthew Knepley     ierr = VecSet(counter,zero);CHKERRQ(ierr);
7812dcb1b2aSMatthew Knepley     ierr = VecSet(is->x,one);CHKERRQ(ierr);
782ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
783ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
784ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
785ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7866bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
787b4319ba4SBarry Smith   }
788958c9bccSBarry Smith   if (!n) {
789b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
790b4319ba4SBarry Smith   } else {
791b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
7922205254eSKarl Rupp 
793b4319ba4SBarry Smith     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
7942b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
795b4319ba4SBarry Smith     for (i=0; i<n; i++) {
796f4df32b1SMatthew Knepley       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
797b4319ba4SBarry Smith     }
798b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
799b4319ba4SBarry Smith     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
800b4319ba4SBarry Smith     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
801b4319ba4SBarry Smith   }
802b4319ba4SBarry Smith   PetscFunctionReturn(0);
803b4319ba4SBarry Smith }
804b4319ba4SBarry Smith 
805b4319ba4SBarry Smith #undef __FUNCT__
806b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
807dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
808b4319ba4SBarry Smith {
809b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
810dfbe8321SBarry Smith   PetscErrorCode ierr;
811b4319ba4SBarry Smith 
812b4319ba4SBarry Smith   PetscFunctionBegin;
813b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
814b4319ba4SBarry Smith   PetscFunctionReturn(0);
815b4319ba4SBarry Smith }
816b4319ba4SBarry Smith 
817b4319ba4SBarry Smith #undef __FUNCT__
818b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
819dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
820b4319ba4SBarry Smith {
821b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
822dfbe8321SBarry Smith   PetscErrorCode ierr;
823b4319ba4SBarry Smith 
824b4319ba4SBarry Smith   PetscFunctionBegin;
825b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
826b4319ba4SBarry Smith   PetscFunctionReturn(0);
827b4319ba4SBarry Smith }
828b4319ba4SBarry Smith 
829b4319ba4SBarry Smith #undef __FUNCT__
830b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
8317087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
832b4319ba4SBarry Smith {
833b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
834b4319ba4SBarry Smith 
835b4319ba4SBarry Smith   PetscFunctionBegin;
836b4319ba4SBarry Smith   *local = is->A;
837b4319ba4SBarry Smith   PetscFunctionReturn(0);
838b4319ba4SBarry Smith }
839b4319ba4SBarry Smith 
840b4319ba4SBarry Smith #undef __FUNCT__
841b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
842b4319ba4SBarry Smith /*@
843b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
844b4319ba4SBarry Smith 
845b4319ba4SBarry Smith   Input Parameter:
846b4319ba4SBarry Smith .  mat - the matrix
847b4319ba4SBarry Smith 
848b4319ba4SBarry Smith   Output Parameter:
849b4319ba4SBarry Smith .  local - the local matrix usually MATSEQAIJ
850b4319ba4SBarry Smith 
851b4319ba4SBarry Smith   Level: advanced
852b4319ba4SBarry Smith 
853b4319ba4SBarry Smith   Notes:
854b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
855b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
856b4319ba4SBarry Smith   of the MatSetValues() operation.
857b4319ba4SBarry Smith 
858b4319ba4SBarry Smith .seealso: MATIS
859b4319ba4SBarry Smith @*/
8607087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
861b4319ba4SBarry Smith {
8624ac538c5SBarry Smith   PetscErrorCode ierr;
863b4319ba4SBarry Smith 
864b4319ba4SBarry Smith   PetscFunctionBegin;
8650700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
866b4319ba4SBarry Smith   PetscValidPointer(local,2);
8674ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
868b4319ba4SBarry Smith   PetscFunctionReturn(0);
869b4319ba4SBarry Smith }
870b4319ba4SBarry Smith 
8713b03a366Sstefano_zampini #undef __FUNCT__
8723b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
8733b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
8743b03a366Sstefano_zampini {
8753b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
8763b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
8773b03a366Sstefano_zampini   PetscErrorCode ierr;
8783b03a366Sstefano_zampini 
8793b03a366Sstefano_zampini   PetscFunctionBegin;
8804e4c7dbeSStefano Zampini   if (is->A) {
8813b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
8823b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
883f23aa3ddSBarry 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);
8844e4c7dbeSStefano Zampini   }
8853b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
8863b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
8873b03a366Sstefano_zampini   is->A = local;
8883b03a366Sstefano_zampini   PetscFunctionReturn(0);
8893b03a366Sstefano_zampini }
8903b03a366Sstefano_zampini 
8913b03a366Sstefano_zampini #undef __FUNCT__
8923b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
8933b03a366Sstefano_zampini /*@
8943b03a366Sstefano_zampini     MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix.
8953b03a366Sstefano_zampini 
8963b03a366Sstefano_zampini   Input Parameter:
8973b03a366Sstefano_zampini .  mat - the matrix
8983b03a366Sstefano_zampini .  local - the local matrix usually MATSEQAIJ
8993b03a366Sstefano_zampini 
9003b03a366Sstefano_zampini   Output Parameter:
9013b03a366Sstefano_zampini 
9023b03a366Sstefano_zampini   Level: advanced
9033b03a366Sstefano_zampini 
9043b03a366Sstefano_zampini   Notes:
9053b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
9063b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
9073b03a366Sstefano_zampini 
9083b03a366Sstefano_zampini .seealso: MATIS
9093b03a366Sstefano_zampini @*/
9103b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
9113b03a366Sstefano_zampini {
9123b03a366Sstefano_zampini   PetscErrorCode ierr;
9133b03a366Sstefano_zampini 
9143b03a366Sstefano_zampini   PetscFunctionBegin;
9153b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
916b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
9173b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
9183b03a366Sstefano_zampini   PetscFunctionReturn(0);
9193b03a366Sstefano_zampini }
9203b03a366Sstefano_zampini 
9216726f965SBarry Smith #undef __FUNCT__
9226726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
9236726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
9246726f965SBarry Smith {
9256726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
9266726f965SBarry Smith   PetscErrorCode ierr;
9276726f965SBarry Smith 
9286726f965SBarry Smith   PetscFunctionBegin;
9296726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
9306726f965SBarry Smith   PetscFunctionReturn(0);
9316726f965SBarry Smith }
9326726f965SBarry Smith 
9336726f965SBarry Smith #undef __FUNCT__
9342e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
9352e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
9362e74eeadSLisandro Dalcin {
9372e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9382e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9392e74eeadSLisandro Dalcin 
9402e74eeadSLisandro Dalcin   PetscFunctionBegin;
9412e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
9422e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9432e74eeadSLisandro Dalcin }
9442e74eeadSLisandro Dalcin 
9452e74eeadSLisandro Dalcin #undef __FUNCT__
9462e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
9472e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
9482e74eeadSLisandro Dalcin {
9492e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9502e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9512e74eeadSLisandro Dalcin 
9522e74eeadSLisandro Dalcin   PetscFunctionBegin;
9532e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
9542e74eeadSLisandro Dalcin   ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr);
9552e74eeadSLisandro Dalcin 
9562e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
9572e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
958ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
959ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9602e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9612e74eeadSLisandro Dalcin }
9622e74eeadSLisandro Dalcin 
9632e74eeadSLisandro Dalcin #undef __FUNCT__
9646726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
965ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
9666726f965SBarry Smith {
9676726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
9686726f965SBarry Smith   PetscErrorCode ierr;
9696726f965SBarry Smith 
9706726f965SBarry Smith   PetscFunctionBegin;
9714e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
9726726f965SBarry Smith   PetscFunctionReturn(0);
9736726f965SBarry Smith }
9746726f965SBarry Smith 
975284134d9SBarry Smith #undef __FUNCT__
976284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
977284134d9SBarry Smith /*@
978284134d9SBarry Smith     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
979284134d9SBarry Smith        process but not across processes.
980284134d9SBarry Smith 
981284134d9SBarry Smith    Input Parameters:
982284134d9SBarry Smith +     comm - MPI communicator that will share the matrix
9834e4c7dbeSStefano Zampini .     bs - local and global block size of the matrix
984df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
985284134d9SBarry Smith -     map - mapping that defines the global number for each local number
986284134d9SBarry Smith 
987284134d9SBarry Smith    Output Parameter:
988284134d9SBarry Smith .    A - the resulting matrix
989284134d9SBarry Smith 
9908e6c10adSSatish Balay    Level: advanced
9918e6c10adSSatish Balay 
992284134d9SBarry Smith    Notes: See MATIS for more details
9938cda6cd7SBarry Smith           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
9948cda6cd7SBarry Smith           by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
9958cda6cd7SBarry Smith           plus the ghost points to global indices.
996284134d9SBarry Smith 
997284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
998284134d9SBarry Smith @*/
9994e4c7dbeSStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
1000284134d9SBarry Smith {
1001284134d9SBarry Smith   PetscErrorCode ierr;
1002284134d9SBarry Smith 
1003284134d9SBarry Smith   PetscFunctionBegin;
1004284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1005d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
1006284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1007284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
10083b03a366Sstefano_zampini   ierr = MatSetUp(*A);CHKERRQ(ierr);
1009784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr);
1010284134d9SBarry Smith   PetscFunctionReturn(0);
1011284134d9SBarry Smith }
1012284134d9SBarry Smith 
1013b4319ba4SBarry Smith /*MC
10146a818285SBarry Smith    MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners, see PCBDDC.
1015b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
1016b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
1017b4319ba4SBarry Smith    product is handled "implicitly".
1018b4319ba4SBarry Smith 
1019b4319ba4SBarry Smith    Operations Provided:
10206726f965SBarry Smith +  MatMult()
10212e74eeadSLisandro Dalcin .  MatMultAdd()
10222e74eeadSLisandro Dalcin .  MatMultTranspose()
10232e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
10246726f965SBarry Smith .  MatZeroEntries()
10256726f965SBarry Smith .  MatSetOption()
10262e74eeadSLisandro Dalcin .  MatZeroRows()
10276726f965SBarry Smith .  MatZeroRowsLocal()
10282e74eeadSLisandro Dalcin .  MatSetValues()
10296726f965SBarry Smith .  MatSetValuesLocal()
10302e74eeadSLisandro Dalcin .  MatScale()
10312e74eeadSLisandro Dalcin .  MatGetDiagonal()
10326726f965SBarry Smith -  MatSetLocalToGlobalMapping()
1033b4319ba4SBarry Smith 
1034b4319ba4SBarry Smith    Options Database Keys:
1035b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
1036b4319ba4SBarry Smith 
1037b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
1038b4319ba4SBarry Smith 
1039b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
1040b4319ba4SBarry Smith 
1041b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
1042b4319ba4SBarry Smith           MatISGetLocalMat()
1043b4319ba4SBarry Smith 
1044b4319ba4SBarry Smith   Level: advanced
1045b4319ba4SBarry Smith 
10466a818285SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), PCBDDC
1047b4319ba4SBarry Smith 
1048b4319ba4SBarry Smith M*/
1049b4319ba4SBarry Smith 
1050b4319ba4SBarry Smith #undef __FUNCT__
1051b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
10528cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
1053b4319ba4SBarry Smith {
1054dfbe8321SBarry Smith   PetscErrorCode ierr;
1055b4319ba4SBarry Smith   Mat_IS         *b;
1056b4319ba4SBarry Smith 
1057b4319ba4SBarry Smith   PetscFunctionBegin;
1058b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1059b4319ba4SBarry Smith   A->data = (void*)b;
1060b4319ba4SBarry Smith   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1061b4319ba4SBarry Smith 
1062b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
10632e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
10642e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
10652e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1066b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1067b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
10682e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
1069b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1070f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
10712e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1072b4319ba4SBarry Smith   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
1073b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1074b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1075b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
10766726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
10772e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
10782e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
10796726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
108069796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
108169796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1082ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
1083b4319ba4SBarry Smith 
108426283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
108526283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1086b4319ba4SBarry Smith 
1087b7ce53b6SStefano Zampini   /* zeros pointer */
1088b4319ba4SBarry Smith   b->A           = 0;
1089b4319ba4SBarry Smith   b->ctx         = 0;
1090b4319ba4SBarry Smith   b->x           = 0;
1091b4319ba4SBarry Smith   b->y           = 0;
1092b09366fdSStefano Zampini   b->mapping     = 0;
109328f4e0baSStefano Zampini   b->sf          = 0;
109428f4e0baSStefano Zampini   b->sf_rootdata = 0;
109528f4e0baSStefano Zampini   b->sf_leafdata = 0;
1096b7ce53b6SStefano Zampini 
1097b7ce53b6SStefano Zampini   /* special MATIS functions */
1098bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1099bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1100b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
11012e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
110217667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1103b4319ba4SBarry Smith   PetscFunctionReturn(0);
1104b4319ba4SBarry Smith }
1105