xref: /petsc/src/mat/impls/is/matis.c (revision 6bd84002eb5dfb07ec36d884cbfaa6037897c44f)
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__
17*6bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS"
18*6bd84002SStefano Zampini PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool  *missing,PetscInt *d)
19*6bd84002SStefano Zampini {
20*6bd84002SStefano Zampini   Mat_IS         *a = (Mat_IS*)A->data;
21*6bd84002SStefano Zampini   PetscErrorCode ierr;
22*6bd84002SStefano Zampini 
23*6bd84002SStefano Zampini   PetscFunctionBegin;
24*6bd84002SStefano Zampini   if (d) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Need to be implemented");
25*6bd84002SStefano Zampini   ierr = MatMissingDiagonal(a->A,missing,NULL);CHKERRQ(ierr);
26*6bd84002SStefano Zampini   PetscFunctionReturn(0);
27*6bd84002SStefano Zampini }
28*6bd84002SStefano Zampini 
29*6bd84002SStefano Zampini #undef __FUNCT__
3028f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private"
31a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B)
3228f4e0baSStefano Zampini {
3328f4e0baSStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
3428f4e0baSStefano Zampini   const PetscInt *gidxs;
3528f4e0baSStefano Zampini   PetscErrorCode ierr;
3628f4e0baSStefano Zampini 
3728f4e0baSStefano Zampini   PetscFunctionBegin;
3828f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr);
3928f4e0baSStefano Zampini   ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr);
4028f4e0baSStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr);
413bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
4228f4e0baSStefano Zampini   /* PETSC_OWN_POINTER refers to ilocal which is NULL */
4328f4e0baSStefano Zampini   ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
443bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
4528f4e0baSStefano Zampini   ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr);
4628f4e0baSStefano Zampini   PetscFunctionReturn(0);
4728f4e0baSStefano Zampini }
482e1947a5SStefano Zampini 
492e1947a5SStefano Zampini #undef __FUNCT__
502e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation"
51eb82efa4SStefano Zampini /*@
52a88811baSStefano Zampini    MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix.
53a88811baSStefano Zampini 
54a88811baSStefano Zampini    Collective on MPI_Comm
55a88811baSStefano Zampini 
56a88811baSStefano Zampini    Input Parameters:
57a88811baSStefano Zampini +  B - the matrix
58a88811baSStefano Zampini .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
59a88811baSStefano Zampini            (same value is used for all local rows)
60a88811baSStefano Zampini .  d_nnz - array containing the number of nonzeros in the various rows of the
61a88811baSStefano Zampini            DIAGONAL portion of the local submatrix (possibly different for each row)
62a88811baSStefano Zampini            or NULL, if d_nz is used to specify the nonzero structure.
63a88811baSStefano Zampini            The size of this array is equal to the number of local rows, i.e 'm'.
64a88811baSStefano Zampini            For matrices that will be factored, you must leave room for (and set)
65a88811baSStefano Zampini            the diagonal entry even if it is zero.
66a88811baSStefano Zampini .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
67a88811baSStefano Zampini            submatrix (same value is used for all local rows).
68a88811baSStefano Zampini -  o_nnz - array containing the number of nonzeros in the various rows of the
69a88811baSStefano Zampini            OFF-DIAGONAL portion of the local submatrix (possibly different for
70a88811baSStefano Zampini            each row) or NULL, if o_nz is used to specify the nonzero
71a88811baSStefano Zampini            structure. The size of this array is equal to the number
72a88811baSStefano Zampini            of local rows, i.e 'm'.
73a88811baSStefano Zampini 
74a88811baSStefano Zampini    If the *_nnz parameter is given then the *_nz parameter is ignored
75a88811baSStefano Zampini 
76a88811baSStefano Zampini    Level: intermediate
77a88811baSStefano Zampini 
78a88811baSStefano Zampini    Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition
79a88811baSStefano Zampini           from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local
80a88811baSStefano Zampini           matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects.
81a88811baSStefano Zampini 
82a88811baSStefano Zampini .keywords: matrix
83a88811baSStefano Zampini 
843c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS
85a88811baSStefano Zampini @*/
862e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
872e1947a5SStefano Zampini {
882e1947a5SStefano Zampini   PetscErrorCode ierr;
892e1947a5SStefano Zampini 
902e1947a5SStefano Zampini   PetscFunctionBegin;
912e1947a5SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
922e1947a5SStefano Zampini   PetscValidType(B,1);
932e1947a5SStefano Zampini   ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
942e1947a5SStefano Zampini   PetscFunctionReturn(0);
952e1947a5SStefano Zampini }
962e1947a5SStefano Zampini 
972e1947a5SStefano Zampini #undef __FUNCT__
982e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS"
997230de76SStefano Zampini static PetscErrorCode  MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
1002e1947a5SStefano Zampini {
1012e1947a5SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
10228f4e0baSStefano Zampini   PetscInt       bs,i,nlocalcols;
1032e1947a5SStefano Zampini   PetscErrorCode ierr;
1042e1947a5SStefano Zampini 
1052e1947a5SStefano Zampini   PetscFunctionBegin;
1066c4ed002SBarry Smith   if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
10728f4e0baSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
10828f4e0baSStefano Zampini     ierr = MatISComputeSF_Private(B);CHKERRQ(ierr);
10928f4e0baSStefano Zampini   }
1102e1947a5SStefano Zampini   if (!d_nnz) {
11128f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz;
1122e1947a5SStefano Zampini   } else {
11328f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i];
1142e1947a5SStefano Zampini   }
1152e1947a5SStefano Zampini   if (!o_nnz) {
11628f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz;
1172e1947a5SStefano Zampini   } else {
11828f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i];
1192e1947a5SStefano Zampini   }
12028f4e0baSStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
12128f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr);
12228f4e0baSStefano Zampini   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
12328f4e0baSStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
12428f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves;i++) {
12528f4e0baSStefano Zampini     matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols);
1262e1947a5SStefano Zampini   }
12728f4e0baSStefano Zampini   ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr);
12828f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
12928f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs;
1302e1947a5SStefano Zampini   }
13128f4e0baSStefano Zampini   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
13228f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
13328f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i]-i;
1342e1947a5SStefano Zampini   }
13528f4e0baSStefano Zampini   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
1362e1947a5SStefano Zampini   PetscFunctionReturn(0);
1372e1947a5SStefano Zampini }
138b4319ba4SBarry Smith 
139b4319ba4SBarry Smith #undef __FUNCT__
1403927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private"
1413927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
1423927de2eSStefano Zampini {
1433927de2eSStefano Zampini   Mat_IS          *matis = (Mat_IS*)(A->data);
1443927de2eSStefano Zampini   PetscInt        *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership;
145ecf5a873SStefano Zampini   const PetscInt  *global_indices_r,*global_indices_c;
1463927de2eSStefano Zampini   PetscInt        i,j,bs,rows,cols;
1473927de2eSStefano Zampini   PetscInt        lrows,lcols;
1483927de2eSStefano Zampini   PetscInt        local_rows,local_cols;
1493927de2eSStefano Zampini   PetscMPIInt     nsubdomains;
1503927de2eSStefano Zampini   PetscBool       isdense,issbaij;
1513927de2eSStefano Zampini   PetscErrorCode  ierr;
1523927de2eSStefano Zampini 
1533927de2eSStefano Zampini   PetscFunctionBegin;
1543927de2eSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr);
1553927de2eSStefano Zampini   ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr);
1563927de2eSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1573927de2eSStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
1583927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
1593927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
160ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
161ecf5a873SStefano Zampini   if (A->rmap->mapping != A->cmap->mapping) {
1627230de76SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
163ecf5a873SStefano Zampini   } else {
164ecf5a873SStefano Zampini     global_indices_c = global_indices_r;
165ecf5a873SStefano Zampini   }
166ecf5a873SStefano Zampini 
1673927de2eSStefano Zampini   if (issbaij) {
1683927de2eSStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
1693927de2eSStefano Zampini   }
1703927de2eSStefano Zampini   /*
171ecf5a873SStefano Zampini      An SF reduce is needed to sum up properly on shared rows.
1723927de2eSStefano Zampini      Note that generally preallocation is not exact, since it overestimates nonzeros
1733927de2eSStefano Zampini   */
1743927de2eSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
1753927de2eSStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
1763927de2eSStefano Zampini   }
1773927de2eSStefano Zampini   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
1783927de2eSStefano Zampini   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
1793927de2eSStefano Zampini   /* All processes need to compute entire row ownership */
1803927de2eSStefano Zampini   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
1813927de2eSStefano Zampini   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
1823927de2eSStefano Zampini   for (i=0;i<nsubdomains;i++) {
1833927de2eSStefano Zampini     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
1843927de2eSStefano Zampini       row_ownership[j] = i;
1853927de2eSStefano Zampini     }
1863927de2eSStefano Zampini   }
1877230de76SStefano Zampini   ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
1883927de2eSStefano Zampini 
1893927de2eSStefano Zampini   /*
1903927de2eSStefano Zampini      my_dnz and my_onz contains exact contribution to preallocation from each local mat
1913927de2eSStefano Zampini      then, they will be summed up properly. This way, preallocation is always sufficient
1923927de2eSStefano Zampini   */
1933927de2eSStefano Zampini   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
1943927de2eSStefano Zampini   /* preallocation as a MATAIJ */
1953927de2eSStefano Zampini   if (isdense) { /* special case for dense local matrices */
1963927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
197ecf5a873SStefano Zampini       PetscInt index_row = global_indices_r[i];
1983927de2eSStefano Zampini       for (j=i;j<local_rows;j++) {
1993927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
200ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[j];
2013927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
2023927de2eSStefano Zampini           my_dnz[i] += 1;
2033927de2eSStefano Zampini         } else { /* offdiag block */
2043927de2eSStefano Zampini           my_onz[i] += 1;
2053927de2eSStefano Zampini         }
2063927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
2073927de2eSStefano Zampini         if (i != j) {
2083927de2eSStefano Zampini           owner = row_ownership[index_col];
2093927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
2103927de2eSStefano Zampini             my_dnz[j] += 1;
2113927de2eSStefano Zampini           } else {
2123927de2eSStefano Zampini             my_onz[j] += 1;
2133927de2eSStefano Zampini           }
2143927de2eSStefano Zampini         }
2153927de2eSStefano Zampini       }
2163927de2eSStefano Zampini     }
217ecf5a873SStefano Zampini   } else { /* TODO: this could be optimized using MatGetRowIJ */
2183927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
2193927de2eSStefano Zampini       const PetscInt *cols;
220ecf5a873SStefano Zampini       PetscInt       ncols,index_row = global_indices_r[i];
2213927de2eSStefano Zampini       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2223927de2eSStefano Zampini       for (j=0;j<ncols;j++) {
2233927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
224ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[cols[j]];
2253927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
2263927de2eSStefano Zampini           my_dnz[i] += 1;
2273927de2eSStefano Zampini         } else { /* offdiag block */
2283927de2eSStefano Zampini           my_onz[i] += 1;
2293927de2eSStefano Zampini         }
2303927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
231d9a9e74cSStefano Zampini         if (issbaij && index_col != index_row) {
2323927de2eSStefano Zampini           owner = row_ownership[index_col];
2333927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
234d9a9e74cSStefano Zampini             my_dnz[cols[j]] += 1;
2353927de2eSStefano Zampini           } else {
236d9a9e74cSStefano Zampini             my_onz[cols[j]] += 1;
2373927de2eSStefano Zampini           }
2383927de2eSStefano Zampini         }
2393927de2eSStefano Zampini       }
2403927de2eSStefano Zampini       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2413927de2eSStefano Zampini     }
2423927de2eSStefano Zampini   }
243ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
244ecf5a873SStefano Zampini   if (global_indices_c != global_indices_r) {
2457230de76SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
246ecf5a873SStefano Zampini   }
2473927de2eSStefano Zampini   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
248ecf5a873SStefano Zampini 
249ecf5a873SStefano Zampini   /* Reduce my_dnz and my_onz */
2503927de2eSStefano Zampini   if (maxreduce) {
2513927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2523927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2533927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2543927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2553927de2eSStefano Zampini   } else {
2563927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2573927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2583927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2593927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2603927de2eSStefano Zampini   }
2613927de2eSStefano Zampini   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
2623927de2eSStefano Zampini 
2633927de2eSStefano Zampini   /* Resize preallocation if overestimated */
2643927de2eSStefano Zampini   for (i=0;i<lrows;i++) {
2653927de2eSStefano Zampini     dnz[i] = PetscMin(dnz[i],lcols);
2663927de2eSStefano Zampini     onz[i] = PetscMin(onz[i],cols-lcols);
2673927de2eSStefano Zampini   }
2683927de2eSStefano Zampini   /* set preallocation */
2693927de2eSStefano Zampini   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
2703927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
2713927de2eSStefano Zampini     dnz[i] = dnz[i*bs]/bs;
2723927de2eSStefano Zampini     onz[i] = onz[i*bs]/bs;
2733927de2eSStefano Zampini   }
2743927de2eSStefano Zampini   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2753927de2eSStefano Zampini   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2763927de2eSStefano Zampini   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
2773927de2eSStefano Zampini   if (issbaij) {
2783927de2eSStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
2793927de2eSStefano Zampini   }
2803927de2eSStefano Zampini   PetscFunctionReturn(0);
2813927de2eSStefano Zampini }
2823927de2eSStefano Zampini 
2833927de2eSStefano Zampini #undef __FUNCT__
284b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
2857230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
286b7ce53b6SStefano Zampini {
287b7ce53b6SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
288d9a9e74cSStefano Zampini   Mat            local_mat;
289b7ce53b6SStefano Zampini   /* info on mat */
2903cfa4ea4SStefano Zampini   PetscInt       bs,rows,cols,lrows,lcols;
291b7ce53b6SStefano Zampini   PetscInt       local_rows,local_cols;
292686e3a49SStefano Zampini   PetscBool      isdense,issbaij,isseqaij;
2937c03b4e8SStefano Zampini   PetscMPIInt    nsubdomains;
294b7ce53b6SStefano Zampini   /* values insertion */
295b7ce53b6SStefano Zampini   PetscScalar    *array;
296b7ce53b6SStefano Zampini   /* work */
297b7ce53b6SStefano Zampini   PetscErrorCode ierr;
298b7ce53b6SStefano Zampini 
299b7ce53b6SStefano Zampini   PetscFunctionBegin;
300b7ce53b6SStefano Zampini   /* get info from mat */
3017c03b4e8SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
3027c03b4e8SStefano Zampini   if (nsubdomains == 1) {
3037c03b4e8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
3047c03b4e8SStefano Zampini       ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr);
3057c03b4e8SStefano Zampini     } else {
3067c03b4e8SStefano Zampini       ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
3077c03b4e8SStefano Zampini     }
3087c03b4e8SStefano Zampini     PetscFunctionReturn(0);
3097c03b4e8SStefano Zampini   }
310b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
311b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
3123cfa4ea4SStefano Zampini   ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr);
313b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
314b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
315686e3a49SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
316b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
317b7ce53b6SStefano Zampini 
318b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
319b7ce53b6SStefano Zampini     MatType     new_mat_type;
3203927de2eSStefano Zampini     PetscBool   issbaij_red;
321b7ce53b6SStefano Zampini 
322b7ce53b6SStefano Zampini     /* determining new matrix type */
323b2566f29SBarry Smith     ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
324b7ce53b6SStefano Zampini     if (issbaij_red) {
325b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
326b7ce53b6SStefano Zampini     } else {
327b7ce53b6SStefano Zampini       if (bs>1) {
328b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
329b7ce53b6SStefano Zampini       } else {
330b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
331b7ce53b6SStefano Zampini       }
332b7ce53b6SStefano Zampini     }
333b7ce53b6SStefano Zampini 
3343927de2eSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
3353cfa4ea4SStefano Zampini     ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr);
3363927de2eSStefano Zampini     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
3373927de2eSStefano Zampini     ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr);
3383927de2eSStefano Zampini     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
339b7ce53b6SStefano Zampini   } else {
3403cfa4ea4SStefano Zampini     PetscInt mbs,mrows,mcols,mlrows,mlcols;
341b7ce53b6SStefano Zampini     /* some checks */
342b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
343b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
3443cfa4ea4SStefano Zampini     ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr);
3456c4ed002SBarry Smith     if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
3466c4ed002SBarry Smith     if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
3476c4ed002SBarry Smith     if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows);
3486c4ed002SBarry Smith     if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols);
3496c4ed002SBarry Smith     if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
350b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
351b7ce53b6SStefano Zampini   }
352d9a9e74cSStefano Zampini 
353d9a9e74cSStefano Zampini   if (issbaij) {
354d9a9e74cSStefano Zampini     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr);
355d9a9e74cSStefano Zampini   } else {
356d9a9e74cSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
357d9a9e74cSStefano Zampini     local_mat = matis->A;
358d9a9e74cSStefano Zampini   }
359686e3a49SStefano Zampini 
360b7ce53b6SStefano Zampini   /* Set values */
361ecf5a873SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr);
362b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
363ecf5a873SStefano Zampini     PetscInt i,*dummy_rows,*dummy_cols;
364ecf5a873SStefano Zampini 
365ecf5a873SStefano Zampini     if (local_rows != local_cols) {
366ecf5a873SStefano Zampini       ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr);
367ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
368ecf5a873SStefano Zampini       for (i=0;i<local_cols;i++) dummy_cols[i] = i;
369ecf5a873SStefano Zampini     } else {
370ecf5a873SStefano Zampini       ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr);
371ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
372ecf5a873SStefano Zampini       dummy_cols = dummy_rows;
373ecf5a873SStefano Zampini     }
374b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
375d9a9e74cSStefano Zampini     ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr);
376ecf5a873SStefano Zampini     ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr);
377d9a9e74cSStefano Zampini     ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr);
378ecf5a873SStefano Zampini     if (dummy_rows != dummy_cols) {
379ecf5a873SStefano Zampini       ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr);
380ecf5a873SStefano Zampini     } else {
381ecf5a873SStefano Zampini       ierr = PetscFree(dummy_rows);CHKERRQ(ierr);
382ecf5a873SStefano Zampini     }
383686e3a49SStefano Zampini   } else if (isseqaij) {
384ecf5a873SStefano Zampini     PetscInt  i,nvtxs,*xadj,*adjncy;
385686e3a49SStefano Zampini     PetscBool done;
386686e3a49SStefano Zampini 
387d9a9e74cSStefano Zampini     ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
3886c4ed002SBarry Smith     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
389d9a9e74cSStefano Zampini     ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr);
390686e3a49SStefano Zampini     for (i=0;i<nvtxs;i++) {
391ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr);
392686e3a49SStefano Zampini     }
393d9a9e74cSStefano Zampini     ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
3946c4ed002SBarry Smith     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
395d9a9e74cSStefano Zampini     ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr);
396686e3a49SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
397ecf5a873SStefano Zampini     PetscInt i;
398c0962df8SStefano Zampini 
399686e3a49SStefano Zampini     for (i=0;i<local_rows;i++) {
400686e3a49SStefano Zampini       PetscInt       j;
401ecf5a873SStefano Zampini       const PetscInt *local_indices_cols;
402686e3a49SStefano Zampini 
403ecf5a873SStefano Zampini       ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
404ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr);
405ecf5a873SStefano Zampini       ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
406686e3a49SStefano Zampini     }
407b7ce53b6SStefano Zampini   }
408b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
409d9a9e74cSStefano Zampini   ierr = MatDestroy(&local_mat);CHKERRQ(ierr);
410b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
411b7ce53b6SStefano Zampini   if (isdense) {
412b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
413b7ce53b6SStefano Zampini   }
414b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
415b7ce53b6SStefano Zampini }
416b7ce53b6SStefano Zampini 
417b7ce53b6SStefano Zampini #undef __FUNCT__
418b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
419b7ce53b6SStefano Zampini /*@
420b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
421b7ce53b6SStefano Zampini 
422b7ce53b6SStefano Zampini   Input Parameter:
423b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
424b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
425b7ce53b6SStefano Zampini 
426b7ce53b6SStefano Zampini   Output Parameter:
427b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
428b7ce53b6SStefano Zampini 
429b7ce53b6SStefano Zampini   Level: developer
430b7ce53b6SStefano Zampini 
431eb82efa4SStefano Zampini   Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested.
432b7ce53b6SStefano Zampini 
433b7ce53b6SStefano Zampini .seealso: MATIS
434b7ce53b6SStefano Zampini @*/
435b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
436b7ce53b6SStefano Zampini {
437b7ce53b6SStefano Zampini   PetscErrorCode ierr;
438b7ce53b6SStefano Zampini 
439b7ce53b6SStefano Zampini   PetscFunctionBegin;
440b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
441b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
442b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
443b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
444b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
445b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
4466c4ed002SBarry Smith     if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix");
447b7ce53b6SStefano Zampini   }
448b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
449b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
450b7ce53b6SStefano Zampini }
451b7ce53b6SStefano Zampini 
452b7ce53b6SStefano Zampini #undef __FUNCT__
453ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
454ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
455ad6194a2SStefano Zampini {
456ad6194a2SStefano Zampini   PetscErrorCode ierr;
457ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
458ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
459ad6194a2SStefano Zampini   Mat            B,localmat;
460ad6194a2SStefano Zampini 
461ad6194a2SStefano Zampini   PetscFunctionBegin;
462ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
463ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
464ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
465e176bc59SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr);
466ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
467ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
468b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
469ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
470ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
471ad6194a2SStefano Zampini   *newmat = B;
472ad6194a2SStefano Zampini   PetscFunctionReturn(0);
473ad6194a2SStefano Zampini }
474ad6194a2SStefano Zampini 
475ad6194a2SStefano Zampini #undef __FUNCT__
47669796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
47769796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
47869796d55SStefano Zampini {
47969796d55SStefano Zampini   PetscErrorCode ierr;
48069796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
48169796d55SStefano Zampini   PetscBool      local_sym;
48269796d55SStefano Zampini 
48369796d55SStefano Zampini   PetscFunctionBegin;
48469796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
485b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
48669796d55SStefano Zampini   PetscFunctionReturn(0);
48769796d55SStefano Zampini }
48869796d55SStefano Zampini 
48969796d55SStefano Zampini #undef __FUNCT__
49069796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
49169796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
49269796d55SStefano Zampini {
49369796d55SStefano Zampini   PetscErrorCode ierr;
49469796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
49569796d55SStefano Zampini   PetscBool      local_sym;
49669796d55SStefano Zampini 
49769796d55SStefano Zampini   PetscFunctionBegin;
49869796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
499b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
50069796d55SStefano Zampini   PetscFunctionReturn(0);
50169796d55SStefano Zampini }
50269796d55SStefano Zampini 
50369796d55SStefano Zampini #undef __FUNCT__
504b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
505dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
506b4319ba4SBarry Smith {
507dfbe8321SBarry Smith   PetscErrorCode ierr;
508b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
509b4319ba4SBarry Smith 
510b4319ba4SBarry Smith   PetscFunctionBegin;
5116bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
512e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr);
513e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr);
5146bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
5156bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
51628f4e0baSStefano Zampini   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
51728f4e0baSStefano Zampini   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
518bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
519dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
520bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
521b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
522b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
5232e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
524b4319ba4SBarry Smith   PetscFunctionReturn(0);
525b4319ba4SBarry Smith }
526b4319ba4SBarry Smith 
527b4319ba4SBarry Smith #undef __FUNCT__
528b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
529dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
530b4319ba4SBarry Smith {
531dfbe8321SBarry Smith   PetscErrorCode ierr;
532b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
533b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
534b4319ba4SBarry Smith 
535b4319ba4SBarry Smith   PetscFunctionBegin;
536b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
537e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
538e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
539b4319ba4SBarry Smith 
540b4319ba4SBarry Smith   /* multiply the local matrix */
541b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
542b4319ba4SBarry Smith 
543b4319ba4SBarry Smith   /* scatter product back into global memory */
5442dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
545e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
546e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
547b4319ba4SBarry Smith   PetscFunctionReturn(0);
548b4319ba4SBarry Smith }
549b4319ba4SBarry Smith 
550b4319ba4SBarry Smith #undef __FUNCT__
5512e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
552b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5532e74eeadSLisandro Dalcin {
554650997f4SStefano Zampini   Vec            temp_vec;
5552e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5562e74eeadSLisandro Dalcin 
5572e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
558650997f4SStefano Zampini   if (v3 != v2) {
559650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
560650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
561650997f4SStefano Zampini   } else {
562650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
563650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
564650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
565650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
566650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
567650997f4SStefano Zampini   }
5682e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5692e74eeadSLisandro Dalcin }
5702e74eeadSLisandro Dalcin 
5712e74eeadSLisandro Dalcin #undef __FUNCT__
5722e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
573e176bc59SStefano Zampini PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x)
5742e74eeadSLisandro Dalcin {
5752e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
5762e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5772e74eeadSLisandro Dalcin 
578e176bc59SStefano Zampini   PetscFunctionBegin;
5792e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
580e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
581e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5822e74eeadSLisandro Dalcin 
5832e74eeadSLisandro Dalcin   /* multiply the local matrix */
584e176bc59SStefano Zampini   ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr);
5852e74eeadSLisandro Dalcin 
5862e74eeadSLisandro Dalcin   /* scatter product back into global vector */
587e176bc59SStefano Zampini   ierr = VecSet(x,0);CHKERRQ(ierr);
588e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
589e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5902e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5912e74eeadSLisandro Dalcin }
5922e74eeadSLisandro Dalcin 
5932e74eeadSLisandro Dalcin #undef __FUNCT__
5942e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
5952e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5962e74eeadSLisandro Dalcin {
597650997f4SStefano Zampini   Vec            temp_vec;
5982e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5992e74eeadSLisandro Dalcin 
6002e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
601650997f4SStefano Zampini   if (v3 != v2) {
602650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
603650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
604650997f4SStefano Zampini   } else {
605650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
606650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
607650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
608650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
609650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
610650997f4SStefano Zampini   }
6112e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
6122e74eeadSLisandro Dalcin }
6132e74eeadSLisandro Dalcin 
6142e74eeadSLisandro Dalcin #undef __FUNCT__
615b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
616dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
617b4319ba4SBarry Smith {
618b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
619dfbe8321SBarry Smith   PetscErrorCode ierr;
620b4319ba4SBarry Smith   PetscViewer    sviewer;
621b4319ba4SBarry Smith 
622b4319ba4SBarry Smith   PetscFunctionBegin;
6233f08860eSBarry Smith   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
624b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
6253f08860eSBarry Smith   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
6266e520ac8SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
627b4319ba4SBarry Smith   PetscFunctionReturn(0);
628b4319ba4SBarry Smith }
629b4319ba4SBarry Smith 
630b4319ba4SBarry Smith #undef __FUNCT__
631b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
632784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
633b4319ba4SBarry Smith {
634dfbe8321SBarry Smith   PetscErrorCode ierr;
635e176bc59SStefano Zampini   PetscInt       nr,rbs,nc,cbs;
636b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
637b4319ba4SBarry Smith   IS             from,to;
638e176bc59SStefano Zampini   Vec            cglobal,rglobal;
639b4319ba4SBarry Smith 
640b4319ba4SBarry Smith   PetscFunctionBegin;
641784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
642e176bc59SStefano Zampini   PetscCheckSameComm(A,1,cmapping,3);
6433bbff08aSStefano Zampini   /* Destroy any previous data */
64470cf5478SStefano Zampini   ierr = VecDestroy(&is->x);CHKERRQ(ierr);
64570cf5478SStefano Zampini   ierr = VecDestroy(&is->y);CHKERRQ(ierr);
646e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr);
647e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr);
6481c47cb0fSStefano Zampini   ierr = MatDestroy(&is->A);CHKERRQ(ierr);
64928f4e0baSStefano Zampini   ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
65028f4e0baSStefano Zampini   ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
6513bbff08aSStefano Zampini 
6523bbff08aSStefano Zampini   /* Setup Layout and set local to global maps */
653fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
654fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
655fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
656fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
657b4319ba4SBarry Smith 
658b4319ba4SBarry Smith   /* Create the local matrix A */
659e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr);
660e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr);
661e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr);
662e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr);
663f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
664e176bc59SStefano Zampini   ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr);
665e176bc59SStefano Zampini   ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr);
666e176bc59SStefano Zampini   ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr);
667ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
668ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
669b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
670b4319ba4SBarry Smith 
671b4319ba4SBarry Smith   /* Create the local work vectors */
6723bbff08aSStefano Zampini   ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr);
673b4319ba4SBarry Smith 
674e176bc59SStefano Zampini   /* setup the global to local scatters */
675e176bc59SStefano Zampini   ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr);
676e176bc59SStefano Zampini   ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr);
677784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
678e176bc59SStefano Zampini   ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr);
679e176bc59SStefano Zampini   if (rmapping != cmapping) {
680e176bc59SStefano Zampini     ierr = ISDestroy(&to);CHKERRQ(ierr);
681e176bc59SStefano Zampini     ierr = ISDestroy(&from);CHKERRQ(ierr);
682e176bc59SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr);
683e176bc59SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr);
684e176bc59SStefano Zampini     ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr);
685e176bc59SStefano Zampini   } else {
686e176bc59SStefano Zampini     ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr);
687e176bc59SStefano Zampini     is->cctx = is->rctx;
688e176bc59SStefano Zampini   }
689e176bc59SStefano Zampini   ierr = VecDestroy(&rglobal);CHKERRQ(ierr);
690e176bc59SStefano Zampini   ierr = VecDestroy(&cglobal);CHKERRQ(ierr);
6916bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
6926bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
69348ff6bf3SStefano Zampini   ierr = MatSetUp(A);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   PetscErrorCode ierr;
70397563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
70497563a80SStefano Zampini   PetscInt       i,zm,zn;
70597563a80SStefano Zampini #endif
70697563a80SStefano Zampini   PetscInt       rows_l[2048],cols_l[2048];
7072e74eeadSLisandro Dalcin 
7082e74eeadSLisandro Dalcin   PetscFunctionBegin;
7092e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
710f23aa3ddSBarry 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);
71197563a80SStefano Zampini   /* count negative indices */
71297563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
71397563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
7142e74eeadSLisandro Dalcin #endif
71597563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
71697563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
71797563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
71897563a80SStefano Zampini   /* count negative indices (should be the same as before) */
71997563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
72097563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
72197563a80SStefano 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");
72297563a80SStefano 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");
72397563a80SStefano Zampini #endif
7242e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
7252e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7262e74eeadSLisandro Dalcin }
7272e74eeadSLisandro Dalcin 
728b4319ba4SBarry Smith #undef __FUNCT__
72997563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS"
73097563a80SStefano Zampini PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
73197563a80SStefano Zampini {
73297563a80SStefano Zampini   Mat_IS         *is = (Mat_IS*)mat->data;
73397563a80SStefano Zampini   PetscErrorCode ierr;
73497563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
73597563a80SStefano Zampini   PetscInt       i,zm,zn;
73697563a80SStefano Zampini #endif
73797563a80SStefano Zampini   PetscInt       rows_l[2048],cols_l[2048];
73897563a80SStefano Zampini 
73997563a80SStefano Zampini   PetscFunctionBegin;
74097563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
74197563a80SStefano 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);
74297563a80SStefano Zampini   /* count negative indices */
74397563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
74497563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
74597563a80SStefano Zampini #endif
74697563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
74797563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
74897563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
74997563a80SStefano Zampini   /* count negative indices (should be the same as before) */
75097563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
75197563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
75297563a80SStefano 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");
75397563a80SStefano 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");
75497563a80SStefano Zampini #endif
75597563a80SStefano Zampini   ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
75697563a80SStefano Zampini   PetscFunctionReturn(0);
75797563a80SStefano Zampini }
75897563a80SStefano Zampini 
75997563a80SStefano Zampini #undef __FUNCT__
760b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
76113f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
762b4319ba4SBarry Smith {
763dfbe8321SBarry Smith   PetscErrorCode ierr;
764b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
765b4319ba4SBarry Smith 
766b4319ba4SBarry Smith   PetscFunctionBegin;
767b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
768b4319ba4SBarry Smith   PetscFunctionReturn(0);
769b4319ba4SBarry Smith }
770b4319ba4SBarry Smith 
771b4319ba4SBarry Smith #undef __FUNCT__
772f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
773f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
774f0006bf2SLisandro Dalcin {
775f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
776f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
777f0006bf2SLisandro Dalcin 
778f0006bf2SLisandro Dalcin   PetscFunctionBegin;
779f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
780f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
781f0006bf2SLisandro Dalcin }
782f0006bf2SLisandro Dalcin 
783f0006bf2SLisandro Dalcin #undef __FUNCT__
7842e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
7852b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
7862e74eeadSLisandro Dalcin {
7876e520ac8SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
7886e520ac8SStefano Zampini   PetscInt       nr,nl,len,i;
7896e520ac8SStefano Zampini   PetscInt       *lrows;
7902e74eeadSLisandro Dalcin   PetscErrorCode ierr;
7912e74eeadSLisandro Dalcin 
7922e74eeadSLisandro Dalcin   PetscFunctionBegin;
7936e520ac8SStefano Zampini   /* get locally owned rows */
7946e520ac8SStefano Zampini   ierr = MatZeroRowsMapLocal_Private(A,n,rows,&len,&lrows);CHKERRQ(ierr);
7956e520ac8SStefano Zampini   /* fix right hand side if needed */
7966e520ac8SStefano Zampini   if (x && b) {
7976e520ac8SStefano Zampini     const PetscScalar *xx;
7986e520ac8SStefano Zampini     PetscScalar       *bb;
7996e520ac8SStefano Zampini 
8006e520ac8SStefano Zampini     ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr);
8016e520ac8SStefano Zampini     ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
8026e520ac8SStefano Zampini     for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]];
8036e520ac8SStefano Zampini     ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr);
8046e520ac8SStefano Zampini     ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr);
8052e74eeadSLisandro Dalcin   }
8066e520ac8SStefano Zampini   /* get rows associated to the local matrices */
8076e520ac8SStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
8086e520ac8SStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
8096e520ac8SStefano Zampini   }
8106e520ac8SStefano Zampini   ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr);
8116e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr);
8126e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
8136e520ac8SStefano Zampini   for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1;
8146e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
8156e520ac8SStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
8166e520ac8SStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
8176e520ac8SStefano Zampini   ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr);
8186e520ac8SStefano Zampini   for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i;
8197230de76SStefano Zampini   ierr = MatISZeroRowsLocal_Private(A,nr,lrows,diag);CHKERRQ(ierr);
8206e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
8212e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
8222e74eeadSLisandro Dalcin }
8232e74eeadSLisandro Dalcin 
8242e74eeadSLisandro Dalcin #undef __FUNCT__
8257230de76SStefano Zampini #define __FUNCT__ "MatISZeroRowsLocal_Private"
8267230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag)
827b4319ba4SBarry Smith {
828b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
829dfbe8321SBarry Smith   PetscErrorCode ierr;
830b4319ba4SBarry Smith 
831b4319ba4SBarry Smith   PetscFunctionBegin;
8326e520ac8SStefano Zampini   if (diag != 0.) {
833b4319ba4SBarry Smith     /*
834b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
835b4319ba4SBarry Smith        work properly in the interface nodes.
836b4319ba4SBarry Smith     */
837b4319ba4SBarry Smith     Vec counter;
838e176bc59SStefano Zampini     ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr);
839e176bc59SStefano Zampini     ierr = VecSet(counter,0.);CHKERRQ(ierr);
840e176bc59SStefano Zampini     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
841e176bc59SStefano Zampini     ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
842e176bc59SStefano Zampini     ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
843e176bc59SStefano Zampini     ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
844e176bc59SStefano Zampini     ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8456bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
846b4319ba4SBarry Smith   }
847958c9bccSBarry Smith   if (!n) {
848b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
849b4319ba4SBarry Smith   } else {
8507230de76SStefano Zampini     PetscInt i;
851b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
8522205254eSKarl Rupp 
8532b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
8546e520ac8SStefano Zampini     if (diag != 0.) {
8556e520ac8SStefano Zampini       const PetscScalar *array;
8566e520ac8SStefano Zampini       ierr = VecGetArrayRead(is->y,&array);CHKERRQ(ierr);
857b4319ba4SBarry Smith       for (i=0; i<n; i++) {
858f4df32b1SMatthew Knepley         ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
859b4319ba4SBarry Smith       }
8606e520ac8SStefano Zampini       ierr = VecRestoreArrayRead(is->y,&array);CHKERRQ(ierr);
8616e520ac8SStefano Zampini     }
862b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
863b4319ba4SBarry Smith     ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
864b4319ba4SBarry Smith   }
865b4319ba4SBarry Smith   PetscFunctionReturn(0);
866b4319ba4SBarry Smith }
867b4319ba4SBarry Smith 
868b4319ba4SBarry Smith #undef __FUNCT__
869b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
870dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
871b4319ba4SBarry Smith {
872b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
873dfbe8321SBarry Smith   PetscErrorCode ierr;
874b4319ba4SBarry Smith 
875b4319ba4SBarry Smith   PetscFunctionBegin;
876b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
877b4319ba4SBarry Smith   PetscFunctionReturn(0);
878b4319ba4SBarry Smith }
879b4319ba4SBarry Smith 
880b4319ba4SBarry Smith #undef __FUNCT__
881b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
882dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
883b4319ba4SBarry Smith {
884b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
885dfbe8321SBarry Smith   PetscErrorCode ierr;
886b4319ba4SBarry Smith 
887b4319ba4SBarry Smith   PetscFunctionBegin;
888b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
889b4319ba4SBarry Smith   PetscFunctionReturn(0);
890b4319ba4SBarry Smith }
891b4319ba4SBarry Smith 
892b4319ba4SBarry Smith #undef __FUNCT__
893b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
8947087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
895b4319ba4SBarry Smith {
896b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
897b4319ba4SBarry Smith 
898b4319ba4SBarry Smith   PetscFunctionBegin;
899b4319ba4SBarry Smith   *local = is->A;
900b4319ba4SBarry Smith   PetscFunctionReturn(0);
901b4319ba4SBarry Smith }
902b4319ba4SBarry Smith 
903b4319ba4SBarry Smith #undef __FUNCT__
904b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
905b4319ba4SBarry Smith /*@
906b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
907b4319ba4SBarry Smith 
908b4319ba4SBarry Smith   Input Parameter:
909b4319ba4SBarry Smith .  mat - the matrix
910b4319ba4SBarry Smith 
911b4319ba4SBarry Smith   Output Parameter:
912eb82efa4SStefano Zampini .  local - the local matrix
913b4319ba4SBarry Smith 
914b4319ba4SBarry Smith   Level: advanced
915b4319ba4SBarry Smith 
916b4319ba4SBarry Smith   Notes:
917b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
918b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
919b4319ba4SBarry Smith   of the MatSetValues() operation.
920b4319ba4SBarry Smith 
921b4319ba4SBarry Smith .seealso: MATIS
922b4319ba4SBarry Smith @*/
9237087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
924b4319ba4SBarry Smith {
9254ac538c5SBarry Smith   PetscErrorCode ierr;
926b4319ba4SBarry Smith 
927b4319ba4SBarry Smith   PetscFunctionBegin;
9280700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
929b4319ba4SBarry Smith   PetscValidPointer(local,2);
9304ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
931b4319ba4SBarry Smith   PetscFunctionReturn(0);
932b4319ba4SBarry Smith }
933b4319ba4SBarry Smith 
9343b03a366Sstefano_zampini #undef __FUNCT__
9353b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
9363b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
9373b03a366Sstefano_zampini {
9383b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
9393b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
9403b03a366Sstefano_zampini   PetscErrorCode ierr;
9413b03a366Sstefano_zampini 
9423b03a366Sstefano_zampini   PetscFunctionBegin;
9434e4c7dbeSStefano Zampini   if (is->A) {
9443b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
9453b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
946f23aa3ddSBarry 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);
9474e4c7dbeSStefano Zampini   }
9483b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
9493b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
9503b03a366Sstefano_zampini   is->A = local;
9513b03a366Sstefano_zampini   PetscFunctionReturn(0);
9523b03a366Sstefano_zampini }
9533b03a366Sstefano_zampini 
9543b03a366Sstefano_zampini #undef __FUNCT__
9553b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
9563b03a366Sstefano_zampini /*@
957eb82efa4SStefano Zampini     MatISSetLocalMat - Replace the local matrix stored inside a MATIS object.
9583b03a366Sstefano_zampini 
9593b03a366Sstefano_zampini   Input Parameter:
9603b03a366Sstefano_zampini .  mat - the matrix
961eb82efa4SStefano Zampini .  local - the local matrix
9623b03a366Sstefano_zampini 
9633b03a366Sstefano_zampini   Output Parameter:
9643b03a366Sstefano_zampini 
9653b03a366Sstefano_zampini   Level: advanced
9663b03a366Sstefano_zampini 
9673b03a366Sstefano_zampini   Notes:
9683b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
9693b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
9703b03a366Sstefano_zampini 
9713b03a366Sstefano_zampini .seealso: MATIS
9723b03a366Sstefano_zampini @*/
9733b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
9743b03a366Sstefano_zampini {
9753b03a366Sstefano_zampini   PetscErrorCode ierr;
9763b03a366Sstefano_zampini 
9773b03a366Sstefano_zampini   PetscFunctionBegin;
9783b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
979b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
9803b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
9813b03a366Sstefano_zampini   PetscFunctionReturn(0);
9823b03a366Sstefano_zampini }
9833b03a366Sstefano_zampini 
9846726f965SBarry Smith #undef __FUNCT__
9856726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
9866726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
9876726f965SBarry Smith {
9886726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
9896726f965SBarry Smith   PetscErrorCode ierr;
9906726f965SBarry Smith 
9916726f965SBarry Smith   PetscFunctionBegin;
9926726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
9936726f965SBarry Smith   PetscFunctionReturn(0);
9946726f965SBarry Smith }
9956726f965SBarry Smith 
9966726f965SBarry Smith #undef __FUNCT__
9972e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
9982e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
9992e74eeadSLisandro Dalcin {
10002e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
10012e74eeadSLisandro Dalcin   PetscErrorCode ierr;
10022e74eeadSLisandro Dalcin 
10032e74eeadSLisandro Dalcin   PetscFunctionBegin;
10042e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
10052e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
10062e74eeadSLisandro Dalcin }
10072e74eeadSLisandro Dalcin 
10082e74eeadSLisandro Dalcin #undef __FUNCT__
10092e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
10102e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
10112e74eeadSLisandro Dalcin {
10122e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
10132e74eeadSLisandro Dalcin   PetscErrorCode ierr;
10142e74eeadSLisandro Dalcin 
10152e74eeadSLisandro Dalcin   PetscFunctionBegin;
10162e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
1017e176bc59SStefano Zampini   ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr);
10182e74eeadSLisandro Dalcin 
10192e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
10202e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
1021e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1022e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10232e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
10242e74eeadSLisandro Dalcin }
10252e74eeadSLisandro Dalcin 
10262e74eeadSLisandro Dalcin #undef __FUNCT__
10276726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
1028ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
10296726f965SBarry Smith {
10306726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
10316726f965SBarry Smith   PetscErrorCode ierr;
10326726f965SBarry Smith 
10336726f965SBarry Smith   PetscFunctionBegin;
10344e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
10356726f965SBarry Smith   PetscFunctionReturn(0);
10366726f965SBarry Smith }
10376726f965SBarry Smith 
1038284134d9SBarry Smith #undef __FUNCT__
1039284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
1040284134d9SBarry Smith /*@
10413c212e90SHong Zhang     MatCreateIS - Creates a "process" unassembled matrix, assembled on each
1042284134d9SBarry Smith        process but not across processes.
1043284134d9SBarry Smith 
1044284134d9SBarry Smith    Input Parameters:
1045284134d9SBarry Smith +     comm    - MPI communicator that will share the matrix
1046e176bc59SStefano Zampini .     bs      - block size of the matrix
1047df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
1048e176bc59SStefano Zampini .     rmap    - local to global map for rows
1049e176bc59SStefano Zampini -     cmap    - local to global map for cols
1050284134d9SBarry Smith 
1051284134d9SBarry Smith    Output Parameter:
1052284134d9SBarry Smith .    A - the resulting matrix
1053284134d9SBarry Smith 
10548e6c10adSSatish Balay    Level: advanced
10558e6c10adSSatish Balay 
10563c212e90SHong Zhang    Notes: See MATIS for more details.
10573c212e90SHong Zhang           m and n are NOT related to the size of the map; they are the size of the part of the vector owned
1058e176bc59SStefano Zampini           by that process. The sizes of rmap and cmap define the size of the local matrices.
10593c212e90SHong Zhang           If either rmap or cmap are NULL, then the matrix is assumed to be square.
1060284134d9SBarry Smith 
1061284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
1062284134d9SBarry Smith @*/
1063e176bc59SStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A)
1064284134d9SBarry Smith {
1065284134d9SBarry Smith   PetscErrorCode ierr;
1066284134d9SBarry Smith 
1067284134d9SBarry Smith   PetscFunctionBegin;
1068e176bc59SStefano Zampini   if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping");
1069284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1070d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
1071284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1072284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
1073e176bc59SStefano Zampini   if (rmap && cmap) {
1074e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr);
1075e176bc59SStefano Zampini   } else if (!rmap) {
1076e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr);
1077e176bc59SStefano Zampini   } else {
1078e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr);
1079e176bc59SStefano Zampini   }
1080284134d9SBarry Smith   PetscFunctionReturn(0);
1081284134d9SBarry Smith }
1082284134d9SBarry Smith 
1083b4319ba4SBarry Smith /*MC
1084eb82efa4SStefano Zampini    MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC).
1085b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
1086b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
1087b4319ba4SBarry Smith    product is handled "implicitly".
1088b4319ba4SBarry Smith 
1089b4319ba4SBarry Smith    Operations Provided:
10906726f965SBarry Smith +  MatMult()
10912e74eeadSLisandro Dalcin .  MatMultAdd()
10922e74eeadSLisandro Dalcin .  MatMultTranspose()
10932e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
10946726f965SBarry Smith .  MatZeroEntries()
10956726f965SBarry Smith .  MatSetOption()
10962e74eeadSLisandro Dalcin .  MatZeroRows()
10972e74eeadSLisandro Dalcin .  MatSetValues()
109897563a80SStefano Zampini .  MatSetValuesBlocked()
10996726f965SBarry Smith .  MatSetValuesLocal()
110097563a80SStefano Zampini .  MatSetValuesBlockedLocal()
11012e74eeadSLisandro Dalcin .  MatScale()
11022e74eeadSLisandro Dalcin .  MatGetDiagonal()
11036726f965SBarry Smith -  MatSetLocalToGlobalMapping()
1104b4319ba4SBarry Smith 
1105b4319ba4SBarry Smith    Options Database Keys:
1106b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
1107b4319ba4SBarry Smith 
1108b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
1109b4319ba4SBarry Smith 
1110b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
1111b4319ba4SBarry Smith 
1112b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
1113eb82efa4SStefano Zampini           MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation()
1114b4319ba4SBarry Smith 
1115b4319ba4SBarry Smith   Level: advanced
1116b4319ba4SBarry Smith 
11173c212e90SHong Zhang .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC, MatCreateIS()
1118b4319ba4SBarry Smith 
1119b4319ba4SBarry Smith M*/
1120b4319ba4SBarry Smith 
1121b4319ba4SBarry Smith #undef __FUNCT__
1122b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
11238cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
1124b4319ba4SBarry Smith {
1125dfbe8321SBarry Smith   PetscErrorCode ierr;
1126b4319ba4SBarry Smith   Mat_IS         *b;
1127b4319ba4SBarry Smith 
1128b4319ba4SBarry Smith   PetscFunctionBegin;
1129b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1130b4319ba4SBarry Smith   A->data = (void*)b;
1131b4319ba4SBarry Smith 
1132e176bc59SStefano Zampini   /* matrix ops */
1133e176bc59SStefano Zampini   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1134b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
11352e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
11362e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
11372e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1138b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1139b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
11402e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
114197563a80SStefano Zampini   A->ops->setvalues               = MatSetValuesBlocked_IS;
1142b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1143f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
11442e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1145b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1146b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1147b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
11486726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
11492e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
11502e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
11516726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
115269796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
115369796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1154ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
1155*6bd84002SStefano Zampini   A->ops->missingdiagonal         = MatMissingDiagonal_IS;
1156b4319ba4SBarry Smith 
1157b7ce53b6SStefano Zampini   /* special MATIS functions */
1158bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1159bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1160b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
11612e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
116217667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1163b4319ba4SBarry Smith   PetscFunctionReturn(0);
1164b4319ba4SBarry Smith }
1165