xref: /petsc/src/mat/impls/is/matis.c (revision 2b4041121a3593a003f52cc04b9a09ee62be1886)
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*2b404112SStefano Zampini #define __FUNCT__ "MatCopy_IS"
18*2b404112SStefano Zampini PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str)
19*2b404112SStefano Zampini {
20*2b404112SStefano Zampini   Mat_IS         *a = (Mat_IS*)A->data,*b;
21*2b404112SStefano Zampini   PetscBool      ismatis;
22*2b404112SStefano Zampini   PetscErrorCode ierr;
23*2b404112SStefano Zampini 
24*2b404112SStefano Zampini   PetscFunctionBegin;
25*2b404112SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr);
26*2b404112SStefano Zampini   if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented");
27*2b404112SStefano Zampini   b = (Mat_IS*)B->data;
28*2b404112SStefano Zampini   ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr);
29*2b404112SStefano Zampini   PetscFunctionReturn(0);
30*2b404112SStefano Zampini }
31*2b404112SStefano Zampini 
32*2b404112SStefano Zampini #undef __FUNCT__
336bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS"
346bd84002SStefano Zampini PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool  *missing,PetscInt *d)
356bd84002SStefano Zampini {
366bd84002SStefano Zampini   Mat_IS         *a = (Mat_IS*)A->data;
376bd84002SStefano Zampini   PetscErrorCode ierr;
386bd84002SStefano Zampini 
396bd84002SStefano Zampini   PetscFunctionBegin;
406bd84002SStefano Zampini   if (d) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Need to be implemented");
416bd84002SStefano Zampini   ierr = MatMissingDiagonal(a->A,missing,NULL);CHKERRQ(ierr);
426bd84002SStefano Zampini   PetscFunctionReturn(0);
436bd84002SStefano Zampini }
446bd84002SStefano Zampini 
456bd84002SStefano Zampini #undef __FUNCT__
4628f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private"
47a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B)
4828f4e0baSStefano Zampini {
4928f4e0baSStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
5028f4e0baSStefano Zampini   const PetscInt *gidxs;
5128f4e0baSStefano Zampini   PetscErrorCode ierr;
5228f4e0baSStefano Zampini 
5328f4e0baSStefano Zampini   PetscFunctionBegin;
5428f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr);
5528f4e0baSStefano Zampini   ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr);
5628f4e0baSStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr);
573bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
5828f4e0baSStefano Zampini   /* PETSC_OWN_POINTER refers to ilocal which is NULL */
5928f4e0baSStefano Zampini   ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
603bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
6128f4e0baSStefano Zampini   ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr);
6228f4e0baSStefano Zampini   PetscFunctionReturn(0);
6328f4e0baSStefano Zampini }
642e1947a5SStefano Zampini 
652e1947a5SStefano Zampini #undef __FUNCT__
662e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation"
67eb82efa4SStefano Zampini /*@
68a88811baSStefano Zampini    MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix.
69a88811baSStefano Zampini 
70a88811baSStefano Zampini    Collective on MPI_Comm
71a88811baSStefano Zampini 
72a88811baSStefano Zampini    Input Parameters:
73a88811baSStefano Zampini +  B - the matrix
74a88811baSStefano Zampini .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
75a88811baSStefano Zampini            (same value is used for all local rows)
76a88811baSStefano Zampini .  d_nnz - array containing the number of nonzeros in the various rows of the
77a88811baSStefano Zampini            DIAGONAL portion of the local submatrix (possibly different for each row)
78a88811baSStefano Zampini            or NULL, if d_nz is used to specify the nonzero structure.
79a88811baSStefano Zampini            The size of this array is equal to the number of local rows, i.e 'm'.
80a88811baSStefano Zampini            For matrices that will be factored, you must leave room for (and set)
81a88811baSStefano Zampini            the diagonal entry even if it is zero.
82a88811baSStefano Zampini .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
83a88811baSStefano Zampini            submatrix (same value is used for all local rows).
84a88811baSStefano Zampini -  o_nnz - array containing the number of nonzeros in the various rows of the
85a88811baSStefano Zampini            OFF-DIAGONAL portion of the local submatrix (possibly different for
86a88811baSStefano Zampini            each row) or NULL, if o_nz is used to specify the nonzero
87a88811baSStefano Zampini            structure. The size of this array is equal to the number
88a88811baSStefano Zampini            of local rows, i.e 'm'.
89a88811baSStefano Zampini 
90a88811baSStefano Zampini    If the *_nnz parameter is given then the *_nz parameter is ignored
91a88811baSStefano Zampini 
92a88811baSStefano Zampini    Level: intermediate
93a88811baSStefano Zampini 
94a88811baSStefano Zampini    Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition
95a88811baSStefano Zampini           from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local
96a88811baSStefano Zampini           matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects.
97a88811baSStefano Zampini 
98a88811baSStefano Zampini .keywords: matrix
99a88811baSStefano Zampini 
1003c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS
101a88811baSStefano Zampini @*/
1022e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
1032e1947a5SStefano Zampini {
1042e1947a5SStefano Zampini   PetscErrorCode ierr;
1052e1947a5SStefano Zampini 
1062e1947a5SStefano Zampini   PetscFunctionBegin;
1072e1947a5SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
1082e1947a5SStefano Zampini   PetscValidType(B,1);
1092e1947a5SStefano Zampini   ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
1102e1947a5SStefano Zampini   PetscFunctionReturn(0);
1112e1947a5SStefano Zampini }
1122e1947a5SStefano Zampini 
1132e1947a5SStefano Zampini #undef __FUNCT__
1142e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS"
1157230de76SStefano Zampini static PetscErrorCode  MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
1162e1947a5SStefano Zampini {
1172e1947a5SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
11828f4e0baSStefano Zampini   PetscInt       bs,i,nlocalcols;
1192e1947a5SStefano Zampini   PetscErrorCode ierr;
1202e1947a5SStefano Zampini 
1212e1947a5SStefano Zampini   PetscFunctionBegin;
1226c4ed002SBarry Smith   if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
12328f4e0baSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
12428f4e0baSStefano Zampini     ierr = MatISComputeSF_Private(B);CHKERRQ(ierr);
12528f4e0baSStefano Zampini   }
1262e1947a5SStefano Zampini   if (!d_nnz) {
12728f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz;
1282e1947a5SStefano Zampini   } else {
12928f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i];
1302e1947a5SStefano Zampini   }
1312e1947a5SStefano Zampini   if (!o_nnz) {
13228f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz;
1332e1947a5SStefano Zampini   } else {
13428f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i];
1352e1947a5SStefano Zampini   }
13628f4e0baSStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
13728f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr);
13828f4e0baSStefano Zampini   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
13928f4e0baSStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
14028f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves;i++) {
14128f4e0baSStefano Zampini     matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols);
1422e1947a5SStefano Zampini   }
14328f4e0baSStefano Zampini   ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr);
14428f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
14528f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs;
1462e1947a5SStefano Zampini   }
14728f4e0baSStefano Zampini   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
14828f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
14928f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i]-i;
1502e1947a5SStefano Zampini   }
15128f4e0baSStefano Zampini   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
1522e1947a5SStefano Zampini   PetscFunctionReturn(0);
1532e1947a5SStefano Zampini }
154b4319ba4SBarry Smith 
155b4319ba4SBarry Smith #undef __FUNCT__
1563927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private"
1573927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
1583927de2eSStefano Zampini {
1593927de2eSStefano Zampini   Mat_IS          *matis = (Mat_IS*)(A->data);
1603927de2eSStefano Zampini   PetscInt        *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership;
161ecf5a873SStefano Zampini   const PetscInt  *global_indices_r,*global_indices_c;
1623927de2eSStefano Zampini   PetscInt        i,j,bs,rows,cols;
1633927de2eSStefano Zampini   PetscInt        lrows,lcols;
1643927de2eSStefano Zampini   PetscInt        local_rows,local_cols;
1653927de2eSStefano Zampini   PetscMPIInt     nsubdomains;
1663927de2eSStefano Zampini   PetscBool       isdense,issbaij;
1673927de2eSStefano Zampini   PetscErrorCode  ierr;
1683927de2eSStefano Zampini 
1693927de2eSStefano Zampini   PetscFunctionBegin;
1703927de2eSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr);
1713927de2eSStefano Zampini   ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr);
1723927de2eSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1733927de2eSStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
1743927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
1753927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
176ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
177ecf5a873SStefano Zampini   if (A->rmap->mapping != A->cmap->mapping) {
1787230de76SStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
179ecf5a873SStefano Zampini   } else {
180ecf5a873SStefano Zampini     global_indices_c = global_indices_r;
181ecf5a873SStefano Zampini   }
182ecf5a873SStefano Zampini 
1833927de2eSStefano Zampini   if (issbaij) {
1843927de2eSStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
1853927de2eSStefano Zampini   }
1863927de2eSStefano Zampini   /*
187ecf5a873SStefano Zampini      An SF reduce is needed to sum up properly on shared rows.
1883927de2eSStefano Zampini      Note that generally preallocation is not exact, since it overestimates nonzeros
1893927de2eSStefano Zampini   */
1903927de2eSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
1913927de2eSStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
1923927de2eSStefano Zampini   }
1933927de2eSStefano Zampini   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
1943927de2eSStefano Zampini   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
1953927de2eSStefano Zampini   /* All processes need to compute entire row ownership */
1963927de2eSStefano Zampini   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
1973927de2eSStefano Zampini   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
1983927de2eSStefano Zampini   for (i=0;i<nsubdomains;i++) {
1993927de2eSStefano Zampini     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
2003927de2eSStefano Zampini       row_ownership[j] = i;
2013927de2eSStefano Zampini     }
2023927de2eSStefano Zampini   }
2037230de76SStefano Zampini   ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
2043927de2eSStefano Zampini 
2053927de2eSStefano Zampini   /*
2063927de2eSStefano Zampini      my_dnz and my_onz contains exact contribution to preallocation from each local mat
2073927de2eSStefano Zampini      then, they will be summed up properly. This way, preallocation is always sufficient
2083927de2eSStefano Zampini   */
2093927de2eSStefano Zampini   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
2103927de2eSStefano Zampini   /* preallocation as a MATAIJ */
2113927de2eSStefano Zampini   if (isdense) { /* special case for dense local matrices */
2123927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
213ecf5a873SStefano Zampini       PetscInt index_row = global_indices_r[i];
2143927de2eSStefano Zampini       for (j=i;j<local_rows;j++) {
2153927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
216ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[j];
2173927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
2183927de2eSStefano Zampini           my_dnz[i] += 1;
2193927de2eSStefano Zampini         } else { /* offdiag block */
2203927de2eSStefano Zampini           my_onz[i] += 1;
2213927de2eSStefano Zampini         }
2223927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
2233927de2eSStefano Zampini         if (i != j) {
2243927de2eSStefano Zampini           owner = row_ownership[index_col];
2253927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
2263927de2eSStefano Zampini             my_dnz[j] += 1;
2273927de2eSStefano Zampini           } else {
2283927de2eSStefano Zampini             my_onz[j] += 1;
2293927de2eSStefano Zampini           }
2303927de2eSStefano Zampini         }
2313927de2eSStefano Zampini       }
2323927de2eSStefano Zampini     }
233ecf5a873SStefano Zampini   } else { /* TODO: this could be optimized using MatGetRowIJ */
2343927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
2353927de2eSStefano Zampini       const PetscInt *cols;
236ecf5a873SStefano Zampini       PetscInt       ncols,index_row = global_indices_r[i];
2373927de2eSStefano Zampini       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2383927de2eSStefano Zampini       for (j=0;j<ncols;j++) {
2393927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
240ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[cols[j]];
2413927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
2423927de2eSStefano Zampini           my_dnz[i] += 1;
2433927de2eSStefano Zampini         } else { /* offdiag block */
2443927de2eSStefano Zampini           my_onz[i] += 1;
2453927de2eSStefano Zampini         }
2463927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
247d9a9e74cSStefano Zampini         if (issbaij && index_col != index_row) {
2483927de2eSStefano Zampini           owner = row_ownership[index_col];
2493927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
250d9a9e74cSStefano Zampini             my_dnz[cols[j]] += 1;
2513927de2eSStefano Zampini           } else {
252d9a9e74cSStefano Zampini             my_onz[cols[j]] += 1;
2533927de2eSStefano Zampini           }
2543927de2eSStefano Zampini         }
2553927de2eSStefano Zampini       }
2563927de2eSStefano Zampini       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2573927de2eSStefano Zampini     }
2583927de2eSStefano Zampini   }
259ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
260ecf5a873SStefano Zampini   if (global_indices_c != global_indices_r) {
2617230de76SStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
262ecf5a873SStefano Zampini   }
2633927de2eSStefano Zampini   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
264ecf5a873SStefano Zampini 
265ecf5a873SStefano Zampini   /* Reduce my_dnz and my_onz */
2663927de2eSStefano Zampini   if (maxreduce) {
2673927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2683927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2693927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2703927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2713927de2eSStefano Zampini   } else {
2723927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2733927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2743927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2753927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2763927de2eSStefano Zampini   }
2773927de2eSStefano Zampini   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
2783927de2eSStefano Zampini 
2793927de2eSStefano Zampini   /* Resize preallocation if overestimated */
2803927de2eSStefano Zampini   for (i=0;i<lrows;i++) {
2813927de2eSStefano Zampini     dnz[i] = PetscMin(dnz[i],lcols);
2823927de2eSStefano Zampini     onz[i] = PetscMin(onz[i],cols-lcols);
2833927de2eSStefano Zampini   }
2843927de2eSStefano Zampini   /* set preallocation */
2853927de2eSStefano Zampini   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
2863927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
2873927de2eSStefano Zampini     dnz[i] = dnz[i*bs]/bs;
2883927de2eSStefano Zampini     onz[i] = onz[i*bs]/bs;
2893927de2eSStefano Zampini   }
2903927de2eSStefano Zampini   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2913927de2eSStefano Zampini   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2923927de2eSStefano Zampini   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
2933927de2eSStefano Zampini   if (issbaij) {
2943927de2eSStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
2953927de2eSStefano Zampini   }
2963927de2eSStefano Zampini   PetscFunctionReturn(0);
2973927de2eSStefano Zampini }
2983927de2eSStefano Zampini 
2993927de2eSStefano Zampini #undef __FUNCT__
300b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
3017230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
302b7ce53b6SStefano Zampini {
303b7ce53b6SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
304d9a9e74cSStefano Zampini   Mat            local_mat;
305b7ce53b6SStefano Zampini   /* info on mat */
3063cfa4ea4SStefano Zampini   PetscInt       bs,rows,cols,lrows,lcols;
307b7ce53b6SStefano Zampini   PetscInt       local_rows,local_cols;
308686e3a49SStefano Zampini   PetscBool      isdense,issbaij,isseqaij;
3097c03b4e8SStefano Zampini   PetscMPIInt    nsubdomains;
310b7ce53b6SStefano Zampini   /* values insertion */
311b7ce53b6SStefano Zampini   PetscScalar    *array;
312b7ce53b6SStefano Zampini   /* work */
313b7ce53b6SStefano Zampini   PetscErrorCode ierr;
314b7ce53b6SStefano Zampini 
315b7ce53b6SStefano Zampini   PetscFunctionBegin;
316b7ce53b6SStefano Zampini   /* get info from mat */
3177c03b4e8SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
3187c03b4e8SStefano Zampini   if (nsubdomains == 1) {
3197c03b4e8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
3207c03b4e8SStefano Zampini       ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr);
3217c03b4e8SStefano Zampini     } else {
3227c03b4e8SStefano Zampini       ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
3237c03b4e8SStefano Zampini     }
3247c03b4e8SStefano Zampini     PetscFunctionReturn(0);
3257c03b4e8SStefano Zampini   }
326b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
327b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
3283cfa4ea4SStefano Zampini   ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr);
329b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
330b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
331686e3a49SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
332b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
333b7ce53b6SStefano Zampini 
334b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
335b7ce53b6SStefano Zampini     MatType     new_mat_type;
3363927de2eSStefano Zampini     PetscBool   issbaij_red;
337b7ce53b6SStefano Zampini 
338b7ce53b6SStefano Zampini     /* determining new matrix type */
339b2566f29SBarry Smith     ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
340b7ce53b6SStefano Zampini     if (issbaij_red) {
341b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
342b7ce53b6SStefano Zampini     } else {
343b7ce53b6SStefano Zampini       if (bs>1) {
344b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
345b7ce53b6SStefano Zampini       } else {
346b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
347b7ce53b6SStefano Zampini       }
348b7ce53b6SStefano Zampini     }
349b7ce53b6SStefano Zampini 
3503927de2eSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
3513cfa4ea4SStefano Zampini     ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr);
3523927de2eSStefano Zampini     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
3533927de2eSStefano Zampini     ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr);
3543927de2eSStefano Zampini     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
355b7ce53b6SStefano Zampini   } else {
3563cfa4ea4SStefano Zampini     PetscInt mbs,mrows,mcols,mlrows,mlcols;
357b7ce53b6SStefano Zampini     /* some checks */
358b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
359b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
3603cfa4ea4SStefano Zampini     ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr);
3616c4ed002SBarry Smith     if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
3626c4ed002SBarry Smith     if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
3636c4ed002SBarry Smith     if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows);
3646c4ed002SBarry Smith     if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols);
3656c4ed002SBarry Smith     if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
366b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
367b7ce53b6SStefano Zampini   }
368d9a9e74cSStefano Zampini 
369d9a9e74cSStefano Zampini   if (issbaij) {
370d9a9e74cSStefano Zampini     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr);
371d9a9e74cSStefano Zampini   } else {
372d9a9e74cSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
373d9a9e74cSStefano Zampini     local_mat = matis->A;
374d9a9e74cSStefano Zampini   }
375686e3a49SStefano Zampini 
376b7ce53b6SStefano Zampini   /* Set values */
377ecf5a873SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr);
378b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
379ecf5a873SStefano Zampini     PetscInt i,*dummy_rows,*dummy_cols;
380ecf5a873SStefano Zampini 
381ecf5a873SStefano Zampini     if (local_rows != local_cols) {
382ecf5a873SStefano Zampini       ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr);
383ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
384ecf5a873SStefano Zampini       for (i=0;i<local_cols;i++) dummy_cols[i] = i;
385ecf5a873SStefano Zampini     } else {
386ecf5a873SStefano Zampini       ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr);
387ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
388ecf5a873SStefano Zampini       dummy_cols = dummy_rows;
389ecf5a873SStefano Zampini     }
390b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
391d9a9e74cSStefano Zampini     ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr);
392ecf5a873SStefano Zampini     ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr);
393d9a9e74cSStefano Zampini     ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr);
394ecf5a873SStefano Zampini     if (dummy_rows != dummy_cols) {
395ecf5a873SStefano Zampini       ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr);
396ecf5a873SStefano Zampini     } else {
397ecf5a873SStefano Zampini       ierr = PetscFree(dummy_rows);CHKERRQ(ierr);
398ecf5a873SStefano Zampini     }
399686e3a49SStefano Zampini   } else if (isseqaij) {
400ecf5a873SStefano Zampini     PetscInt  i,nvtxs,*xadj,*adjncy;
401686e3a49SStefano Zampini     PetscBool done;
402686e3a49SStefano Zampini 
403d9a9e74cSStefano Zampini     ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
4046c4ed002SBarry Smith     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
405d9a9e74cSStefano Zampini     ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr);
406686e3a49SStefano Zampini     for (i=0;i<nvtxs;i++) {
407ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr);
408686e3a49SStefano Zampini     }
409d9a9e74cSStefano Zampini     ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
4106c4ed002SBarry Smith     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
411d9a9e74cSStefano Zampini     ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr);
412686e3a49SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
413ecf5a873SStefano Zampini     PetscInt i;
414c0962df8SStefano Zampini 
415686e3a49SStefano Zampini     for (i=0;i<local_rows;i++) {
416686e3a49SStefano Zampini       PetscInt       j;
417ecf5a873SStefano Zampini       const PetscInt *local_indices_cols;
418686e3a49SStefano Zampini 
419ecf5a873SStefano Zampini       ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
420ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr);
421ecf5a873SStefano Zampini       ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
422686e3a49SStefano Zampini     }
423b7ce53b6SStefano Zampini   }
424b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
425d9a9e74cSStefano Zampini   ierr = MatDestroy(&local_mat);CHKERRQ(ierr);
426b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
427b7ce53b6SStefano Zampini   if (isdense) {
428b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
429b7ce53b6SStefano Zampini   }
430b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
431b7ce53b6SStefano Zampini }
432b7ce53b6SStefano Zampini 
433b7ce53b6SStefano Zampini #undef __FUNCT__
434b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
435b7ce53b6SStefano Zampini /*@
436b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
437b7ce53b6SStefano Zampini 
438b7ce53b6SStefano Zampini   Input Parameter:
439b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
440b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
441b7ce53b6SStefano Zampini 
442b7ce53b6SStefano Zampini   Output Parameter:
443b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
444b7ce53b6SStefano Zampini 
445b7ce53b6SStefano Zampini   Level: developer
446b7ce53b6SStefano Zampini 
447eb82efa4SStefano Zampini   Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested.
448b7ce53b6SStefano Zampini 
449b7ce53b6SStefano Zampini .seealso: MATIS
450b7ce53b6SStefano Zampini @*/
451b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
452b7ce53b6SStefano Zampini {
453b7ce53b6SStefano Zampini   PetscErrorCode ierr;
454b7ce53b6SStefano Zampini 
455b7ce53b6SStefano Zampini   PetscFunctionBegin;
456b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
457b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
458b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
459b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
460b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
461b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
4626c4ed002SBarry Smith     if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix");
463b7ce53b6SStefano Zampini   }
464b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
465b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
466b7ce53b6SStefano Zampini }
467b7ce53b6SStefano Zampini 
468b7ce53b6SStefano Zampini #undef __FUNCT__
469ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
470ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
471ad6194a2SStefano Zampini {
472ad6194a2SStefano Zampini   PetscErrorCode ierr;
473ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
474ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
475ad6194a2SStefano Zampini   Mat            B,localmat;
476ad6194a2SStefano Zampini 
477ad6194a2SStefano Zampini   PetscFunctionBegin;
478ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
479ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
480ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
481e176bc59SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr);
482ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
483ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
484b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
485ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
486ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
487ad6194a2SStefano Zampini   *newmat = B;
488ad6194a2SStefano Zampini   PetscFunctionReturn(0);
489ad6194a2SStefano Zampini }
490ad6194a2SStefano Zampini 
491ad6194a2SStefano Zampini #undef __FUNCT__
49269796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
49369796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
49469796d55SStefano Zampini {
49569796d55SStefano Zampini   PetscErrorCode ierr;
49669796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
49769796d55SStefano Zampini   PetscBool      local_sym;
49869796d55SStefano Zampini 
49969796d55SStefano Zampini   PetscFunctionBegin;
50069796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
501b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
50269796d55SStefano Zampini   PetscFunctionReturn(0);
50369796d55SStefano Zampini }
50469796d55SStefano Zampini 
50569796d55SStefano Zampini #undef __FUNCT__
50669796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
50769796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
50869796d55SStefano Zampini {
50969796d55SStefano Zampini   PetscErrorCode ierr;
51069796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
51169796d55SStefano Zampini   PetscBool      local_sym;
51269796d55SStefano Zampini 
51369796d55SStefano Zampini   PetscFunctionBegin;
51469796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
515b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
51669796d55SStefano Zampini   PetscFunctionReturn(0);
51769796d55SStefano Zampini }
51869796d55SStefano Zampini 
51969796d55SStefano Zampini #undef __FUNCT__
520b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
521dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
522b4319ba4SBarry Smith {
523dfbe8321SBarry Smith   PetscErrorCode ierr;
524b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
525b4319ba4SBarry Smith 
526b4319ba4SBarry Smith   PetscFunctionBegin;
5276bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
528e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr);
529e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr);
5306bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
5316bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
53228f4e0baSStefano Zampini   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
53328f4e0baSStefano Zampini   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
534bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
535dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
536bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
537b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
538b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
5392e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
540b4319ba4SBarry Smith   PetscFunctionReturn(0);
541b4319ba4SBarry Smith }
542b4319ba4SBarry Smith 
543b4319ba4SBarry Smith #undef __FUNCT__
544b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
545dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
546b4319ba4SBarry Smith {
547dfbe8321SBarry Smith   PetscErrorCode ierr;
548b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
549b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
550b4319ba4SBarry Smith 
551b4319ba4SBarry Smith   PetscFunctionBegin;
552b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
553e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
554e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
555b4319ba4SBarry Smith 
556b4319ba4SBarry Smith   /* multiply the local matrix */
557b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
558b4319ba4SBarry Smith 
559b4319ba4SBarry Smith   /* scatter product back into global memory */
5602dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
561e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
562e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
563b4319ba4SBarry Smith   PetscFunctionReturn(0);
564b4319ba4SBarry Smith }
565b4319ba4SBarry Smith 
566b4319ba4SBarry Smith #undef __FUNCT__
5672e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
568b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5692e74eeadSLisandro Dalcin {
570650997f4SStefano Zampini   Vec            temp_vec;
5712e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5722e74eeadSLisandro Dalcin 
5732e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
574650997f4SStefano Zampini   if (v3 != v2) {
575650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
576650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
577650997f4SStefano Zampini   } else {
578650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
579650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
580650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
581650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
582650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
583650997f4SStefano Zampini   }
5842e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5852e74eeadSLisandro Dalcin }
5862e74eeadSLisandro Dalcin 
5872e74eeadSLisandro Dalcin #undef __FUNCT__
5882e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
589e176bc59SStefano Zampini PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x)
5902e74eeadSLisandro Dalcin {
5912e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
5922e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5932e74eeadSLisandro Dalcin 
594e176bc59SStefano Zampini   PetscFunctionBegin;
5952e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
596e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
597e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5982e74eeadSLisandro Dalcin 
5992e74eeadSLisandro Dalcin   /* multiply the local matrix */
600e176bc59SStefano Zampini   ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr);
6012e74eeadSLisandro Dalcin 
6022e74eeadSLisandro Dalcin   /* scatter product back into global vector */
603e176bc59SStefano Zampini   ierr = VecSet(x,0);CHKERRQ(ierr);
604e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
605e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6062e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
6072e74eeadSLisandro Dalcin }
6082e74eeadSLisandro Dalcin 
6092e74eeadSLisandro Dalcin #undef __FUNCT__
6102e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
6112e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
6122e74eeadSLisandro Dalcin {
613650997f4SStefano Zampini   Vec            temp_vec;
6142e74eeadSLisandro Dalcin   PetscErrorCode ierr;
6152e74eeadSLisandro Dalcin 
6162e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
617650997f4SStefano Zampini   if (v3 != v2) {
618650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
619650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
620650997f4SStefano Zampini   } else {
621650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
622650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
623650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
624650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
625650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
626650997f4SStefano Zampini   }
6272e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
6282e74eeadSLisandro Dalcin }
6292e74eeadSLisandro Dalcin 
6302e74eeadSLisandro Dalcin #undef __FUNCT__
631b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
632dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
633b4319ba4SBarry Smith {
634b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
635dfbe8321SBarry Smith   PetscErrorCode ierr;
636b4319ba4SBarry Smith   PetscViewer    sviewer;
637b4319ba4SBarry Smith 
638b4319ba4SBarry Smith   PetscFunctionBegin;
6393f08860eSBarry Smith   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
640b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
6413f08860eSBarry Smith   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
6426e520ac8SStefano Zampini   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
643b4319ba4SBarry Smith   PetscFunctionReturn(0);
644b4319ba4SBarry Smith }
645b4319ba4SBarry Smith 
646b4319ba4SBarry Smith #undef __FUNCT__
647b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
648784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
649b4319ba4SBarry Smith {
650dfbe8321SBarry Smith   PetscErrorCode ierr;
651e176bc59SStefano Zampini   PetscInt       nr,rbs,nc,cbs;
652b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
653b4319ba4SBarry Smith   IS             from,to;
654e176bc59SStefano Zampini   Vec            cglobal,rglobal;
655b4319ba4SBarry Smith 
656b4319ba4SBarry Smith   PetscFunctionBegin;
657784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
658e176bc59SStefano Zampini   PetscCheckSameComm(A,1,cmapping,3);
6593bbff08aSStefano Zampini   /* Destroy any previous data */
66070cf5478SStefano Zampini   ierr = VecDestroy(&is->x);CHKERRQ(ierr);
66170cf5478SStefano Zampini   ierr = VecDestroy(&is->y);CHKERRQ(ierr);
662e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr);
663e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr);
6641c47cb0fSStefano Zampini   ierr = MatDestroy(&is->A);CHKERRQ(ierr);
66528f4e0baSStefano Zampini   ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
66628f4e0baSStefano Zampini   ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
6673bbff08aSStefano Zampini 
6683bbff08aSStefano Zampini   /* Setup Layout and set local to global maps */
669fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
670fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
671fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
672fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
673b4319ba4SBarry Smith 
674b4319ba4SBarry Smith   /* Create the local matrix A */
675e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr);
676e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr);
677e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr);
678e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr);
679f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
680e176bc59SStefano Zampini   ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr);
681e176bc59SStefano Zampini   ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr);
682e176bc59SStefano Zampini   ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr);
683ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
684ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
685b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
686b4319ba4SBarry Smith 
687b4319ba4SBarry Smith   /* Create the local work vectors */
6883bbff08aSStefano Zampini   ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr);
689b4319ba4SBarry Smith 
690e176bc59SStefano Zampini   /* setup the global to local scatters */
691e176bc59SStefano Zampini   ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr);
692e176bc59SStefano Zampini   ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr);
693784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
694e176bc59SStefano Zampini   ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr);
695e176bc59SStefano Zampini   if (rmapping != cmapping) {
696e176bc59SStefano Zampini     ierr = ISDestroy(&to);CHKERRQ(ierr);
697e176bc59SStefano Zampini     ierr = ISDestroy(&from);CHKERRQ(ierr);
698e176bc59SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr);
699e176bc59SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr);
700e176bc59SStefano Zampini     ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr);
701e176bc59SStefano Zampini   } else {
702e176bc59SStefano Zampini     ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr);
703e176bc59SStefano Zampini     is->cctx = is->rctx;
704e176bc59SStefano Zampini   }
705e176bc59SStefano Zampini   ierr = VecDestroy(&rglobal);CHKERRQ(ierr);
706e176bc59SStefano Zampini   ierr = VecDestroy(&cglobal);CHKERRQ(ierr);
7076bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
7086bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
70948ff6bf3SStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
710b4319ba4SBarry Smith   PetscFunctionReturn(0);
711b4319ba4SBarry Smith }
712b4319ba4SBarry Smith 
7132e74eeadSLisandro Dalcin #undef __FUNCT__
7142e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
7152e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
7162e74eeadSLisandro Dalcin {
7172e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
7182e74eeadSLisandro Dalcin   PetscErrorCode ierr;
71997563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
72097563a80SStefano Zampini   PetscInt       i,zm,zn;
72197563a80SStefano Zampini #endif
72297563a80SStefano Zampini   PetscInt       rows_l[2048],cols_l[2048];
7232e74eeadSLisandro Dalcin 
7242e74eeadSLisandro Dalcin   PetscFunctionBegin;
7252e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
726f23aa3ddSBarry 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);
72797563a80SStefano Zampini   /* count negative indices */
72897563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
72997563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
7302e74eeadSLisandro Dalcin #endif
73197563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
73297563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
73397563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
73497563a80SStefano Zampini   /* count negative indices (should be the same as before) */
73597563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
73697563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
73797563a80SStefano 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");
73897563a80SStefano 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");
73997563a80SStefano Zampini #endif
7402e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
7412e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7422e74eeadSLisandro Dalcin }
7432e74eeadSLisandro Dalcin 
744b4319ba4SBarry Smith #undef __FUNCT__
74597563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS"
74697563a80SStefano Zampini PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
74797563a80SStefano Zampini {
74897563a80SStefano Zampini   Mat_IS         *is = (Mat_IS*)mat->data;
74997563a80SStefano Zampini   PetscErrorCode ierr;
75097563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
75197563a80SStefano Zampini   PetscInt       i,zm,zn;
75297563a80SStefano Zampini #endif
75397563a80SStefano Zampini   PetscInt       rows_l[2048],cols_l[2048];
75497563a80SStefano Zampini 
75597563a80SStefano Zampini   PetscFunctionBegin;
75697563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
75797563a80SStefano 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);
75897563a80SStefano Zampini   /* count negative indices */
75997563a80SStefano Zampini   for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++;
76097563a80SStefano Zampini   for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++;
76197563a80SStefano Zampini #endif
76297563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr);
76397563a80SStefano Zampini   ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr);
76497563a80SStefano Zampini #if defined(PETSC_USE_DEBUG)
76597563a80SStefano Zampini   /* count negative indices (should be the same as before) */
76697563a80SStefano Zampini   for (i=0;i<m;i++) if (rows_l[i] < 0) zm--;
76797563a80SStefano Zampini   for (i=0;i<n;i++) if (cols_l[i] < 0) zn--;
76897563a80SStefano 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");
76997563a80SStefano 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");
77097563a80SStefano Zampini #endif
77197563a80SStefano Zampini   ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
77297563a80SStefano Zampini   PetscFunctionReturn(0);
77397563a80SStefano Zampini }
77497563a80SStefano Zampini 
77597563a80SStefano Zampini #undef __FUNCT__
776b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
77713f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
778b4319ba4SBarry Smith {
779dfbe8321SBarry Smith   PetscErrorCode ierr;
780b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
781b4319ba4SBarry Smith 
782b4319ba4SBarry Smith   PetscFunctionBegin;
783b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
784b4319ba4SBarry Smith   PetscFunctionReturn(0);
785b4319ba4SBarry Smith }
786b4319ba4SBarry Smith 
787b4319ba4SBarry Smith #undef __FUNCT__
788f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
789f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
790f0006bf2SLisandro Dalcin {
791f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
792f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
793f0006bf2SLisandro Dalcin 
794f0006bf2SLisandro Dalcin   PetscFunctionBegin;
795f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
796f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
797f0006bf2SLisandro Dalcin }
798f0006bf2SLisandro Dalcin 
799f0006bf2SLisandro Dalcin #undef __FUNCT__
8002e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
8012b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
8022e74eeadSLisandro Dalcin {
8036e520ac8SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
8046e520ac8SStefano Zampini   PetscInt       nr,nl,len,i;
8056e520ac8SStefano Zampini   PetscInt       *lrows;
8062e74eeadSLisandro Dalcin   PetscErrorCode ierr;
8072e74eeadSLisandro Dalcin 
8082e74eeadSLisandro Dalcin   PetscFunctionBegin;
8096e520ac8SStefano Zampini   /* get locally owned rows */
8106e520ac8SStefano Zampini   ierr = MatZeroRowsMapLocal_Private(A,n,rows,&len,&lrows);CHKERRQ(ierr);
8116e520ac8SStefano Zampini   /* fix right hand side if needed */
8126e520ac8SStefano Zampini   if (x && b) {
8136e520ac8SStefano Zampini     const PetscScalar *xx;
8146e520ac8SStefano Zampini     PetscScalar       *bb;
8156e520ac8SStefano Zampini 
8166e520ac8SStefano Zampini     ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr);
8176e520ac8SStefano Zampini     ierr = VecGetArray(b, &bb);CHKERRQ(ierr);
8186e520ac8SStefano Zampini     for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]];
8196e520ac8SStefano Zampini     ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr);
8206e520ac8SStefano Zampini     ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr);
8212e74eeadSLisandro Dalcin   }
8226e520ac8SStefano Zampini   /* get rows associated to the local matrices */
8236e520ac8SStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
8246e520ac8SStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
8256e520ac8SStefano Zampini   }
8266e520ac8SStefano Zampini   ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr);
8276e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr);
8286e520ac8SStefano Zampini   ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr);
8296e520ac8SStefano Zampini   for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1;
8306e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
8316e520ac8SStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
8326e520ac8SStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
8336e520ac8SStefano Zampini   ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr);
8346e520ac8SStefano Zampini   for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i;
8357230de76SStefano Zampini   ierr = MatISZeroRowsLocal_Private(A,nr,lrows,diag);CHKERRQ(ierr);
8366e520ac8SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
8372e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
8382e74eeadSLisandro Dalcin }
8392e74eeadSLisandro Dalcin 
8402e74eeadSLisandro Dalcin #undef __FUNCT__
8417230de76SStefano Zampini #define __FUNCT__ "MatISZeroRowsLocal_Private"
8427230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag)
843b4319ba4SBarry Smith {
844b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
845dfbe8321SBarry Smith   PetscErrorCode ierr;
846b4319ba4SBarry Smith 
847b4319ba4SBarry Smith   PetscFunctionBegin;
8486e520ac8SStefano Zampini   if (diag != 0.) {
849b4319ba4SBarry Smith     /*
850b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
851b4319ba4SBarry Smith        work properly in the interface nodes.
852b4319ba4SBarry Smith     */
853b4319ba4SBarry Smith     Vec counter;
854e176bc59SStefano Zampini     ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr);
855e176bc59SStefano Zampini     ierr = VecSet(counter,0.);CHKERRQ(ierr);
856e176bc59SStefano Zampini     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
857e176bc59SStefano Zampini     ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
858e176bc59SStefano Zampini     ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
859e176bc59SStefano Zampini     ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
860e176bc59SStefano Zampini     ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
8616bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
862b4319ba4SBarry Smith   }
863*2b404112SStefano Zampini 
864958c9bccSBarry Smith   if (!n) {
865b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
866b4319ba4SBarry Smith   } else {
8677230de76SStefano Zampini     PetscInt i;
868b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
8692205254eSKarl Rupp 
8702b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
8716e520ac8SStefano Zampini     if (diag != 0.) {
8726e520ac8SStefano Zampini       const PetscScalar *array;
8736e520ac8SStefano Zampini       ierr = VecGetArrayRead(is->y,&array);CHKERRQ(ierr);
874b4319ba4SBarry Smith       for (i=0; i<n; i++) {
875f4df32b1SMatthew Knepley         ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
876b4319ba4SBarry Smith       }
8776e520ac8SStefano Zampini       ierr = VecRestoreArrayRead(is->y,&array);CHKERRQ(ierr);
8786e520ac8SStefano Zampini     }
879b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
880b4319ba4SBarry Smith     ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
881b4319ba4SBarry Smith   }
882b4319ba4SBarry Smith   PetscFunctionReturn(0);
883b4319ba4SBarry Smith }
884b4319ba4SBarry Smith 
885b4319ba4SBarry Smith #undef __FUNCT__
886b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
887dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
888b4319ba4SBarry Smith {
889b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
890dfbe8321SBarry Smith   PetscErrorCode ierr;
891b4319ba4SBarry Smith 
892b4319ba4SBarry Smith   PetscFunctionBegin;
893b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
894b4319ba4SBarry Smith   PetscFunctionReturn(0);
895b4319ba4SBarry Smith }
896b4319ba4SBarry Smith 
897b4319ba4SBarry Smith #undef __FUNCT__
898b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
899dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
900b4319ba4SBarry Smith {
901b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
902dfbe8321SBarry Smith   PetscErrorCode ierr;
903b4319ba4SBarry Smith 
904b4319ba4SBarry Smith   PetscFunctionBegin;
905b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
906b4319ba4SBarry Smith   PetscFunctionReturn(0);
907b4319ba4SBarry Smith }
908b4319ba4SBarry Smith 
909b4319ba4SBarry Smith #undef __FUNCT__
910b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
9117087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
912b4319ba4SBarry Smith {
913b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
914b4319ba4SBarry Smith 
915b4319ba4SBarry Smith   PetscFunctionBegin;
916b4319ba4SBarry Smith   *local = is->A;
917b4319ba4SBarry Smith   PetscFunctionReturn(0);
918b4319ba4SBarry Smith }
919b4319ba4SBarry Smith 
920b4319ba4SBarry Smith #undef __FUNCT__
921b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
922b4319ba4SBarry Smith /*@
923b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
924b4319ba4SBarry Smith 
925b4319ba4SBarry Smith   Input Parameter:
926b4319ba4SBarry Smith .  mat - the matrix
927b4319ba4SBarry Smith 
928b4319ba4SBarry Smith   Output Parameter:
929eb82efa4SStefano Zampini .  local - the local matrix
930b4319ba4SBarry Smith 
931b4319ba4SBarry Smith   Level: advanced
932b4319ba4SBarry Smith 
933b4319ba4SBarry Smith   Notes:
934b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
935b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
936b4319ba4SBarry Smith   of the MatSetValues() operation.
937b4319ba4SBarry Smith 
938b4319ba4SBarry Smith .seealso: MATIS
939b4319ba4SBarry Smith @*/
9407087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
941b4319ba4SBarry Smith {
9424ac538c5SBarry Smith   PetscErrorCode ierr;
943b4319ba4SBarry Smith 
944b4319ba4SBarry Smith   PetscFunctionBegin;
9450700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
946b4319ba4SBarry Smith   PetscValidPointer(local,2);
9474ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
948b4319ba4SBarry Smith   PetscFunctionReturn(0);
949b4319ba4SBarry Smith }
950b4319ba4SBarry Smith 
9513b03a366Sstefano_zampini #undef __FUNCT__
9523b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
9533b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
9543b03a366Sstefano_zampini {
9553b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
9563b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
9573b03a366Sstefano_zampini   PetscErrorCode ierr;
9583b03a366Sstefano_zampini 
9593b03a366Sstefano_zampini   PetscFunctionBegin;
9604e4c7dbeSStefano Zampini   if (is->A) {
9613b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
9623b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
963f23aa3ddSBarry 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);
9644e4c7dbeSStefano Zampini   }
9653b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
9663b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
9673b03a366Sstefano_zampini   is->A = local;
9683b03a366Sstefano_zampini   PetscFunctionReturn(0);
9693b03a366Sstefano_zampini }
9703b03a366Sstefano_zampini 
9713b03a366Sstefano_zampini #undef __FUNCT__
9723b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
9733b03a366Sstefano_zampini /*@
974eb82efa4SStefano Zampini     MatISSetLocalMat - Replace the local matrix stored inside a MATIS object.
9753b03a366Sstefano_zampini 
9763b03a366Sstefano_zampini   Input Parameter:
9773b03a366Sstefano_zampini .  mat - the matrix
978eb82efa4SStefano Zampini .  local - the local matrix
9793b03a366Sstefano_zampini 
9803b03a366Sstefano_zampini   Output Parameter:
9813b03a366Sstefano_zampini 
9823b03a366Sstefano_zampini   Level: advanced
9833b03a366Sstefano_zampini 
9843b03a366Sstefano_zampini   Notes:
9853b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
9863b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
9873b03a366Sstefano_zampini 
9883b03a366Sstefano_zampini .seealso: MATIS
9893b03a366Sstefano_zampini @*/
9903b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
9913b03a366Sstefano_zampini {
9923b03a366Sstefano_zampini   PetscErrorCode ierr;
9933b03a366Sstefano_zampini 
9943b03a366Sstefano_zampini   PetscFunctionBegin;
9953b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
996b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
9973b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
9983b03a366Sstefano_zampini   PetscFunctionReturn(0);
9993b03a366Sstefano_zampini }
10003b03a366Sstefano_zampini 
10016726f965SBarry Smith #undef __FUNCT__
10026726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
10036726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
10046726f965SBarry Smith {
10056726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
10066726f965SBarry Smith   PetscErrorCode ierr;
10076726f965SBarry Smith 
10086726f965SBarry Smith   PetscFunctionBegin;
10096726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
10106726f965SBarry Smith   PetscFunctionReturn(0);
10116726f965SBarry Smith }
10126726f965SBarry Smith 
10136726f965SBarry Smith #undef __FUNCT__
10142e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
10152e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
10162e74eeadSLisandro Dalcin {
10172e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
10182e74eeadSLisandro Dalcin   PetscErrorCode ierr;
10192e74eeadSLisandro Dalcin 
10202e74eeadSLisandro Dalcin   PetscFunctionBegin;
10212e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
10222e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
10232e74eeadSLisandro Dalcin }
10242e74eeadSLisandro Dalcin 
10252e74eeadSLisandro Dalcin #undef __FUNCT__
10262e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
10272e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
10282e74eeadSLisandro Dalcin {
10292e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
10302e74eeadSLisandro Dalcin   PetscErrorCode ierr;
10312e74eeadSLisandro Dalcin 
10322e74eeadSLisandro Dalcin   PetscFunctionBegin;
10332e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
1034e176bc59SStefano Zampini   ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr);
10352e74eeadSLisandro Dalcin 
10362e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
10372e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
1038e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
1039e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
10402e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
10412e74eeadSLisandro Dalcin }
10422e74eeadSLisandro Dalcin 
10432e74eeadSLisandro Dalcin #undef __FUNCT__
10446726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
1045ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
10466726f965SBarry Smith {
10476726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
10486726f965SBarry Smith   PetscErrorCode ierr;
10496726f965SBarry Smith 
10506726f965SBarry Smith   PetscFunctionBegin;
10514e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
10526726f965SBarry Smith   PetscFunctionReturn(0);
10536726f965SBarry Smith }
10546726f965SBarry Smith 
1055284134d9SBarry Smith #undef __FUNCT__
1056284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
1057284134d9SBarry Smith /*@
10583c212e90SHong Zhang     MatCreateIS - Creates a "process" unassembled matrix, assembled on each
1059284134d9SBarry Smith        process but not across processes.
1060284134d9SBarry Smith 
1061284134d9SBarry Smith    Input Parameters:
1062284134d9SBarry Smith +     comm    - MPI communicator that will share the matrix
1063e176bc59SStefano Zampini .     bs      - block size of the matrix
1064df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
1065e176bc59SStefano Zampini .     rmap    - local to global map for rows
1066e176bc59SStefano Zampini -     cmap    - local to global map for cols
1067284134d9SBarry Smith 
1068284134d9SBarry Smith    Output Parameter:
1069284134d9SBarry Smith .    A - the resulting matrix
1070284134d9SBarry Smith 
10718e6c10adSSatish Balay    Level: advanced
10728e6c10adSSatish Balay 
10733c212e90SHong Zhang    Notes: See MATIS for more details.
10743c212e90SHong Zhang           m and n are NOT related to the size of the map; they are the size of the part of the vector owned
1075e176bc59SStefano Zampini           by that process. The sizes of rmap and cmap define the size of the local matrices.
10763c212e90SHong Zhang           If either rmap or cmap are NULL, then the matrix is assumed to be square.
1077284134d9SBarry Smith 
1078284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
1079284134d9SBarry Smith @*/
1080e176bc59SStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A)
1081284134d9SBarry Smith {
1082284134d9SBarry Smith   PetscErrorCode ierr;
1083284134d9SBarry Smith 
1084284134d9SBarry Smith   PetscFunctionBegin;
1085e176bc59SStefano Zampini   if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping");
1086284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1087d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
1088284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1089284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
1090e176bc59SStefano Zampini   if (rmap && cmap) {
1091e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr);
1092e176bc59SStefano Zampini   } else if (!rmap) {
1093e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr);
1094e176bc59SStefano Zampini   } else {
1095e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr);
1096e176bc59SStefano Zampini   }
1097284134d9SBarry Smith   PetscFunctionReturn(0);
1098284134d9SBarry Smith }
1099284134d9SBarry Smith 
1100b4319ba4SBarry Smith /*MC
1101eb82efa4SStefano Zampini    MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC).
1102b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
1103b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
1104b4319ba4SBarry Smith    product is handled "implicitly".
1105b4319ba4SBarry Smith 
1106b4319ba4SBarry Smith    Operations Provided:
11076726f965SBarry Smith +  MatMult()
11082e74eeadSLisandro Dalcin .  MatMultAdd()
11092e74eeadSLisandro Dalcin .  MatMultTranspose()
11102e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
11116726f965SBarry Smith .  MatZeroEntries()
11126726f965SBarry Smith .  MatSetOption()
11132e74eeadSLisandro Dalcin .  MatZeroRows()
11142e74eeadSLisandro Dalcin .  MatSetValues()
111597563a80SStefano Zampini .  MatSetValuesBlocked()
11166726f965SBarry Smith .  MatSetValuesLocal()
111797563a80SStefano Zampini .  MatSetValuesBlockedLocal()
11182e74eeadSLisandro Dalcin .  MatScale()
11192e74eeadSLisandro Dalcin .  MatGetDiagonal()
1120*2b404112SStefano Zampini .  MatMissingDiagonal()
1121*2b404112SStefano Zampini .  MatDuplicate()
1122*2b404112SStefano Zampini .  MatCopy()
11236726f965SBarry Smith -  MatSetLocalToGlobalMapping()
1124b4319ba4SBarry Smith 
1125b4319ba4SBarry Smith    Options Database Keys:
1126b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
1127b4319ba4SBarry Smith 
1128b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
1129b4319ba4SBarry Smith 
1130b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
1131b4319ba4SBarry Smith 
1132b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
1133eb82efa4SStefano Zampini           MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation()
1134b4319ba4SBarry Smith 
1135b4319ba4SBarry Smith   Level: advanced
1136b4319ba4SBarry Smith 
11373c212e90SHong Zhang .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC, MatCreateIS()
1138b4319ba4SBarry Smith 
1139b4319ba4SBarry Smith M*/
1140b4319ba4SBarry Smith 
1141b4319ba4SBarry Smith #undef __FUNCT__
1142b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
11438cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
1144b4319ba4SBarry Smith {
1145dfbe8321SBarry Smith   PetscErrorCode ierr;
1146b4319ba4SBarry Smith   Mat_IS         *b;
1147b4319ba4SBarry Smith 
1148b4319ba4SBarry Smith   PetscFunctionBegin;
1149b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1150b4319ba4SBarry Smith   A->data = (void*)b;
1151b4319ba4SBarry Smith 
1152e176bc59SStefano Zampini   /* matrix ops */
1153e176bc59SStefano Zampini   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1154b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
11552e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
11562e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
11572e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1158b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1159b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
11602e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
116197563a80SStefano Zampini   A->ops->setvalues               = MatSetValuesBlocked_IS;
1162b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1163f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
11642e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1165b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1166b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1167b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
11686726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
11692e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
11702e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
11716726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
117269796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
117369796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1174ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
11756bd84002SStefano Zampini   A->ops->missingdiagonal         = MatMissingDiagonal_IS;
1176*2b404112SStefano Zampini   A->ops->copy                    = MatCopy_IS;
1177b4319ba4SBarry Smith 
1178b7ce53b6SStefano Zampini   /* special MATIS functions */
1179bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1180bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1181b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
11822e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
118317667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1184b4319ba4SBarry Smith   PetscFunctionReturn(0);
1185b4319ba4SBarry Smith }
1186