xref: /petsc/src/mat/impls/is/matis.c (revision 97563a804b114d3565bf1d112fdcd9e3b1c38d9b)
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     Currently this allows for only one subdomain per processor.
9b4319ba4SBarry Smith */
10b4319ba4SBarry Smith 
11c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h>      /*I "petscmat.h" I*/
1228f4e0baSStefano Zampini 
13a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat);
147230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat,PetscInt,const PetscInt[],PetscScalar);
15a72627d2SStefano Zampini 
1628f4e0baSStefano Zampini #undef __FUNCT__
1728f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private"
18a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B)
1928f4e0baSStefano Zampini {
2028f4e0baSStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
2128f4e0baSStefano Zampini   const PetscInt *gidxs;
2228f4e0baSStefano Zampini   PetscErrorCode ierr;
2328f4e0baSStefano Zampini 
2428f4e0baSStefano Zampini   PetscFunctionBegin;
2528f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr);
2628f4e0baSStefano Zampini   ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr);
2728f4e0baSStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr);
283bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
2928f4e0baSStefano Zampini   /* PETSC_OWN_POINTER refers to ilocal which is NULL */
3028f4e0baSStefano Zampini   ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
313bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
3228f4e0baSStefano Zampini   ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr);
3328f4e0baSStefano Zampini   PetscFunctionReturn(0);
3428f4e0baSStefano Zampini }
352e1947a5SStefano Zampini 
362e1947a5SStefano Zampini #undef __FUNCT__
372e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation"
38eb82efa4SStefano Zampini /*@
39a88811baSStefano Zampini    MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix.
40a88811baSStefano Zampini 
41a88811baSStefano Zampini    Collective on MPI_Comm
42a88811baSStefano Zampini 
43a88811baSStefano Zampini    Input Parameters:
44a88811baSStefano Zampini +  B - the matrix
45a88811baSStefano Zampini .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
46a88811baSStefano Zampini            (same value is used for all local rows)
47a88811baSStefano Zampini .  d_nnz - array containing the number of nonzeros in the various rows of the
48a88811baSStefano Zampini            DIAGONAL portion of the local submatrix (possibly different for each row)
49a88811baSStefano Zampini            or NULL, if d_nz is used to specify the nonzero structure.
50a88811baSStefano Zampini            The size of this array is equal to the number of local rows, i.e 'm'.
51a88811baSStefano Zampini            For matrices that will be factored, you must leave room for (and set)
52a88811baSStefano Zampini            the diagonal entry even if it is zero.
53a88811baSStefano Zampini .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
54a88811baSStefano Zampini            submatrix (same value is used for all local rows).
55a88811baSStefano Zampini -  o_nnz - array containing the number of nonzeros in the various rows of the
56a88811baSStefano Zampini            OFF-DIAGONAL portion of the local submatrix (possibly different for
57a88811baSStefano Zampini            each row) or NULL, if o_nz is used to specify the nonzero
58a88811baSStefano Zampini            structure. The size of this array is equal to the number
59a88811baSStefano Zampini            of local rows, i.e 'm'.
60a88811baSStefano Zampini 
61a88811baSStefano Zampini    If the *_nnz parameter is given then the *_nz parameter is ignored
62a88811baSStefano Zampini 
63a88811baSStefano Zampini    Level: intermediate
64a88811baSStefano Zampini 
65a88811baSStefano Zampini    Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition
66a88811baSStefano Zampini           from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local
67a88811baSStefano Zampini           matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects.
68a88811baSStefano Zampini 
69a88811baSStefano Zampini .keywords: matrix
70a88811baSStefano Zampini 
713c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS
72a88811baSStefano Zampini @*/
732e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
742e1947a5SStefano Zampini {
752e1947a5SStefano Zampini   PetscErrorCode ierr;
762e1947a5SStefano Zampini 
772e1947a5SStefano Zampini   PetscFunctionBegin;
782e1947a5SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
792e1947a5SStefano Zampini   PetscValidType(B,1);
802e1947a5SStefano Zampini   ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
812e1947a5SStefano Zampini   PetscFunctionReturn(0);
822e1947a5SStefano Zampini }
832e1947a5SStefano Zampini 
842e1947a5SStefano Zampini #undef __FUNCT__
852e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS"
867230de76SStefano Zampini static PetscErrorCode  MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
872e1947a5SStefano Zampini {
882e1947a5SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
8928f4e0baSStefano Zampini   PetscInt       bs,i,nlocalcols;
902e1947a5SStefano Zampini   PetscErrorCode ierr;
912e1947a5SStefano Zampini 
922e1947a5SStefano Zampini   PetscFunctionBegin;
936c4ed002SBarry Smith   if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
9428f4e0baSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
9528f4e0baSStefano Zampini     ierr = MatISComputeSF_Private(B);CHKERRQ(ierr);
9628f4e0baSStefano Zampini   }
972e1947a5SStefano Zampini   if (!d_nnz) {
9828f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz;
992e1947a5SStefano Zampini   } else {
10028f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i];
1012e1947a5SStefano Zampini   }
1022e1947a5SStefano Zampini   if (!o_nnz) {
10328f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz;
1042e1947a5SStefano Zampini   } else {
10528f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i];
1062e1947a5SStefano Zampini   }
10728f4e0baSStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
10828f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr);
10928f4e0baSStefano Zampini   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
11028f4e0baSStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
11128f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves;i++) {
11228f4e0baSStefano Zampini     matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols);
1132e1947a5SStefano Zampini   }
11428f4e0baSStefano Zampini   ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr);
11528f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
11628f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs;
1172e1947a5SStefano Zampini   }
11828f4e0baSStefano Zampini   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
11928f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
12028f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i]-i;
1212e1947a5SStefano Zampini   }
12228f4e0baSStefano Zampini   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
1232e1947a5SStefano Zampini   PetscFunctionReturn(0);
1242e1947a5SStefano Zampini }
125b4319ba4SBarry Smith 
126b4319ba4SBarry Smith #undef __FUNCT__
1273927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private"
1283927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
1293927de2eSStefano Zampini {
1303927de2eSStefano Zampini   Mat_IS          *matis = (Mat_IS*)(A->data);
1313927de2eSStefano Zampini   PetscInt        *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership;
132ecf5a873SStefano Zampini   const PetscInt  *global_indices_r,*global_indices_c;
1333927de2eSStefano Zampini   PetscInt        i,j,bs,rows,cols;
1343927de2eSStefano Zampini   PetscInt        lrows,lcols;
1353927de2eSStefano Zampini   PetscInt        local_rows,local_cols;
1363927de2eSStefano Zampini   PetscMPIInt     nsubdomains;
1373927de2eSStefano Zampini   PetscBool       isdense,issbaij;
1383927de2eSStefano Zampini   PetscErrorCode  ierr;
1393927de2eSStefano Zampini 
1403927de2eSStefano Zampini   PetscFunctionBegin;
1413927de2eSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr);
1423927de2eSStefano Zampini   ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr);
1433927de2eSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1443927de2eSStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
1453927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
1463927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
147ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
148ecf5a873SStefano Zampini   if (A->rmap->mapping != A->cmap->mapping) {
1497230de76SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
150ecf5a873SStefano Zampini   } else {
151ecf5a873SStefano Zampini     global_indices_c = global_indices_r;
152ecf5a873SStefano Zampini   }
153ecf5a873SStefano Zampini 
1543927de2eSStefano Zampini   if (issbaij) {
1553927de2eSStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
1563927de2eSStefano Zampini   }
1573927de2eSStefano Zampini   /*
158ecf5a873SStefano Zampini      An SF reduce is needed to sum up properly on shared rows.
1593927de2eSStefano Zampini      Note that generally preallocation is not exact, since it overestimates nonzeros
1603927de2eSStefano Zampini   */
1613927de2eSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
1623927de2eSStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
1633927de2eSStefano Zampini   }
1643927de2eSStefano Zampini   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
1653927de2eSStefano Zampini   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
1663927de2eSStefano Zampini   /* All processes need to compute entire row ownership */
1673927de2eSStefano Zampini   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
1683927de2eSStefano Zampini   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
1693927de2eSStefano Zampini   for (i=0;i<nsubdomains;i++) {
1703927de2eSStefano Zampini     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
1713927de2eSStefano Zampini       row_ownership[j] = i;
1723927de2eSStefano Zampini     }
1733927de2eSStefano Zampini   }
1747230de76SStefano Zampini   ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
1753927de2eSStefano Zampini 
1763927de2eSStefano Zampini   /*
1773927de2eSStefano Zampini      my_dnz and my_onz contains exact contribution to preallocation from each local mat
1783927de2eSStefano Zampini      then, they will be summed up properly. This way, preallocation is always sufficient
1793927de2eSStefano Zampini   */
1803927de2eSStefano Zampini   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
1813927de2eSStefano Zampini   /* preallocation as a MATAIJ */
1823927de2eSStefano Zampini   if (isdense) { /* special case for dense local matrices */
1833927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
184ecf5a873SStefano Zampini       PetscInt index_row = global_indices_r[i];
1853927de2eSStefano Zampini       for (j=i;j<local_rows;j++) {
1863927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
187ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[j];
1883927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
1893927de2eSStefano Zampini           my_dnz[i] += 1;
1903927de2eSStefano Zampini         } else { /* offdiag block */
1913927de2eSStefano Zampini           my_onz[i] += 1;
1923927de2eSStefano Zampini         }
1933927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
1943927de2eSStefano Zampini         if (i != j) {
1953927de2eSStefano Zampini           owner = row_ownership[index_col];
1963927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
1973927de2eSStefano Zampini             my_dnz[j] += 1;
1983927de2eSStefano Zampini           } else {
1993927de2eSStefano Zampini             my_onz[j] += 1;
2003927de2eSStefano Zampini           }
2013927de2eSStefano Zampini         }
2023927de2eSStefano Zampini       }
2033927de2eSStefano Zampini     }
204ecf5a873SStefano Zampini   } else { /* TODO: this could be optimized using MatGetRowIJ */
2053927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
2063927de2eSStefano Zampini       const PetscInt *cols;
207ecf5a873SStefano Zampini       PetscInt       ncols,index_row = global_indices_r[i];
2083927de2eSStefano Zampini       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2093927de2eSStefano Zampini       for (j=0;j<ncols;j++) {
2103927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
211ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[cols[j]];
2123927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
2133927de2eSStefano Zampini           my_dnz[i] += 1;
2143927de2eSStefano Zampini         } else { /* offdiag block */
2153927de2eSStefano Zampini           my_onz[i] += 1;
2163927de2eSStefano Zampini         }
2173927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
218d9a9e74cSStefano Zampini         if (issbaij && index_col != index_row) {
2193927de2eSStefano Zampini           owner = row_ownership[index_col];
2203927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
221d9a9e74cSStefano Zampini             my_dnz[cols[j]] += 1;
2223927de2eSStefano Zampini           } else {
223d9a9e74cSStefano Zampini             my_onz[cols[j]] += 1;
2243927de2eSStefano Zampini           }
2253927de2eSStefano Zampini         }
2263927de2eSStefano Zampini       }
2273927de2eSStefano Zampini       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2283927de2eSStefano Zampini     }
2293927de2eSStefano Zampini   }
230ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
231ecf5a873SStefano Zampini   if (global_indices_c != global_indices_r) {
2327230de76SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
233ecf5a873SStefano Zampini   }
2343927de2eSStefano Zampini   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
235ecf5a873SStefano Zampini 
236ecf5a873SStefano Zampini   /* Reduce my_dnz and my_onz */
2373927de2eSStefano Zampini   if (maxreduce) {
2383927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2393927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2403927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2413927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2423927de2eSStefano Zampini   } else {
2433927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2443927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2453927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2463927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2473927de2eSStefano Zampini   }
2483927de2eSStefano Zampini   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
2493927de2eSStefano Zampini 
2503927de2eSStefano Zampini   /* Resize preallocation if overestimated */
2513927de2eSStefano Zampini   for (i=0;i<lrows;i++) {
2523927de2eSStefano Zampini     dnz[i] = PetscMin(dnz[i],lcols);
2533927de2eSStefano Zampini     onz[i] = PetscMin(onz[i],cols-lcols);
2543927de2eSStefano Zampini   }
2553927de2eSStefano Zampini   /* set preallocation */
2563927de2eSStefano Zampini   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
2573927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
2583927de2eSStefano Zampini     dnz[i] = dnz[i*bs]/bs;
2593927de2eSStefano Zampini     onz[i] = onz[i*bs]/bs;
2603927de2eSStefano Zampini   }
2613927de2eSStefano Zampini   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2623927de2eSStefano Zampini   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2633927de2eSStefano Zampini   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
2643927de2eSStefano Zampini   if (issbaij) {
2653927de2eSStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
2663927de2eSStefano Zampini   }
2673927de2eSStefano Zampini   PetscFunctionReturn(0);
2683927de2eSStefano Zampini }
2693927de2eSStefano Zampini 
2703927de2eSStefano Zampini #undef __FUNCT__
271b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
2727230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
273b7ce53b6SStefano Zampini {
274b7ce53b6SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
275d9a9e74cSStefano Zampini   Mat            local_mat;
276b7ce53b6SStefano Zampini   /* info on mat */
2773cfa4ea4SStefano Zampini   PetscInt       bs,rows,cols,lrows,lcols;
278b7ce53b6SStefano Zampini   PetscInt       local_rows,local_cols;
279686e3a49SStefano Zampini   PetscBool      isdense,issbaij,isseqaij;
2807c03b4e8SStefano Zampini   PetscMPIInt    nsubdomains;
281b7ce53b6SStefano Zampini   /* values insertion */
282b7ce53b6SStefano Zampini   PetscScalar    *array;
283b7ce53b6SStefano Zampini   /* work */
284b7ce53b6SStefano Zampini   PetscErrorCode ierr;
285b7ce53b6SStefano Zampini 
286b7ce53b6SStefano Zampini   PetscFunctionBegin;
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 = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
298b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
2993cfa4ea4SStefano Zampini   ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr);
300b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
301b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
302686e3a49SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
303b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
304b7ce53b6SStefano Zampini 
305b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
306b7ce53b6SStefano Zampini     MatType     new_mat_type;
3073927de2eSStefano Zampini     PetscBool   issbaij_red;
308b7ce53b6SStefano Zampini 
309b7ce53b6SStefano Zampini     /* determining new matrix type */
310b2566f29SBarry Smith     ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
311b7ce53b6SStefano Zampini     if (issbaij_red) {
312b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
313b7ce53b6SStefano Zampini     } else {
314b7ce53b6SStefano Zampini       if (bs>1) {
315b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
316b7ce53b6SStefano Zampini       } else {
317b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
318b7ce53b6SStefano Zampini       }
319b7ce53b6SStefano Zampini     }
320b7ce53b6SStefano Zampini 
3213927de2eSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
3223cfa4ea4SStefano Zampini     ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr);
3233927de2eSStefano Zampini     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
3243927de2eSStefano Zampini     ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr);
3253927de2eSStefano Zampini     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
326b7ce53b6SStefano Zampini   } else {
3273cfa4ea4SStefano Zampini     PetscInt mbs,mrows,mcols,mlrows,mlcols;
328b7ce53b6SStefano Zampini     /* some checks */
329b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
330b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
3313cfa4ea4SStefano Zampini     ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr);
3326c4ed002SBarry Smith     if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
3336c4ed002SBarry Smith     if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
3346c4ed002SBarry Smith     if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows);
3356c4ed002SBarry Smith     if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols);
3366c4ed002SBarry Smith     if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
337b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
338b7ce53b6SStefano Zampini   }
339d9a9e74cSStefano Zampini 
340d9a9e74cSStefano Zampini   if (issbaij) {
341d9a9e74cSStefano Zampini     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr);
342d9a9e74cSStefano Zampini   } else {
343d9a9e74cSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
344d9a9e74cSStefano Zampini     local_mat = matis->A;
345d9a9e74cSStefano Zampini   }
346686e3a49SStefano Zampini 
347b7ce53b6SStefano Zampini   /* Set values */
348ecf5a873SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr);
349b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
350ecf5a873SStefano Zampini     PetscInt i,*dummy_rows,*dummy_cols;
351ecf5a873SStefano Zampini 
352ecf5a873SStefano Zampini     if (local_rows != local_cols) {
353ecf5a873SStefano Zampini       ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr);
354ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
355ecf5a873SStefano Zampini       for (i=0;i<local_cols;i++) dummy_cols[i] = i;
356ecf5a873SStefano Zampini     } else {
357ecf5a873SStefano Zampini       ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr);
358ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
359ecf5a873SStefano Zampini       dummy_cols = dummy_rows;
360ecf5a873SStefano Zampini     }
361b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
362d9a9e74cSStefano Zampini     ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr);
363ecf5a873SStefano Zampini     ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr);
364d9a9e74cSStefano Zampini     ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr);
365ecf5a873SStefano Zampini     if (dummy_rows != dummy_cols) {
366ecf5a873SStefano Zampini       ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr);
367ecf5a873SStefano Zampini     } else {
368ecf5a873SStefano Zampini       ierr = PetscFree(dummy_rows);CHKERRQ(ierr);
369ecf5a873SStefano Zampini     }
370686e3a49SStefano Zampini   } else if (isseqaij) {
371ecf5a873SStefano Zampini     PetscInt  i,nvtxs,*xadj,*adjncy;
372686e3a49SStefano Zampini     PetscBool done;
373686e3a49SStefano Zampini 
374d9a9e74cSStefano Zampini     ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
3756c4ed002SBarry Smith     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
376d9a9e74cSStefano Zampini     ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr);
377686e3a49SStefano Zampini     for (i=0;i<nvtxs;i++) {
378ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr);
379686e3a49SStefano Zampini     }
380d9a9e74cSStefano Zampini     ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
3816c4ed002SBarry Smith     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
382d9a9e74cSStefano Zampini     ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr);
383686e3a49SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
384ecf5a873SStefano Zampini     PetscInt i;
385c0962df8SStefano Zampini 
386686e3a49SStefano Zampini     for (i=0;i<local_rows;i++) {
387686e3a49SStefano Zampini       PetscInt       j;
388ecf5a873SStefano Zampini       const PetscInt *local_indices_cols;
389686e3a49SStefano Zampini 
390ecf5a873SStefano Zampini       ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
391ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr);
392ecf5a873SStefano Zampini       ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
393686e3a49SStefano Zampini     }
394b7ce53b6SStefano Zampini   }
395b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
396d9a9e74cSStefano Zampini   ierr = MatDestroy(&local_mat);CHKERRQ(ierr);
397b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
398b7ce53b6SStefano Zampini   if (isdense) {
399b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
400b7ce53b6SStefano Zampini   }
401b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
402b7ce53b6SStefano Zampini }
403b7ce53b6SStefano Zampini 
404b7ce53b6SStefano Zampini #undef __FUNCT__
405b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
406b7ce53b6SStefano Zampini /*@
407b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
408b7ce53b6SStefano Zampini 
409b7ce53b6SStefano Zampini   Input Parameter:
410b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
411b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
412b7ce53b6SStefano Zampini 
413b7ce53b6SStefano Zampini   Output Parameter:
414b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
415b7ce53b6SStefano Zampini 
416b7ce53b6SStefano Zampini   Level: developer
417b7ce53b6SStefano Zampini 
418eb82efa4SStefano Zampini   Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested.
419b7ce53b6SStefano Zampini 
420b7ce53b6SStefano Zampini .seealso: MATIS
421b7ce53b6SStefano Zampini @*/
422b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
423b7ce53b6SStefano Zampini {
424b7ce53b6SStefano Zampini   PetscErrorCode ierr;
425b7ce53b6SStefano Zampini 
426b7ce53b6SStefano Zampini   PetscFunctionBegin;
427b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
428b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
429b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
430b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
431b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
432b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
4336c4ed002SBarry Smith     if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix");
434b7ce53b6SStefano Zampini   }
435b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
436b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
437b7ce53b6SStefano Zampini }
438b7ce53b6SStefano Zampini 
439b7ce53b6SStefano Zampini #undef __FUNCT__
440ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
441ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
442ad6194a2SStefano Zampini {
443ad6194a2SStefano Zampini   PetscErrorCode ierr;
444ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
445ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
446ad6194a2SStefano Zampini   Mat            B,localmat;
447ad6194a2SStefano Zampini 
448ad6194a2SStefano Zampini   PetscFunctionBegin;
449ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
450ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
451ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
452e176bc59SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr);
453ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
454ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
455b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
456ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
457ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
458ad6194a2SStefano Zampini   *newmat = B;
459ad6194a2SStefano Zampini   PetscFunctionReturn(0);
460ad6194a2SStefano Zampini }
461ad6194a2SStefano Zampini 
462ad6194a2SStefano Zampini #undef __FUNCT__
46369796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
46469796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
46569796d55SStefano Zampini {
46669796d55SStefano Zampini   PetscErrorCode ierr;
46769796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
46869796d55SStefano Zampini   PetscBool      local_sym;
46969796d55SStefano Zampini 
47069796d55SStefano Zampini   PetscFunctionBegin;
47169796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
472b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
47369796d55SStefano Zampini   PetscFunctionReturn(0);
47469796d55SStefano Zampini }
47569796d55SStefano Zampini 
47669796d55SStefano Zampini #undef __FUNCT__
47769796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
47869796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
47969796d55SStefano Zampini {
48069796d55SStefano Zampini   PetscErrorCode ierr;
48169796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
48269796d55SStefano Zampini   PetscBool      local_sym;
48369796d55SStefano Zampini 
48469796d55SStefano Zampini   PetscFunctionBegin;
48569796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
486b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
48769796d55SStefano Zampini   PetscFunctionReturn(0);
48869796d55SStefano Zampini }
48969796d55SStefano Zampini 
49069796d55SStefano Zampini #undef __FUNCT__
491b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
492dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
493b4319ba4SBarry Smith {
494dfbe8321SBarry Smith   PetscErrorCode ierr;
495b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
496b4319ba4SBarry Smith 
497b4319ba4SBarry Smith   PetscFunctionBegin;
4986bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
499e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr);
500e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr);
5016bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
5026bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
50328f4e0baSStefano Zampini   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
50428f4e0baSStefano Zampini   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
505bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
506dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
507bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
508b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
509b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
5102e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
511b4319ba4SBarry Smith   PetscFunctionReturn(0);
512b4319ba4SBarry Smith }
513b4319ba4SBarry Smith 
514b4319ba4SBarry Smith #undef __FUNCT__
515b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
516dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
517b4319ba4SBarry Smith {
518dfbe8321SBarry Smith   PetscErrorCode ierr;
519b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
520b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
521b4319ba4SBarry Smith 
522b4319ba4SBarry Smith   PetscFunctionBegin;
523b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
524e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
525e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
526b4319ba4SBarry Smith 
527b4319ba4SBarry Smith   /* multiply the local matrix */
528b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
529b4319ba4SBarry Smith 
530b4319ba4SBarry Smith   /* scatter product back into global memory */
5312dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
532e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
533e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
534b4319ba4SBarry Smith   PetscFunctionReturn(0);
535b4319ba4SBarry Smith }
536b4319ba4SBarry Smith 
537b4319ba4SBarry Smith #undef __FUNCT__
5382e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
539b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5402e74eeadSLisandro Dalcin {
541650997f4SStefano Zampini   Vec            temp_vec;
5422e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5432e74eeadSLisandro Dalcin 
5442e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
545650997f4SStefano Zampini   if (v3 != v2) {
546650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
547650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
548650997f4SStefano Zampini   } else {
549650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
550650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
551650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
552650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
553650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
554650997f4SStefano Zampini   }
5552e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5562e74eeadSLisandro Dalcin }
5572e74eeadSLisandro Dalcin 
5582e74eeadSLisandro Dalcin #undef __FUNCT__
5592e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
560e176bc59SStefano Zampini PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x)
5612e74eeadSLisandro Dalcin {
5622e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
5632e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5642e74eeadSLisandro Dalcin 
565e176bc59SStefano Zampini   PetscFunctionBegin;
5662e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
567e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
568e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5692e74eeadSLisandro Dalcin 
5702e74eeadSLisandro Dalcin   /* multiply the local matrix */
571e176bc59SStefano Zampini   ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr);
5722e74eeadSLisandro Dalcin 
5732e74eeadSLisandro Dalcin   /* scatter product back into global vector */
574e176bc59SStefano Zampini   ierr = VecSet(x,0);CHKERRQ(ierr);
575e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
576e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5772e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5782e74eeadSLisandro Dalcin }
5792e74eeadSLisandro Dalcin 
5802e74eeadSLisandro Dalcin #undef __FUNCT__
5812e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
5822e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5832e74eeadSLisandro Dalcin {
584650997f4SStefano Zampini   Vec            temp_vec;
5852e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5862e74eeadSLisandro Dalcin 
5872e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
588650997f4SStefano Zampini   if (v3 != v2) {
589650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
590650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
591650997f4SStefano Zampini   } else {
592650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
593650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
594650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
595650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
596650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
597650997f4SStefano Zampini   }
5982e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5992e74eeadSLisandro Dalcin }
6002e74eeadSLisandro Dalcin 
6012e74eeadSLisandro Dalcin #undef __FUNCT__
602b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
603dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
604b4319ba4SBarry Smith {
605b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
606dfbe8321SBarry Smith   PetscErrorCode ierr;
607b4319ba4SBarry Smith   PetscViewer    sviewer;
608b4319ba4SBarry Smith 
609b4319ba4SBarry Smith   PetscFunctionBegin;
6103f08860eSBarry Smith   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
611b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
6123f08860eSBarry Smith   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
6136e520ac8SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
614b4319ba4SBarry Smith   PetscFunctionReturn(0);
615b4319ba4SBarry Smith }
616b4319ba4SBarry Smith 
617b4319ba4SBarry Smith #undef __FUNCT__
618b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
619784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
620b4319ba4SBarry Smith {
621dfbe8321SBarry Smith   PetscErrorCode ierr;
622e176bc59SStefano Zampini   PetscInt       nr,rbs,nc,cbs;
623b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
624b4319ba4SBarry Smith   IS             from,to;
625e176bc59SStefano Zampini   Vec            cglobal,rglobal;
626b4319ba4SBarry Smith 
627b4319ba4SBarry Smith   PetscFunctionBegin;
628784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
629e176bc59SStefano Zampini   PetscCheckSameComm(A,1,cmapping,3);
6303bbff08aSStefano Zampini   /* Destroy any previous data */
63170cf5478SStefano Zampini   ierr = VecDestroy(&is->x);CHKERRQ(ierr);
63270cf5478SStefano Zampini   ierr = VecDestroy(&is->y);CHKERRQ(ierr);
633e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr);
634e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr);
6351c47cb0fSStefano Zampini   ierr = MatDestroy(&is->A);CHKERRQ(ierr);
63628f4e0baSStefano Zampini   ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
63728f4e0baSStefano Zampini   ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
6383bbff08aSStefano Zampini 
6393bbff08aSStefano Zampini   /* Setup Layout and set local to global maps */
640fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
641fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
642fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
643fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
644b4319ba4SBarry Smith 
645b4319ba4SBarry Smith   /* Create the local matrix A */
646e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr);
647e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr);
648e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr);
649e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr);
650f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
651e176bc59SStefano Zampini   ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr);
652e176bc59SStefano Zampini   ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr);
653e176bc59SStefano Zampini   ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr);
654ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
655ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
656b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
657b4319ba4SBarry Smith 
658b4319ba4SBarry Smith   /* Create the local work vectors */
6593bbff08aSStefano Zampini   ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr);
660b4319ba4SBarry Smith 
661e176bc59SStefano Zampini   /* setup the global to local scatters */
662e176bc59SStefano Zampini   ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr);
663e176bc59SStefano Zampini   ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr);
664784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
665e176bc59SStefano Zampini   ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr);
666e176bc59SStefano Zampini   if (rmapping != cmapping) {
667e176bc59SStefano Zampini     ierr = ISDestroy(&to);CHKERRQ(ierr);
668e176bc59SStefano Zampini     ierr = ISDestroy(&from);CHKERRQ(ierr);
669e176bc59SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr);
670e176bc59SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr);
671e176bc59SStefano Zampini     ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr);
672e176bc59SStefano Zampini   } else {
673e176bc59SStefano Zampini     ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr);
674e176bc59SStefano Zampini     is->cctx = is->rctx;
675e176bc59SStefano Zampini   }
676e176bc59SStefano Zampini   ierr = VecDestroy(&rglobal);CHKERRQ(ierr);
677e176bc59SStefano Zampini   ierr = VecDestroy(&cglobal);CHKERRQ(ierr);
6786bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
6796bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
68048ff6bf3SStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
681b4319ba4SBarry Smith   PetscFunctionReturn(0);
682b4319ba4SBarry Smith }
683b4319ba4SBarry Smith 
6842e74eeadSLisandro Dalcin #undef __FUNCT__
6852e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
6862e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
6872e74eeadSLisandro Dalcin {
6882e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
6892e74eeadSLisandro Dalcin   PetscErrorCode ierr;
690*97563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
691*97563a80SStefano Zampini   PetscInt       i,zm,zn;
692*97563a80SStefano Zampini #endif
693*97563a80SStefano Zampini   PetscInt       rows_l[2048],cols_l[2048];
6942e74eeadSLisandro Dalcin 
6952e74eeadSLisandro Dalcin   PetscFunctionBegin;
6962e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
697f23aa3ddSBarry 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);
698*97563a80SStefano Zampini   /* count negative indices */
699*97563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
700*97563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
7012e74eeadSLisandro Dalcin #endif
702*97563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
703*97563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
704*97563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
705*97563a80SStefano Zampini   /* count negative indices (should be the same as before) */
706*97563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
707*97563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
708*97563a80SStefano Zampini   if (zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the row indices can not be mapped! Maybe you should not use MATIS");
709*97563a80SStefano Zampini   if (zn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the column indices can not be mapped! Maybe you should not use MATIS");
710*97563a80SStefano Zampini #endif
7112e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
7122e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7132e74eeadSLisandro Dalcin }
7142e74eeadSLisandro Dalcin 
715b4319ba4SBarry Smith #undef __FUNCT__
716*97563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS"
717*97563a80SStefano Zampini PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
718*97563a80SStefano Zampini {
719*97563a80SStefano Zampini   Mat_IS         *is = (Mat_IS*)mat->data;
720*97563a80SStefano Zampini   PetscErrorCode ierr;
721*97563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
722*97563a80SStefano Zampini   PetscInt       i,zm,zn;
723*97563a80SStefano Zampini #endif
724*97563a80SStefano Zampini   PetscInt       rows_l[2048],cols_l[2048];
725*97563a80SStefano Zampini 
726*97563a80SStefano Zampini   PetscFunctionBegin;
727*97563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
728*97563a80SStefano Zampini   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);
729*97563a80SStefano Zampini   /* count negative indices */
730*97563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
731*97563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
732*97563a80SStefano Zampini #endif
733*97563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
734*97563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
735*97563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
736*97563a80SStefano Zampini   /* count negative indices (should be the same as before) */
737*97563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
738*97563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
739*97563a80SStefano Zampini   if (zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the row indices can not be mapped! Maybe you should not use MATIS");
740*97563a80SStefano Zampini   if (zn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the column indices can not be mapped! Maybe you should not use MATIS");
741*97563a80SStefano Zampini #endif
742*97563a80SStefano Zampini   ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
743*97563a80SStefano Zampini   PetscFunctionReturn(0);
744*97563a80SStefano Zampini }
745*97563a80SStefano Zampini 
746*97563a80SStefano Zampini #undef __FUNCT__
747b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
74813f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
749b4319ba4SBarry Smith {
750dfbe8321SBarry Smith   PetscErrorCode ierr;
751b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
752b4319ba4SBarry Smith 
753b4319ba4SBarry Smith   PetscFunctionBegin;
754b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
755b4319ba4SBarry Smith   PetscFunctionReturn(0);
756b4319ba4SBarry Smith }
757b4319ba4SBarry Smith 
758b4319ba4SBarry Smith #undef __FUNCT__
759f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
760f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
761f0006bf2SLisandro Dalcin {
762f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
763f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
764f0006bf2SLisandro Dalcin 
765f0006bf2SLisandro Dalcin   PetscFunctionBegin;
766f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
767f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
768f0006bf2SLisandro Dalcin }
769f0006bf2SLisandro Dalcin 
770f0006bf2SLisandro Dalcin #undef __FUNCT__
7712e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
7722b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
7732e74eeadSLisandro Dalcin {
7746e520ac8SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
7756e520ac8SStefano Zampini   PetscInt       nr,nl,len,i;
7766e520ac8SStefano Zampini   PetscInt       *lrows;
7772e74eeadSLisandro Dalcin   PetscErrorCode ierr;
7782e74eeadSLisandro Dalcin 
7792e74eeadSLisandro Dalcin   PetscFunctionBegin;
7806e520ac8SStefano Zampini   /* get locally owned rows */
7816e520ac8SStefano Zampini   ierr = MatZeroRowsMapLocal_Private(A,n,rows,&len,&lrows);CHKERRQ(ierr);
7826e520ac8SStefano Zampini   /* fix right hand side if needed */
7836e520ac8SStefano Zampini   if (x && b) {
7846e520ac8SStefano Zampini     const PetscScalar *xx;
7856e520ac8SStefano Zampini     PetscScalar       *bb;
7866e520ac8SStefano Zampini 
7876e520ac8SStefano Zampini     ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr);
7886e520ac8SStefano Zampini     ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
7896e520ac8SStefano Zampini     for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]];
7906e520ac8SStefano Zampini     ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr);
7916e520ac8SStefano Zampini     ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr);
7922e74eeadSLisandro Dalcin   }
7936e520ac8SStefano Zampini   /* get rows associated to the local matrices */
7946e520ac8SStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
7956e520ac8SStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
7966e520ac8SStefano Zampini   }
7976e520ac8SStefano Zampini   ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr);
7986e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr);
7996e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
8006e520ac8SStefano Zampini   for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1;
8016e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
8026e520ac8SStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
8036e520ac8SStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
8046e520ac8SStefano Zampini   ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr);
8056e520ac8SStefano Zampini   for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i;
8067230de76SStefano Zampini   ierr = MatISZeroRowsLocal_Private(A,nr,lrows,diag);CHKERRQ(ierr);
8076e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
8082e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
8092e74eeadSLisandro Dalcin }
8102e74eeadSLisandro Dalcin 
8112e74eeadSLisandro Dalcin #undef __FUNCT__
8127230de76SStefano Zampini #define __FUNCT__ "MatISZeroRowsLocal_Private"
8137230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag)
814b4319ba4SBarry Smith {
815b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
816dfbe8321SBarry Smith   PetscErrorCode ierr;
817b4319ba4SBarry Smith 
818b4319ba4SBarry Smith   PetscFunctionBegin;
8196e520ac8SStefano Zampini   if (diag != 0.) {
820b4319ba4SBarry Smith     /*
821b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
822b4319ba4SBarry Smith        work properly in the interface nodes.
823b4319ba4SBarry Smith     */
824b4319ba4SBarry Smith     Vec counter;
825e176bc59SStefano Zampini     ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr);
826e176bc59SStefano Zampini     ierr = VecSet(counter,0.);CHKERRQ(ierr);
827e176bc59SStefano Zampini     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
828e176bc59SStefano Zampini     ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
829e176bc59SStefano Zampini     ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
830e176bc59SStefano Zampini     ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
831e176bc59SStefano Zampini     ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8326bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
833b4319ba4SBarry Smith   }
834958c9bccSBarry Smith   if (!n) {
835b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
836b4319ba4SBarry Smith   } else {
8377230de76SStefano Zampini     PetscInt i;
838b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
8392205254eSKarl Rupp 
8402b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
8416e520ac8SStefano Zampini     if (diag != 0.) {
8426e520ac8SStefano Zampini       const PetscScalar *array;
8436e520ac8SStefano Zampini       ierr = VecGetArrayRead(is->y,&array);CHKERRQ(ierr);
844b4319ba4SBarry Smith       for (i=0; i<n; i++) {
845f4df32b1SMatthew Knepley         ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
846b4319ba4SBarry Smith       }
8476e520ac8SStefano Zampini       ierr = VecRestoreArrayRead(is->y,&array);CHKERRQ(ierr);
8486e520ac8SStefano Zampini     }
849b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
850b4319ba4SBarry Smith     ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
851b4319ba4SBarry Smith   }
852b4319ba4SBarry Smith   PetscFunctionReturn(0);
853b4319ba4SBarry Smith }
854b4319ba4SBarry Smith 
855b4319ba4SBarry Smith #undef __FUNCT__
856b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
857dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
858b4319ba4SBarry Smith {
859b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
860dfbe8321SBarry Smith   PetscErrorCode ierr;
861b4319ba4SBarry Smith 
862b4319ba4SBarry Smith   PetscFunctionBegin;
863b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
864b4319ba4SBarry Smith   PetscFunctionReturn(0);
865b4319ba4SBarry Smith }
866b4319ba4SBarry Smith 
867b4319ba4SBarry Smith #undef __FUNCT__
868b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
869dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
870b4319ba4SBarry Smith {
871b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
872dfbe8321SBarry Smith   PetscErrorCode ierr;
873b4319ba4SBarry Smith 
874b4319ba4SBarry Smith   PetscFunctionBegin;
875b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
876b4319ba4SBarry Smith   PetscFunctionReturn(0);
877b4319ba4SBarry Smith }
878b4319ba4SBarry Smith 
879b4319ba4SBarry Smith #undef __FUNCT__
880b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
8817087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
882b4319ba4SBarry Smith {
883b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
884b4319ba4SBarry Smith 
885b4319ba4SBarry Smith   PetscFunctionBegin;
886b4319ba4SBarry Smith   *local = is->A;
887b4319ba4SBarry Smith   PetscFunctionReturn(0);
888b4319ba4SBarry Smith }
889b4319ba4SBarry Smith 
890b4319ba4SBarry Smith #undef __FUNCT__
891b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
892b4319ba4SBarry Smith /*@
893b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
894b4319ba4SBarry Smith 
895b4319ba4SBarry Smith   Input Parameter:
896b4319ba4SBarry Smith .  mat - the matrix
897b4319ba4SBarry Smith 
898b4319ba4SBarry Smith   Output Parameter:
899eb82efa4SStefano Zampini .  local - the local matrix
900b4319ba4SBarry Smith 
901b4319ba4SBarry Smith   Level: advanced
902b4319ba4SBarry Smith 
903b4319ba4SBarry Smith   Notes:
904b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
905b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
906b4319ba4SBarry Smith   of the MatSetValues() operation.
907b4319ba4SBarry Smith 
908b4319ba4SBarry Smith .seealso: MATIS
909b4319ba4SBarry Smith @*/
9107087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
911b4319ba4SBarry Smith {
9124ac538c5SBarry Smith   PetscErrorCode ierr;
913b4319ba4SBarry Smith 
914b4319ba4SBarry Smith   PetscFunctionBegin;
9150700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
916b4319ba4SBarry Smith   PetscValidPointer(local,2);
9174ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
918b4319ba4SBarry Smith   PetscFunctionReturn(0);
919b4319ba4SBarry Smith }
920b4319ba4SBarry Smith 
9213b03a366Sstefano_zampini #undef __FUNCT__
9223b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
9233b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
9243b03a366Sstefano_zampini {
9253b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
9263b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
9273b03a366Sstefano_zampini   PetscErrorCode ierr;
9283b03a366Sstefano_zampini 
9293b03a366Sstefano_zampini   PetscFunctionBegin;
9304e4c7dbeSStefano Zampini   if (is->A) {
9313b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
9323b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
933f23aa3ddSBarry 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);
9344e4c7dbeSStefano Zampini   }
9353b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
9363b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
9373b03a366Sstefano_zampini   is->A = local;
9383b03a366Sstefano_zampini   PetscFunctionReturn(0);
9393b03a366Sstefano_zampini }
9403b03a366Sstefano_zampini 
9413b03a366Sstefano_zampini #undef __FUNCT__
9423b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
9433b03a366Sstefano_zampini /*@
944eb82efa4SStefano Zampini     MatISSetLocalMat - Replace the local matrix stored inside a MATIS object.
9453b03a366Sstefano_zampini 
9463b03a366Sstefano_zampini   Input Parameter:
9473b03a366Sstefano_zampini .  mat - the matrix
948eb82efa4SStefano Zampini .  local - the local matrix
9493b03a366Sstefano_zampini 
9503b03a366Sstefano_zampini   Output Parameter:
9513b03a366Sstefano_zampini 
9523b03a366Sstefano_zampini   Level: advanced
9533b03a366Sstefano_zampini 
9543b03a366Sstefano_zampini   Notes:
9553b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
9563b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
9573b03a366Sstefano_zampini 
9583b03a366Sstefano_zampini .seealso: MATIS
9593b03a366Sstefano_zampini @*/
9603b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
9613b03a366Sstefano_zampini {
9623b03a366Sstefano_zampini   PetscErrorCode ierr;
9633b03a366Sstefano_zampini 
9643b03a366Sstefano_zampini   PetscFunctionBegin;
9653b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
966b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
9673b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
9683b03a366Sstefano_zampini   PetscFunctionReturn(0);
9693b03a366Sstefano_zampini }
9703b03a366Sstefano_zampini 
9716726f965SBarry Smith #undef __FUNCT__
9726726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
9736726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
9746726f965SBarry Smith {
9756726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
9766726f965SBarry Smith   PetscErrorCode ierr;
9776726f965SBarry Smith 
9786726f965SBarry Smith   PetscFunctionBegin;
9796726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
9806726f965SBarry Smith   PetscFunctionReturn(0);
9816726f965SBarry Smith }
9826726f965SBarry Smith 
9836726f965SBarry Smith #undef __FUNCT__
9842e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
9852e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
9862e74eeadSLisandro Dalcin {
9872e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9882e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9892e74eeadSLisandro Dalcin 
9902e74eeadSLisandro Dalcin   PetscFunctionBegin;
9912e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
9922e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9932e74eeadSLisandro Dalcin }
9942e74eeadSLisandro Dalcin 
9952e74eeadSLisandro Dalcin #undef __FUNCT__
9962e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
9972e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
9982e74eeadSLisandro Dalcin {
9992e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
10002e74eeadSLisandro Dalcin   PetscErrorCode ierr;
10012e74eeadSLisandro Dalcin 
10022e74eeadSLisandro Dalcin   PetscFunctionBegin;
10032e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
1004e176bc59SStefano Zampini   ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr);
10052e74eeadSLisandro Dalcin 
10062e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
10072e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
1008e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1009e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10102e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
10112e74eeadSLisandro Dalcin }
10122e74eeadSLisandro Dalcin 
10132e74eeadSLisandro Dalcin #undef __FUNCT__
10146726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
1015ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
10166726f965SBarry Smith {
10176726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
10186726f965SBarry Smith   PetscErrorCode ierr;
10196726f965SBarry Smith 
10206726f965SBarry Smith   PetscFunctionBegin;
10214e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
10226726f965SBarry Smith   PetscFunctionReturn(0);
10236726f965SBarry Smith }
10246726f965SBarry Smith 
1025284134d9SBarry Smith #undef __FUNCT__
1026284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
1027284134d9SBarry Smith /*@
10283c212e90SHong Zhang     MatCreateIS - Creates a "process" unassembled matrix, assembled on each
1029284134d9SBarry Smith        process but not across processes.
1030284134d9SBarry Smith 
1031284134d9SBarry Smith    Input Parameters:
1032284134d9SBarry Smith +     comm    - MPI communicator that will share the matrix
1033e176bc59SStefano Zampini .     bs      - block size of the matrix
1034df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
1035e176bc59SStefano Zampini .     rmap    - local to global map for rows
1036e176bc59SStefano Zampini -     cmap    - local to global map for cols
1037284134d9SBarry Smith 
1038284134d9SBarry Smith    Output Parameter:
1039284134d9SBarry Smith .    A - the resulting matrix
1040284134d9SBarry Smith 
10418e6c10adSSatish Balay    Level: advanced
10428e6c10adSSatish Balay 
10433c212e90SHong Zhang    Notes: See MATIS for more details.
10443c212e90SHong Zhang           m and n are NOT related to the size of the map; they are the size of the part of the vector owned
1045e176bc59SStefano Zampini           by that process. The sizes of rmap and cmap define the size of the local matrices.
10463c212e90SHong Zhang           If either rmap or cmap are NULL, then the matrix is assumed to be square.
1047284134d9SBarry Smith 
1048284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
1049284134d9SBarry Smith @*/
1050e176bc59SStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A)
1051284134d9SBarry Smith {
1052284134d9SBarry Smith   PetscErrorCode ierr;
1053284134d9SBarry Smith 
1054284134d9SBarry Smith   PetscFunctionBegin;
1055e176bc59SStefano Zampini   if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping");
1056284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1057d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
1058284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1059284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
1060e176bc59SStefano Zampini   if (rmap && cmap) {
1061e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr);
1062e176bc59SStefano Zampini   } else if (!rmap) {
1063e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr);
1064e176bc59SStefano Zampini   } else {
1065e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr);
1066e176bc59SStefano Zampini   }
1067284134d9SBarry Smith   PetscFunctionReturn(0);
1068284134d9SBarry Smith }
1069284134d9SBarry Smith 
1070b4319ba4SBarry Smith /*MC
1071eb82efa4SStefano Zampini    MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC).
1072b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
1073b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
1074b4319ba4SBarry Smith    product is handled "implicitly".
1075b4319ba4SBarry Smith 
1076b4319ba4SBarry Smith    Operations Provided:
10776726f965SBarry Smith +  MatMult()
10782e74eeadSLisandro Dalcin .  MatMultAdd()
10792e74eeadSLisandro Dalcin .  MatMultTranspose()
10802e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
10816726f965SBarry Smith .  MatZeroEntries()
10826726f965SBarry Smith .  MatSetOption()
10832e74eeadSLisandro Dalcin .  MatZeroRows()
10842e74eeadSLisandro Dalcin .  MatSetValues()
1085*97563a80SStefano Zampini .  MatSetValuesBlocked()
10866726f965SBarry Smith .  MatSetValuesLocal()
1087*97563a80SStefano Zampini .  MatSetValuesBlockedLocal()
10882e74eeadSLisandro Dalcin .  MatScale()
10892e74eeadSLisandro Dalcin .  MatGetDiagonal()
10906726f965SBarry Smith -  MatSetLocalToGlobalMapping()
1091b4319ba4SBarry Smith 
1092b4319ba4SBarry Smith    Options Database Keys:
1093b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
1094b4319ba4SBarry Smith 
1095b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
1096b4319ba4SBarry Smith 
1097b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
1098b4319ba4SBarry Smith 
1099b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
1100eb82efa4SStefano Zampini           MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation()
1101b4319ba4SBarry Smith 
1102b4319ba4SBarry Smith   Level: advanced
1103b4319ba4SBarry Smith 
11043c212e90SHong Zhang .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC, MatCreateIS()
1105b4319ba4SBarry Smith 
1106b4319ba4SBarry Smith M*/
1107b4319ba4SBarry Smith 
1108b4319ba4SBarry Smith #undef __FUNCT__
1109b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
11108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
1111b4319ba4SBarry Smith {
1112dfbe8321SBarry Smith   PetscErrorCode ierr;
1113b4319ba4SBarry Smith   Mat_IS         *b;
1114b4319ba4SBarry Smith 
1115b4319ba4SBarry Smith   PetscFunctionBegin;
1116b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1117b4319ba4SBarry Smith   A->data = (void*)b;
1118b4319ba4SBarry Smith 
1119e176bc59SStefano Zampini   /* matrix ops */
1120e176bc59SStefano Zampini   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1121b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
11222e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
11232e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
11242e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1125b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1126b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
11272e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
1128*97563a80SStefano Zampini   A->ops->setvalues               = MatSetValuesBlocked_IS;
1129b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1130f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
11312e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1132b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1133b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1134b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
11356726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
11362e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
11372e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
11386726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
113969796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
114069796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1141ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
1142b4319ba4SBarry Smith 
1143b7ce53b6SStefano Zampini   /* special MATIS functions */
1144bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1145bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1146b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
11472e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
114817667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1149b4319ba4SBarry Smith   PetscFunctionReturn(0);
1150b4319ba4SBarry Smith }
1151