xref: /petsc/src/mat/impls/is/matis.c (revision 0ea065fb06d751599c4157d36bfe1a1b41348e0b)
1be1d678aSKris Buschelman 
2b4319ba4SBarry Smith /*
3b4319ba4SBarry Smith     Creates a matrix class for using the Neumann-Neumann type preconditioners.
4b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
5b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
6b4319ba4SBarry Smith    product is handled "implicitly".
7b4319ba4SBarry Smith 
8b4319ba4SBarry Smith      We provide:
9b4319ba4SBarry Smith          MatMult()
10b4319ba4SBarry Smith 
11b4319ba4SBarry Smith     Currently this allows for only one subdomain per processor.
12b4319ba4SBarry Smith 
13b4319ba4SBarry Smith */
14b4319ba4SBarry Smith 
15c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h>      /*I "petscmat.h" I*/
1628f4e0baSStefano Zampini 
17a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat);
18a72627d2SStefano Zampini 
1928f4e0baSStefano Zampini #undef __FUNCT__
2028f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private"
21a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B)
2228f4e0baSStefano Zampini {
2328f4e0baSStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
2428f4e0baSStefano Zampini   const PetscInt *gidxs;
2528f4e0baSStefano Zampini   PetscErrorCode ierr;
2628f4e0baSStefano Zampini 
2728f4e0baSStefano Zampini   PetscFunctionBegin;
2828f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr);
2928f4e0baSStefano Zampini   ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr);
3028f4e0baSStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr);
313bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
3228f4e0baSStefano Zampini   /* PETSC_OWN_POINTER refers to ilocal which is NULL */
3328f4e0baSStefano Zampini   ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
343bbff08aSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr);
3528f4e0baSStefano Zampini   ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr);
3628f4e0baSStefano Zampini   PetscFunctionReturn(0);
3728f4e0baSStefano Zampini }
382e1947a5SStefano Zampini 
392e1947a5SStefano Zampini #undef __FUNCT__
402e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation"
41eb82efa4SStefano Zampini /*@
42a88811baSStefano Zampini    MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix.
43a88811baSStefano Zampini 
44a88811baSStefano Zampini    Collective on MPI_Comm
45a88811baSStefano Zampini 
46a88811baSStefano Zampini    Input Parameters:
47a88811baSStefano Zampini +  B - the matrix
48a88811baSStefano Zampini .  d_nz  - number of nonzeros per row in DIAGONAL portion of local submatrix
49a88811baSStefano Zampini            (same value is used for all local rows)
50a88811baSStefano Zampini .  d_nnz - array containing the number of nonzeros in the various rows of the
51a88811baSStefano Zampini            DIAGONAL portion of the local submatrix (possibly different for each row)
52a88811baSStefano Zampini            or NULL, if d_nz is used to specify the nonzero structure.
53a88811baSStefano Zampini            The size of this array is equal to the number of local rows, i.e 'm'.
54a88811baSStefano Zampini            For matrices that will be factored, you must leave room for (and set)
55a88811baSStefano Zampini            the diagonal entry even if it is zero.
56a88811baSStefano Zampini .  o_nz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
57a88811baSStefano Zampini            submatrix (same value is used for all local rows).
58a88811baSStefano Zampini -  o_nnz - array containing the number of nonzeros in the various rows of the
59a88811baSStefano Zampini            OFF-DIAGONAL portion of the local submatrix (possibly different for
60a88811baSStefano Zampini            each row) or NULL, if o_nz is used to specify the nonzero
61a88811baSStefano Zampini            structure. The size of this array is equal to the number
62a88811baSStefano Zampini            of local rows, i.e 'm'.
63a88811baSStefano Zampini 
64a88811baSStefano Zampini    If the *_nnz parameter is given then the *_nz parameter is ignored
65a88811baSStefano Zampini 
66a88811baSStefano Zampini    Level: intermediate
67a88811baSStefano Zampini 
68a88811baSStefano Zampini    Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition
69a88811baSStefano Zampini           from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local
70a88811baSStefano Zampini           matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects.
71a88811baSStefano Zampini 
72a88811baSStefano Zampini .keywords: matrix
73a88811baSStefano Zampini 
74a88811baSStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat()
75a88811baSStefano Zampini @*/
762e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
772e1947a5SStefano Zampini {
782e1947a5SStefano Zampini   PetscErrorCode ierr;
792e1947a5SStefano Zampini 
802e1947a5SStefano Zampini   PetscFunctionBegin;
812e1947a5SStefano Zampini   PetscValidHeaderSpecific(B,MAT_CLASSID,1);
822e1947a5SStefano Zampini   PetscValidType(B,1);
832e1947a5SStefano Zampini   ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr);
842e1947a5SStefano Zampini   PetscFunctionReturn(0);
852e1947a5SStefano Zampini }
862e1947a5SStefano Zampini 
872e1947a5SStefano Zampini #undef __FUNCT__
882e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS"
892e1947a5SStefano Zampini PetscErrorCode  MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
902e1947a5SStefano Zampini {
912e1947a5SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
9228f4e0baSStefano Zampini   PetscInt       bs,i,nlocalcols;
932e1947a5SStefano Zampini   PetscErrorCode ierr;
942e1947a5SStefano Zampini 
952e1947a5SStefano Zampini   PetscFunctionBegin;
966c4ed002SBarry Smith   if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
9728f4e0baSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
9828f4e0baSStefano Zampini     ierr = MatISComputeSF_Private(B);CHKERRQ(ierr);
9928f4e0baSStefano Zampini   }
1002e1947a5SStefano Zampini   if (!d_nnz) {
10128f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz;
1022e1947a5SStefano Zampini   } else {
10328f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i];
1042e1947a5SStefano Zampini   }
1052e1947a5SStefano Zampini   if (!o_nnz) {
10628f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz;
1072e1947a5SStefano Zampini   } else {
10828f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i];
1092e1947a5SStefano Zampini   }
11028f4e0baSStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
11128f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr);
11228f4e0baSStefano Zampini   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
11328f4e0baSStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
11428f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves;i++) {
11528f4e0baSStefano Zampini     matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols);
1162e1947a5SStefano Zampini   }
11728f4e0baSStefano Zampini   ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr);
11828f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
11928f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs;
1202e1947a5SStefano Zampini   }
12128f4e0baSStefano Zampini   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
12228f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
12328f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i]-i;
1242e1947a5SStefano Zampini   }
12528f4e0baSStefano Zampini   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
1262e1947a5SStefano Zampini   PetscFunctionReturn(0);
1272e1947a5SStefano Zampini }
128b4319ba4SBarry Smith 
129b4319ba4SBarry Smith #undef __FUNCT__
1303927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private"
1313927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
1323927de2eSStefano Zampini {
1333927de2eSStefano Zampini   Mat_IS          *matis = (Mat_IS*)(A->data);
1343927de2eSStefano Zampini   PetscInt        *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership;
135ecf5a873SStefano Zampini   const PetscInt  *global_indices_r,*global_indices_c;
1363927de2eSStefano Zampini   PetscInt        i,j,bs,rows,cols;
1373927de2eSStefano Zampini   PetscInt        lrows,lcols;
1383927de2eSStefano Zampini   PetscInt        local_rows,local_cols;
1393927de2eSStefano Zampini   PetscMPIInt     nsubdomains;
1403927de2eSStefano Zampini   PetscBool       isdense,issbaij;
1413927de2eSStefano Zampini   PetscErrorCode  ierr;
1423927de2eSStefano Zampini 
1433927de2eSStefano Zampini   PetscFunctionBegin;
1443927de2eSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr);
1453927de2eSStefano Zampini   ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr);
1463927de2eSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1473927de2eSStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
1483927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
1493927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
150ecf5a873SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
151ecf5a873SStefano Zampini   if (A->rmap->mapping != A->cmap->mapping) {
152*0ea065fbSStefano Zampini     ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
153ecf5a873SStefano Zampini   } else {
154ecf5a873SStefano Zampini     global_indices_c = global_indices_r;
155ecf5a873SStefano Zampini   }
156ecf5a873SStefano Zampini 
1573927de2eSStefano Zampini   if (issbaij) {
1583927de2eSStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
1593927de2eSStefano Zampini   }
1603927de2eSStefano Zampini   /*
161ecf5a873SStefano Zampini      An SF reduce is needed to sum up properly on shared rows.
1623927de2eSStefano Zampini      Note that generally preallocation is not exact, since it overestimates nonzeros
1633927de2eSStefano Zampini   */
1643927de2eSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
1653927de2eSStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
1663927de2eSStefano Zampini   }
1673927de2eSStefano Zampini   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
1683927de2eSStefano Zampini   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
1693927de2eSStefano Zampini   /* All processes need to compute entire row ownership */
1703927de2eSStefano Zampini   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
1713927de2eSStefano Zampini   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
1723927de2eSStefano Zampini   for (i=0;i<nsubdomains;i++) {
1733927de2eSStefano Zampini     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
1743927de2eSStefano Zampini       row_ownership[j] = i;
1753927de2eSStefano Zampini     }
1763927de2eSStefano Zampini   }
177*0ea065fbSStefano Zampini   ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
1783927de2eSStefano Zampini 
1793927de2eSStefano Zampini   /*
1803927de2eSStefano Zampini      my_dnz and my_onz contains exact contribution to preallocation from each local mat
1813927de2eSStefano Zampini      then, they will be summed up properly. This way, preallocation is always sufficient
1823927de2eSStefano Zampini   */
1833927de2eSStefano Zampini   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
1843927de2eSStefano Zampini   /* preallocation as a MATAIJ */
1853927de2eSStefano Zampini   if (isdense) { /* special case for dense local matrices */
1863927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
187ecf5a873SStefano Zampini       PetscInt index_row = global_indices_r[i];
1884c8dd594SStefano Zampini       for (j=i;j<local_cols;j++) {
1893927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
190ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[j];
1913927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
1923927de2eSStefano Zampini           my_dnz[i] += 1;
1933927de2eSStefano Zampini         } else { /* offdiag block */
1943927de2eSStefano Zampini           my_onz[i] += 1;
1953927de2eSStefano Zampini         }
1963927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
1973927de2eSStefano Zampini         if (i != j) {
1983927de2eSStefano Zampini           owner = row_ownership[index_col];
1993927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
2003927de2eSStefano Zampini             my_dnz[j] += 1;
2013927de2eSStefano Zampini           } else {
2023927de2eSStefano Zampini             my_onz[j] += 1;
2033927de2eSStefano Zampini           }
2043927de2eSStefano Zampini         }
2053927de2eSStefano Zampini       }
2063927de2eSStefano Zampini     }
207ecf5a873SStefano Zampini   } else { /* TODO: this could be optimized using MatGetRowIJ */
2083927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
2093927de2eSStefano Zampini       const PetscInt *cols;
210ecf5a873SStefano Zampini       PetscInt       ncols,index_row = global_indices_r[i];
2113927de2eSStefano Zampini       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2123927de2eSStefano Zampini       for (j=0;j<ncols;j++) {
2133927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
214ecf5a873SStefano Zampini         PetscInt index_col = global_indices_c[cols[j]];
2153927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
2163927de2eSStefano Zampini           my_dnz[i] += 1;
2173927de2eSStefano Zampini         } else { /* offdiag block */
2183927de2eSStefano Zampini           my_onz[i] += 1;
2193927de2eSStefano Zampini         }
2203927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
221d9a9e74cSStefano Zampini         if (issbaij && index_col != index_row) {
2223927de2eSStefano Zampini           owner = row_ownership[index_col];
2233927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
224d9a9e74cSStefano Zampini             my_dnz[cols[j]] += 1;
2253927de2eSStefano Zampini           } else {
226d9a9e74cSStefano Zampini             my_onz[cols[j]] += 1;
2273927de2eSStefano Zampini           }
2283927de2eSStefano Zampini         }
2293927de2eSStefano Zampini       }
2303927de2eSStefano Zampini       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2313927de2eSStefano Zampini     }
2323927de2eSStefano Zampini   }
233ecf5a873SStefano Zampini   if (global_indices_c != global_indices_r) {
234*0ea065fbSStefano Zampini     ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr);
235ecf5a873SStefano Zampini   }
2364f619741Sstefano_zampini   ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr);
2373927de2eSStefano Zampini   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
238ecf5a873SStefano Zampini 
239ecf5a873SStefano Zampini   /* Reduce my_dnz and my_onz */
2403927de2eSStefano Zampini   if (maxreduce) {
2413927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2423927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2433927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2443927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2453927de2eSStefano Zampini   } else {
2463927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2473927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2483927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2493927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2503927de2eSStefano Zampini   }
2513927de2eSStefano Zampini   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
2523927de2eSStefano Zampini 
2533927de2eSStefano Zampini   /* Resize preallocation if overestimated */
2543927de2eSStefano Zampini   for (i=0;i<lrows;i++) {
2553927de2eSStefano Zampini     dnz[i] = PetscMin(dnz[i],lcols);
2563927de2eSStefano Zampini     onz[i] = PetscMin(onz[i],cols-lcols);
2573927de2eSStefano Zampini   }
2583927de2eSStefano Zampini   /* set preallocation */
2593927de2eSStefano Zampini   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
2603927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
2613927de2eSStefano Zampini     dnz[i] = dnz[i*bs]/bs;
2623927de2eSStefano Zampini     onz[i] = onz[i*bs]/bs;
2633927de2eSStefano Zampini   }
2643927de2eSStefano Zampini   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2653927de2eSStefano Zampini   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2663927de2eSStefano Zampini   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
2673927de2eSStefano Zampini   if (issbaij) {
2683927de2eSStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
2693927de2eSStefano Zampini   }
2703927de2eSStefano Zampini   PetscFunctionReturn(0);
2713927de2eSStefano Zampini }
2723927de2eSStefano Zampini 
2733927de2eSStefano Zampini #undef __FUNCT__
274b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
275b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
276b7ce53b6SStefano Zampini {
277b7ce53b6SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
278d9a9e74cSStefano Zampini   Mat            local_mat;
279b7ce53b6SStefano Zampini   /* info on mat */
2803cfa4ea4SStefano Zampini   PetscInt       bs,rows,cols,lrows,lcols;
281b7ce53b6SStefano Zampini   PetscInt       local_rows,local_cols;
282686e3a49SStefano Zampini   PetscBool      isdense,issbaij,isseqaij;
2837c03b4e8SStefano Zampini   PetscMPIInt    nsubdomains;
284b7ce53b6SStefano Zampini   /* values insertion */
285b7ce53b6SStefano Zampini   PetscScalar    *array;
286b7ce53b6SStefano Zampini   /* work */
287b7ce53b6SStefano Zampini   PetscErrorCode ierr;
288b7ce53b6SStefano Zampini 
289b7ce53b6SStefano Zampini   PetscFunctionBegin;
290b7ce53b6SStefano Zampini   /* get info from mat */
2917c03b4e8SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
2927c03b4e8SStefano Zampini   if (nsubdomains == 1) {
2937c03b4e8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
2947c03b4e8SStefano Zampini       ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr);
2957c03b4e8SStefano Zampini     } else {
2967c03b4e8SStefano Zampini       ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
2977c03b4e8SStefano Zampini     }
2987c03b4e8SStefano Zampini     PetscFunctionReturn(0);
2997c03b4e8SStefano Zampini   }
300b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
301b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
3023cfa4ea4SStefano Zampini   ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr);
303b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
304b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
305686e3a49SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
306b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
307b7ce53b6SStefano Zampini 
308b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
309b7ce53b6SStefano Zampini     MatType     new_mat_type;
3103927de2eSStefano Zampini     PetscBool   issbaij_red;
311b7ce53b6SStefano Zampini 
312b7ce53b6SStefano Zampini     /* determining new matrix type */
313b2566f29SBarry Smith     ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
314b7ce53b6SStefano Zampini     if (issbaij_red) {
315b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
316b7ce53b6SStefano Zampini     } else {
317b7ce53b6SStefano Zampini       if (bs>1) {
318b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
319b7ce53b6SStefano Zampini       } else {
320b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
321b7ce53b6SStefano Zampini       }
322b7ce53b6SStefano Zampini     }
323b7ce53b6SStefano Zampini 
3243927de2eSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
3253cfa4ea4SStefano Zampini     ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr);
3263927de2eSStefano Zampini     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
3273927de2eSStefano Zampini     ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr);
3283927de2eSStefano Zampini     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
329b7ce53b6SStefano Zampini   } else {
3303cfa4ea4SStefano Zampini     PetscInt mbs,mrows,mcols,mlrows,mlcols;
331b7ce53b6SStefano Zampini     /* some checks */
332b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
333b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
3343cfa4ea4SStefano Zampini     ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr);
3356c4ed002SBarry Smith     if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
3366c4ed002SBarry Smith     if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
3376c4ed002SBarry Smith     if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows);
3386c4ed002SBarry Smith     if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols);
3396c4ed002SBarry Smith     if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
340b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
341b7ce53b6SStefano Zampini   }
342d9a9e74cSStefano Zampini 
343d9a9e74cSStefano Zampini   if (issbaij) {
344d9a9e74cSStefano Zampini     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr);
345d9a9e74cSStefano Zampini   } else {
346d9a9e74cSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
347d9a9e74cSStefano Zampini     local_mat = matis->A;
348d9a9e74cSStefano Zampini   }
349686e3a49SStefano Zampini 
350b7ce53b6SStefano Zampini   /* Set values */
351ecf5a873SStefano Zampini   ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr);
352b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
353ecf5a873SStefano Zampini     PetscInt i,*dummy_rows,*dummy_cols;
354ecf5a873SStefano Zampini 
355ecf5a873SStefano Zampini     if (local_rows != local_cols) {
356ecf5a873SStefano Zampini       ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr);
357ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
358ecf5a873SStefano Zampini       for (i=0;i<local_cols;i++) dummy_cols[i] = i;
359ecf5a873SStefano Zampini     } else {
360ecf5a873SStefano Zampini       ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr);
361ecf5a873SStefano Zampini       for (i=0;i<local_rows;i++) dummy_rows[i] = i;
362ecf5a873SStefano Zampini       dummy_cols = dummy_rows;
363ecf5a873SStefano Zampini     }
364b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
365d9a9e74cSStefano Zampini     ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr);
366ecf5a873SStefano Zampini     ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr);
367d9a9e74cSStefano Zampini     ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr);
368ecf5a873SStefano Zampini     if (dummy_rows != dummy_cols) {
369ecf5a873SStefano Zampini       ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr);
370ecf5a873SStefano Zampini     } else {
371ecf5a873SStefano Zampini       ierr = PetscFree(dummy_rows);CHKERRQ(ierr);
372ecf5a873SStefano Zampini     }
373686e3a49SStefano Zampini   } else if (isseqaij) {
374ecf5a873SStefano Zampini     PetscInt  i,nvtxs,*xadj,*adjncy;
375686e3a49SStefano Zampini     PetscBool done;
376686e3a49SStefano Zampini 
377d9a9e74cSStefano Zampini     ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
3786c4ed002SBarry Smith     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
379d9a9e74cSStefano Zampini     ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr);
380686e3a49SStefano Zampini     for (i=0;i<nvtxs;i++) {
381ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr);
382686e3a49SStefano Zampini     }
383d9a9e74cSStefano Zampini     ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
3846c4ed002SBarry Smith     if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
385d9a9e74cSStefano Zampini     ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr);
386686e3a49SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
387ecf5a873SStefano Zampini     PetscInt i;
388c0962df8SStefano Zampini 
389686e3a49SStefano Zampini     for (i=0;i<local_rows;i++) {
390686e3a49SStefano Zampini       PetscInt       j;
391ecf5a873SStefano Zampini       const PetscInt *local_indices_cols;
392686e3a49SStefano Zampini 
393ecf5a873SStefano Zampini       ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
394ecf5a873SStefano Zampini       ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr);
395ecf5a873SStefano Zampini       ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr);
396686e3a49SStefano Zampini     }
397b7ce53b6SStefano Zampini   }
398b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
399d9a9e74cSStefano Zampini   ierr = MatDestroy(&local_mat);CHKERRQ(ierr);
400b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
401b7ce53b6SStefano Zampini   if (isdense) {
402b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
403b7ce53b6SStefano Zampini   }
404b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
405b7ce53b6SStefano Zampini }
406b7ce53b6SStefano Zampini 
407b7ce53b6SStefano Zampini #undef __FUNCT__
408b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
409b7ce53b6SStefano Zampini /*@
410b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
411b7ce53b6SStefano Zampini 
412b7ce53b6SStefano Zampini   Input Parameter:
413b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
414b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
415b7ce53b6SStefano Zampini 
416b7ce53b6SStefano Zampini   Output Parameter:
417b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
418b7ce53b6SStefano Zampini 
419b7ce53b6SStefano Zampini   Level: developer
420b7ce53b6SStefano Zampini 
421eb82efa4SStefano Zampini   Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested.
422b7ce53b6SStefano Zampini 
423b7ce53b6SStefano Zampini .seealso: MATIS
424b7ce53b6SStefano Zampini @*/
425b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
426b7ce53b6SStefano Zampini {
427b7ce53b6SStefano Zampini   PetscErrorCode ierr;
428b7ce53b6SStefano Zampini 
429b7ce53b6SStefano Zampini   PetscFunctionBegin;
430b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
431b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
432b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
433b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
434b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
435b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
4366c4ed002SBarry Smith     if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix");
437b7ce53b6SStefano Zampini   }
438b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
439b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
440b7ce53b6SStefano Zampini }
441b7ce53b6SStefano Zampini 
442b7ce53b6SStefano Zampini #undef __FUNCT__
443ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
444ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
445ad6194a2SStefano Zampini {
446ad6194a2SStefano Zampini   PetscErrorCode ierr;
447ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
448ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
449ad6194a2SStefano Zampini   Mat            B,localmat;
450ad6194a2SStefano Zampini 
451ad6194a2SStefano Zampini   PetscFunctionBegin;
452ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
453ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
454ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
455e176bc59SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr);
456ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
457ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
458b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
459ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
460ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
461ad6194a2SStefano Zampini   *newmat = B;
462ad6194a2SStefano Zampini   PetscFunctionReturn(0);
463ad6194a2SStefano Zampini }
464ad6194a2SStefano Zampini 
465ad6194a2SStefano Zampini #undef __FUNCT__
46669796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
46769796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
46869796d55SStefano Zampini {
46969796d55SStefano Zampini   PetscErrorCode ierr;
47069796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
47169796d55SStefano Zampini   PetscBool      local_sym;
47269796d55SStefano Zampini 
47369796d55SStefano Zampini   PetscFunctionBegin;
47469796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
475b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
47669796d55SStefano Zampini   PetscFunctionReturn(0);
47769796d55SStefano Zampini }
47869796d55SStefano Zampini 
47969796d55SStefano Zampini #undef __FUNCT__
48069796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
48169796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
48269796d55SStefano Zampini {
48369796d55SStefano Zampini   PetscErrorCode ierr;
48469796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
48569796d55SStefano Zampini   PetscBool      local_sym;
48669796d55SStefano Zampini 
48769796d55SStefano Zampini   PetscFunctionBegin;
48869796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
489b2566f29SBarry Smith   ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
49069796d55SStefano Zampini   PetscFunctionReturn(0);
49169796d55SStefano Zampini }
49269796d55SStefano Zampini 
49369796d55SStefano Zampini #undef __FUNCT__
494b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
495dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
496b4319ba4SBarry Smith {
497dfbe8321SBarry Smith   PetscErrorCode ierr;
498b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
499b4319ba4SBarry Smith 
500b4319ba4SBarry Smith   PetscFunctionBegin;
5016bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
502e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr);
503e176bc59SStefano Zampini   ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr);
5046bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
5056bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
50628f4e0baSStefano Zampini   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
50728f4e0baSStefano Zampini   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
508bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
509dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
510bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
511b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
512b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
5132e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
514b4319ba4SBarry Smith   PetscFunctionReturn(0);
515b4319ba4SBarry Smith }
516b4319ba4SBarry Smith 
517b4319ba4SBarry Smith #undef __FUNCT__
518b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
519dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
520b4319ba4SBarry Smith {
521dfbe8321SBarry Smith   PetscErrorCode ierr;
522b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
523b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
524b4319ba4SBarry Smith 
525b4319ba4SBarry Smith   PetscFunctionBegin;
526b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
527e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
528e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
529b4319ba4SBarry Smith 
530b4319ba4SBarry Smith   /* multiply the local matrix */
531b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
532b4319ba4SBarry Smith 
533b4319ba4SBarry Smith   /* scatter product back into global memory */
5342dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
535e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
536e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
537b4319ba4SBarry Smith   PetscFunctionReturn(0);
538b4319ba4SBarry Smith }
539b4319ba4SBarry Smith 
540b4319ba4SBarry Smith #undef __FUNCT__
5412e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
542b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5432e74eeadSLisandro Dalcin {
544650997f4SStefano Zampini   Vec            temp_vec;
5452e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5462e74eeadSLisandro Dalcin 
5472e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
548650997f4SStefano Zampini   if (v3 != v2) {
549650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
550650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
551650997f4SStefano Zampini   } else {
552650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
553650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
554650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
555650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
556650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
557650997f4SStefano Zampini   }
5582e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5592e74eeadSLisandro Dalcin }
5602e74eeadSLisandro Dalcin 
5612e74eeadSLisandro Dalcin #undef __FUNCT__
5622e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
563e176bc59SStefano Zampini PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x)
5642e74eeadSLisandro Dalcin {
5652e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
5662e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5672e74eeadSLisandro Dalcin 
568e176bc59SStefano Zampini   PetscFunctionBegin;
5692e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
570e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
571e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5722e74eeadSLisandro Dalcin 
5732e74eeadSLisandro Dalcin   /* multiply the local matrix */
574e176bc59SStefano Zampini   ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr);
5752e74eeadSLisandro Dalcin 
5762e74eeadSLisandro Dalcin   /* scatter product back into global vector */
577e176bc59SStefano Zampini   ierr = VecSet(x,0);CHKERRQ(ierr);
578e176bc59SStefano Zampini   ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
579e176bc59SStefano Zampini   ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5802e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5812e74eeadSLisandro Dalcin }
5822e74eeadSLisandro Dalcin 
5832e74eeadSLisandro Dalcin #undef __FUNCT__
5842e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
5852e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5862e74eeadSLisandro Dalcin {
587650997f4SStefano Zampini   Vec            temp_vec;
5882e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5892e74eeadSLisandro Dalcin 
5902e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
591650997f4SStefano Zampini   if (v3 != v2) {
592650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
593650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
594650997f4SStefano Zampini   } else {
595650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
596650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
597650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
598650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
599650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
600650997f4SStefano Zampini   }
6012e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
6022e74eeadSLisandro Dalcin }
6032e74eeadSLisandro Dalcin 
6042e74eeadSLisandro Dalcin #undef __FUNCT__
605b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
606dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
607b4319ba4SBarry Smith {
608b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
609dfbe8321SBarry Smith   PetscErrorCode ierr;
610b4319ba4SBarry Smith   PetscViewer    sviewer;
611b4319ba4SBarry Smith 
612b4319ba4SBarry Smith   PetscFunctionBegin;
6133f08860eSBarry Smith   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
614b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
6153f08860eSBarry Smith   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
616b4319ba4SBarry Smith   PetscFunctionReturn(0);
617b4319ba4SBarry Smith }
618b4319ba4SBarry Smith 
619b4319ba4SBarry Smith #undef __FUNCT__
620b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
621784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
622b4319ba4SBarry Smith {
623dfbe8321SBarry Smith   PetscErrorCode ierr;
624e176bc59SStefano Zampini   PetscInt       nr,rbs,nc,cbs;
625b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
626b4319ba4SBarry Smith   IS             from,to;
627e176bc59SStefano Zampini   Vec            cglobal,rglobal;
628b4319ba4SBarry Smith 
629b4319ba4SBarry Smith   PetscFunctionBegin;
630784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
631e176bc59SStefano Zampini   PetscCheckSameComm(A,1,cmapping,3);
6323bbff08aSStefano Zampini   /* Destroy any previous data */
63370cf5478SStefano Zampini   ierr = VecDestroy(&is->x);CHKERRQ(ierr);
63470cf5478SStefano Zampini   ierr = VecDestroy(&is->y);CHKERRQ(ierr);
635e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr);
636e176bc59SStefano Zampini   ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr);
6371c47cb0fSStefano Zampini   ierr = MatDestroy(&is->A);CHKERRQ(ierr);
63828f4e0baSStefano Zampini   ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
63928f4e0baSStefano Zampini   ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
6403bbff08aSStefano Zampini 
6413bbff08aSStefano Zampini   /* Setup Layout and set local to global maps */
642fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
643fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
644fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
645fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
646b4319ba4SBarry Smith 
647b4319ba4SBarry Smith   /* Create the local matrix A */
648e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr);
649e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr);
650e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr);
651e176bc59SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr);
652f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
653e176bc59SStefano Zampini   ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr);
654e176bc59SStefano Zampini   ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr);
655e176bc59SStefano Zampini   ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr);
656ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
657ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
658b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
659b4319ba4SBarry Smith 
660b4319ba4SBarry Smith   /* Create the local work vectors */
6613bbff08aSStefano Zampini   ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr);
662b4319ba4SBarry Smith 
663e176bc59SStefano Zampini   /* setup the global to local scatters */
664e176bc59SStefano Zampini   ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr);
665e176bc59SStefano Zampini   ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr);
666784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
667e176bc59SStefano Zampini   ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr);
668e176bc59SStefano Zampini   if (rmapping != cmapping) {
669e176bc59SStefano Zampini     ierr = ISDestroy(&to);CHKERRQ(ierr);
670e176bc59SStefano Zampini     ierr = ISDestroy(&from);CHKERRQ(ierr);
671e176bc59SStefano Zampini     ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr);
672e176bc59SStefano Zampini     ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr);
673e176bc59SStefano Zampini     ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr);
674e176bc59SStefano Zampini   } else {
675e176bc59SStefano Zampini     ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr);
676e176bc59SStefano Zampini     is->cctx = is->rctx;
677e176bc59SStefano Zampini   }
678e176bc59SStefano Zampini   ierr = VecDestroy(&rglobal);CHKERRQ(ierr);
679e176bc59SStefano Zampini   ierr = VecDestroy(&cglobal);CHKERRQ(ierr);
6806bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
6816bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
682b4319ba4SBarry Smith   PetscFunctionReturn(0);
683b4319ba4SBarry Smith }
684b4319ba4SBarry Smith 
6852e74eeadSLisandro Dalcin #undef __FUNCT__
6862e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
6872e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
6882e74eeadSLisandro Dalcin {
6892e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
6902e74eeadSLisandro Dalcin   PetscInt       rows_l[2048],cols_l[2048];
6912e74eeadSLisandro Dalcin   PetscErrorCode ierr;
6922e74eeadSLisandro Dalcin 
6932e74eeadSLisandro Dalcin   PetscFunctionBegin;
6942e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
695f23aa3ddSBarry 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);
6962e74eeadSLisandro Dalcin #endif
6973bbff08aSStefano Zampini   ierr = ISG2LMapApply(mat->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr);
6983bbff08aSStefano Zampini   ierr = ISG2LMapApply(mat->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr);
6992e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
7002e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7012e74eeadSLisandro Dalcin }
7022e74eeadSLisandro Dalcin 
703b4319ba4SBarry Smith #undef __FUNCT__
704b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
70513f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
706b4319ba4SBarry Smith {
707dfbe8321SBarry Smith   PetscErrorCode ierr;
708b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
709b4319ba4SBarry Smith 
710b4319ba4SBarry Smith   PetscFunctionBegin;
711b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
712b4319ba4SBarry Smith   PetscFunctionReturn(0);
713b4319ba4SBarry Smith }
714b4319ba4SBarry Smith 
715b4319ba4SBarry Smith #undef __FUNCT__
716f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
717f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
718f0006bf2SLisandro Dalcin {
719f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
720f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
721f0006bf2SLisandro Dalcin 
722f0006bf2SLisandro Dalcin   PetscFunctionBegin;
723f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
724f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
725f0006bf2SLisandro Dalcin }
726f0006bf2SLisandro Dalcin 
727f0006bf2SLisandro Dalcin #undef __FUNCT__
7282e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
7292b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
7302e74eeadSLisandro Dalcin {
7310298fd71SBarry Smith   PetscInt       n_l = 0, *rows_l = NULL;
7322e74eeadSLisandro Dalcin   PetscErrorCode ierr;
7332e74eeadSLisandro Dalcin 
7342e74eeadSLisandro Dalcin   PetscFunctionBegin;
735ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
7362e74eeadSLisandro Dalcin   if (n) {
737785e854fSJed Brown     ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr);
7383bbff08aSStefano Zampini     ierr = ISGlobalToLocalMappingApply(A->rmap->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
7392e74eeadSLisandro Dalcin   }
7402b40b63fSBarry Smith   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
741c31cb41cSBarry Smith   ierr = PetscFree(rows_l);CHKERRQ(ierr);
7422e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7432e74eeadSLisandro Dalcin }
7442e74eeadSLisandro Dalcin 
7452e74eeadSLisandro Dalcin #undef __FUNCT__
746b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
7472b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
748b4319ba4SBarry Smith {
749b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
750dfbe8321SBarry Smith   PetscErrorCode ierr;
751f4df32b1SMatthew Knepley   PetscInt       i;
752b4319ba4SBarry Smith   PetscScalar    *array;
753b4319ba4SBarry Smith 
754b4319ba4SBarry Smith   PetscFunctionBegin;
755ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
756b4319ba4SBarry Smith   {
757b4319ba4SBarry Smith     /*
758b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
759b4319ba4SBarry Smith        work properly in the interface nodes.
760b4319ba4SBarry Smith     */
761b4319ba4SBarry Smith     Vec counter;
762e176bc59SStefano Zampini     ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr);
763e176bc59SStefano Zampini     ierr = VecSet(counter,0.);CHKERRQ(ierr);
764e176bc59SStefano Zampini     ierr = VecSet(is->y,1.);CHKERRQ(ierr);
765e176bc59SStefano Zampini     ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
766e176bc59SStefano Zampini     ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
767e176bc59SStefano Zampini     ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
768e176bc59SStefano Zampini     ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7696bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
770b4319ba4SBarry Smith   }
771958c9bccSBarry Smith   if (!n) {
772b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
773b4319ba4SBarry Smith   } else {
774b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
7752205254eSKarl Rupp 
776e176bc59SStefano Zampini     ierr = VecGetArray(is->y,&array);CHKERRQ(ierr);
7772b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
778b4319ba4SBarry Smith     for (i=0; i<n; i++) {
779f4df32b1SMatthew Knepley       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
780b4319ba4SBarry Smith     }
781b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
782b4319ba4SBarry Smith     ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
783e176bc59SStefano Zampini     ierr = VecRestoreArray(is->y,&array);CHKERRQ(ierr);
784b4319ba4SBarry Smith   }
785b4319ba4SBarry Smith   PetscFunctionReturn(0);
786b4319ba4SBarry Smith }
787b4319ba4SBarry Smith 
788b4319ba4SBarry Smith #undef __FUNCT__
789b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
790dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
791b4319ba4SBarry Smith {
792b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
793dfbe8321SBarry Smith   PetscErrorCode ierr;
794b4319ba4SBarry Smith 
795b4319ba4SBarry Smith   PetscFunctionBegin;
796b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
797b4319ba4SBarry Smith   PetscFunctionReturn(0);
798b4319ba4SBarry Smith }
799b4319ba4SBarry Smith 
800b4319ba4SBarry Smith #undef __FUNCT__
801b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
802dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
803b4319ba4SBarry Smith {
804b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
805dfbe8321SBarry Smith   PetscErrorCode ierr;
806b4319ba4SBarry Smith 
807b4319ba4SBarry Smith   PetscFunctionBegin;
808b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
809b4319ba4SBarry Smith   PetscFunctionReturn(0);
810b4319ba4SBarry Smith }
811b4319ba4SBarry Smith 
812b4319ba4SBarry Smith #undef __FUNCT__
813b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
8147087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
815b4319ba4SBarry Smith {
816b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
817b4319ba4SBarry Smith 
818b4319ba4SBarry Smith   PetscFunctionBegin;
819b4319ba4SBarry Smith   *local = is->A;
820b4319ba4SBarry Smith   PetscFunctionReturn(0);
821b4319ba4SBarry Smith }
822b4319ba4SBarry Smith 
823b4319ba4SBarry Smith #undef __FUNCT__
824b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
825b4319ba4SBarry Smith /*@
826b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
827b4319ba4SBarry Smith 
828b4319ba4SBarry Smith   Input Parameter:
829b4319ba4SBarry Smith .  mat - the matrix
830b4319ba4SBarry Smith 
831b4319ba4SBarry Smith   Output Parameter:
832eb82efa4SStefano Zampini .  local - the local matrix
833b4319ba4SBarry Smith 
834b4319ba4SBarry Smith   Level: advanced
835b4319ba4SBarry Smith 
836b4319ba4SBarry Smith   Notes:
837b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
838b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
839b4319ba4SBarry Smith   of the MatSetValues() operation.
840b4319ba4SBarry Smith 
84196a6f129SJed Brown   This function does not increase the reference count for the local Mat.  Do not destroy it and do not attempt to use
84296a6f129SJed Brown   your reference after destroying the parent mat.
84396a6f129SJed Brown 
844b4319ba4SBarry Smith .seealso: MATIS
845b4319ba4SBarry Smith @*/
8467087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
847b4319ba4SBarry Smith {
8484ac538c5SBarry Smith   PetscErrorCode ierr;
849b4319ba4SBarry Smith 
850b4319ba4SBarry Smith   PetscFunctionBegin;
8510700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
852b4319ba4SBarry Smith   PetscValidPointer(local,2);
8534ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
854b4319ba4SBarry Smith   PetscFunctionReturn(0);
855b4319ba4SBarry Smith }
856b4319ba4SBarry Smith 
8573b03a366Sstefano_zampini #undef __FUNCT__
8583b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
8593b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
8603b03a366Sstefano_zampini {
8613b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
8623b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
8633b03a366Sstefano_zampini   PetscErrorCode ierr;
8643b03a366Sstefano_zampini 
8653b03a366Sstefano_zampini   PetscFunctionBegin;
8664e4c7dbeSStefano Zampini   if (is->A) {
8673b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
8683b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
869f23aa3ddSBarry 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);
8704e4c7dbeSStefano Zampini   }
8713b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
8723b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
8733b03a366Sstefano_zampini   is->A = local;
8743b03a366Sstefano_zampini   PetscFunctionReturn(0);
8753b03a366Sstefano_zampini }
8763b03a366Sstefano_zampini 
8773b03a366Sstefano_zampini #undef __FUNCT__
8783b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
8793b03a366Sstefano_zampini /*@
880eb82efa4SStefano Zampini     MatISSetLocalMat - Replace the local matrix stored inside a MATIS object.
8813b03a366Sstefano_zampini 
8823b03a366Sstefano_zampini   Input Parameter:
8833b03a366Sstefano_zampini .  mat - the matrix
884eb82efa4SStefano Zampini .  local - the local matrix
8853b03a366Sstefano_zampini 
8863b03a366Sstefano_zampini   Output Parameter:
8873b03a366Sstefano_zampini 
8883b03a366Sstefano_zampini   Level: advanced
8893b03a366Sstefano_zampini 
8903b03a366Sstefano_zampini   Notes:
8913b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
8923b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
8933b03a366Sstefano_zampini 
8943b03a366Sstefano_zampini .seealso: MATIS
8953b03a366Sstefano_zampini @*/
8963b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
8973b03a366Sstefano_zampini {
8983b03a366Sstefano_zampini   PetscErrorCode ierr;
8993b03a366Sstefano_zampini 
9003b03a366Sstefano_zampini   PetscFunctionBegin;
9013b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
902b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
9033b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
9043b03a366Sstefano_zampini   PetscFunctionReturn(0);
9053b03a366Sstefano_zampini }
9063b03a366Sstefano_zampini 
9076726f965SBarry Smith #undef __FUNCT__
9086726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
9096726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
9106726f965SBarry Smith {
9116726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
9126726f965SBarry Smith   PetscErrorCode ierr;
9136726f965SBarry Smith 
9146726f965SBarry Smith   PetscFunctionBegin;
9156726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
9166726f965SBarry Smith   PetscFunctionReturn(0);
9176726f965SBarry Smith }
9186726f965SBarry Smith 
9196726f965SBarry Smith #undef __FUNCT__
9202e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
9212e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
9222e74eeadSLisandro Dalcin {
9232e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9242e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9252e74eeadSLisandro Dalcin 
9262e74eeadSLisandro Dalcin   PetscFunctionBegin;
9272e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
9282e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9292e74eeadSLisandro Dalcin }
9302e74eeadSLisandro Dalcin 
9312e74eeadSLisandro Dalcin #undef __FUNCT__
9322e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
9332e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
9342e74eeadSLisandro Dalcin {
9352e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9362e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9372e74eeadSLisandro Dalcin 
9382e74eeadSLisandro Dalcin   PetscFunctionBegin;
9392e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
940e176bc59SStefano Zampini   ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr);
9412e74eeadSLisandro Dalcin 
9422e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
9432e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
944e176bc59SStefano Zampini   ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
945e176bc59SStefano Zampini   ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9462e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9472e74eeadSLisandro Dalcin }
9482e74eeadSLisandro Dalcin 
9492e74eeadSLisandro Dalcin #undef __FUNCT__
9506726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
951ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
9526726f965SBarry Smith {
9536726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
9546726f965SBarry Smith   PetscErrorCode ierr;
9556726f965SBarry Smith 
9566726f965SBarry Smith   PetscFunctionBegin;
9574e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
9586726f965SBarry Smith   PetscFunctionReturn(0);
9596726f965SBarry Smith }
9606726f965SBarry Smith 
961284134d9SBarry Smith #undef __FUNCT__
962284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
963284134d9SBarry Smith /*@
964284134d9SBarry Smith     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
965284134d9SBarry Smith        process but not across processes.
966284134d9SBarry Smith 
967284134d9SBarry Smith    Input Parameters:
968284134d9SBarry Smith +     comm    - MPI communicator that will share the matrix
969e176bc59SStefano Zampini .     bs      - block size of the matrix
970df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
971e176bc59SStefano Zampini .     rmap    - local to global map for rows
972e176bc59SStefano Zampini -     cmap    - local to global map for cols
973284134d9SBarry Smith 
974284134d9SBarry Smith    Output Parameter:
975284134d9SBarry Smith .    A - the resulting matrix
976284134d9SBarry Smith 
9778e6c10adSSatish Balay    Level: advanced
9788e6c10adSSatish Balay 
979284134d9SBarry Smith    Notes: See MATIS for more details
9808cda6cd7SBarry Smith           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
981e176bc59SStefano Zampini           by that process. The sizes of rmap and cmap define the size of the local matrices.
982e176bc59SStefano Zampini           If either rmap or cmap are NULL, than the matrix is assumed to be square
983284134d9SBarry Smith 
984284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
985284134d9SBarry Smith @*/
986e176bc59SStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A)
987284134d9SBarry Smith {
988284134d9SBarry Smith   PetscErrorCode ierr;
989284134d9SBarry Smith 
990284134d9SBarry Smith   PetscFunctionBegin;
991e176bc59SStefano Zampini   if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping");
992284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
993284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
9944f619741Sstefano_zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
995284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
9963b03a366Sstefano_zampini   ierr = MatSetUp(*A);CHKERRQ(ierr);
997e176bc59SStefano Zampini   if (rmap && cmap) {
998e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr);
999e176bc59SStefano Zampini   } else if (!rmap) {
1000e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr);
1001e176bc59SStefano Zampini   } else {
1002e176bc59SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr);
1003e176bc59SStefano Zampini   }
1004284134d9SBarry Smith   PetscFunctionReturn(0);
1005284134d9SBarry Smith }
1006284134d9SBarry Smith 
1007b4319ba4SBarry Smith /*MC
1008eb82efa4SStefano Zampini    MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC).
1009b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
1010b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
1011b4319ba4SBarry Smith    product is handled "implicitly".
1012b4319ba4SBarry Smith 
1013b4319ba4SBarry Smith    Operations Provided:
10146726f965SBarry Smith +  MatMult()
10152e74eeadSLisandro Dalcin .  MatMultAdd()
10162e74eeadSLisandro Dalcin .  MatMultTranspose()
10172e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
10186726f965SBarry Smith .  MatZeroEntries()
10196726f965SBarry Smith .  MatSetOption()
10202e74eeadSLisandro Dalcin .  MatZeroRows()
10216726f965SBarry Smith .  MatZeroRowsLocal()
10222e74eeadSLisandro Dalcin .  MatSetValues()
10236726f965SBarry Smith .  MatSetValuesLocal()
10242e74eeadSLisandro Dalcin .  MatScale()
10252e74eeadSLisandro Dalcin .  MatGetDiagonal()
10266726f965SBarry Smith -  MatSetLocalToGlobalMapping()
1027b4319ba4SBarry Smith 
1028b4319ba4SBarry Smith    Options Database Keys:
1029b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
1030b4319ba4SBarry Smith 
1031b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
1032b4319ba4SBarry Smith 
1033b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
1034b4319ba4SBarry Smith 
1035b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
1036eb82efa4SStefano Zampini           MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation()
1037b4319ba4SBarry Smith 
1038b4319ba4SBarry Smith   Level: advanced
1039b4319ba4SBarry Smith 
1040eb82efa4SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC
1041b4319ba4SBarry Smith 
1042b4319ba4SBarry Smith M*/
1043b4319ba4SBarry Smith 
1044b4319ba4SBarry Smith #undef __FUNCT__
1045b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
10468cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
1047b4319ba4SBarry Smith {
1048dfbe8321SBarry Smith   PetscErrorCode ierr;
1049b4319ba4SBarry Smith   Mat_IS         *b;
1050b4319ba4SBarry Smith 
1051b4319ba4SBarry Smith   PetscFunctionBegin;
1052b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1053b4319ba4SBarry Smith   A->data = (void*)b;
1054b4319ba4SBarry Smith 
1055e176bc59SStefano Zampini   /* matrix ops */
1056e176bc59SStefano Zampini   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1057b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
10582e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
10592e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
10602e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1061b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1062b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
10632e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
1064b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1065f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
10662e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1067b4319ba4SBarry Smith   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
1068b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1069b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1070b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
10716726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
10722e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
10732e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
10746726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
107569796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
107669796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1077ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
1078b4319ba4SBarry Smith 
1079b7ce53b6SStefano Zampini   /* special MATIS functions */
1080bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1081bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1082b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
10832e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
108417667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1085b4319ba4SBarry Smith   PetscFunctionReturn(0);
1086b4319ba4SBarry Smith }
1087