xref: /petsc/src/mat/impls/is/matis.c (revision fc27028ae294fed1ae1964acd0cd129dfdb3cb8e)
1be1d678aSKris Buschelman 
2b4319ba4SBarry Smith /*
3b4319ba4SBarry Smith     Creates a matrix class for using the Neumann-Neumann type preconditioners.
4b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
5b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
6b4319ba4SBarry Smith    product is handled "implicitly".
7b4319ba4SBarry Smith 
8b4319ba4SBarry Smith      We provide:
9b4319ba4SBarry Smith          MatMult()
10b4319ba4SBarry Smith 
11b4319ba4SBarry Smith     Currently this allows for only one subdomain per processor.
12b4319ba4SBarry Smith 
13b4319ba4SBarry Smith */
14b4319ba4SBarry Smith 
15c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h>      /*I "petscmat.h" I*/
1628f4e0baSStefano Zampini 
17a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat);
18a72627d2SStefano Zampini 
1928f4e0baSStefano Zampini #undef __FUNCT__
2028f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private"
21a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B)
2228f4e0baSStefano Zampini {
2328f4e0baSStefano Zampini   Mat_IS         *matis = (Mat_IS*)(B->data);
2428f4e0baSStefano Zampini   const PetscInt *gidxs;
2528f4e0baSStefano Zampini   PetscErrorCode ierr;
2628f4e0baSStefano Zampini 
2728f4e0baSStefano Zampini   PetscFunctionBegin;
2828f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr);
2928f4e0baSStefano Zampini   ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr);
3028f4e0baSStefano Zampini   ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr);
3128f4e0baSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&gidxs);CHKERRQ(ierr);
3228f4e0baSStefano Zampini   /* PETSC_OWN_POINTER refers to ilocal which is NULL */
3328f4e0baSStefano Zampini   ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr);
3428f4e0baSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&gidxs);CHKERRQ(ierr);
3528f4e0baSStefano Zampini   ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr);
3628f4e0baSStefano Zampini   PetscFunctionReturn(0);
3728f4e0baSStefano Zampini }
382e1947a5SStefano Zampini 
392e1947a5SStefano Zampini #undef __FUNCT__
402e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation"
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;
962e1947a5SStefano Zampini   if (!matis->A) {
972e1947a5SStefano Zampini     SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping");
982e1947a5SStefano Zampini   }
9928f4e0baSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
10028f4e0baSStefano Zampini     ierr = MatISComputeSF_Private(B);CHKERRQ(ierr);
10128f4e0baSStefano Zampini   }
1022e1947a5SStefano Zampini   if (!d_nnz) {
10328f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz;
1042e1947a5SStefano Zampini   } else {
10528f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i];
1062e1947a5SStefano Zampini   }
1072e1947a5SStefano Zampini   if (!o_nnz) {
10828f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz;
1092e1947a5SStefano Zampini   } else {
11028f4e0baSStefano Zampini     for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i];
1112e1947a5SStefano Zampini   }
11228f4e0baSStefano Zampini   ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
11328f4e0baSStefano Zampini   ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr);
11428f4e0baSStefano Zampini   ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr);
11528f4e0baSStefano Zampini   ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr);
11628f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves;i++) {
11728f4e0baSStefano Zampini     matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols);
1182e1947a5SStefano Zampini   }
11928f4e0baSStefano Zampini   ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr);
12028f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
12128f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs;
1222e1947a5SStefano Zampini   }
12328f4e0baSStefano Zampini   ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
12428f4e0baSStefano Zampini   for (i=0;i<matis->sf_nleaves/bs;i++) {
12528f4e0baSStefano Zampini     matis->sf_leafdata[i] = matis->sf_leafdata[i]-i;
1262e1947a5SStefano Zampini   }
12728f4e0baSStefano Zampini   ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr);
1282e1947a5SStefano Zampini   PetscFunctionReturn(0);
1292e1947a5SStefano Zampini }
130b4319ba4SBarry Smith 
131b4319ba4SBarry Smith #undef __FUNCT__
1323927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private"
1333927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce)
1343927de2eSStefano Zampini {
1353927de2eSStefano Zampini   Mat_IS          *matis = (Mat_IS*)(A->data);
1363927de2eSStefano Zampini   PetscInt        *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership;
1373927de2eSStefano Zampini   const PetscInt* global_indices;
1383927de2eSStefano Zampini   PetscInt        i,j,bs,rows,cols;
1393927de2eSStefano Zampini   PetscInt        lrows,lcols;
1403927de2eSStefano Zampini   PetscInt        local_rows,local_cols;
1413927de2eSStefano Zampini   PetscMPIInt     nsubdomains;
1423927de2eSStefano Zampini   PetscBool       isdense,issbaij;
1433927de2eSStefano Zampini   PetscErrorCode  ierr;
1443927de2eSStefano Zampini 
1453927de2eSStefano Zampini   PetscFunctionBegin;
1463927de2eSStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr);
1473927de2eSStefano Zampini   ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr);
1483927de2eSStefano Zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1493927de2eSStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
1503927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
1513927de2eSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
1523927de2eSStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&global_indices);CHKERRQ(ierr);
1533927de2eSStefano Zampini   if (issbaij) {
1543927de2eSStefano Zampini     ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr);
1553927de2eSStefano Zampini   }
1563927de2eSStefano Zampini   /*
1573927de2eSStefano Zampini      An SF reduce is needed to sum up properly on shared interface dofs.
1583927de2eSStefano Zampini      Note that generally preallocation is not exact, since it overestimates nonzeros
1593927de2eSStefano Zampini   */
1603927de2eSStefano Zampini   if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */
1613927de2eSStefano Zampini     ierr = MatISComputeSF_Private(A);CHKERRQ(ierr);
1623927de2eSStefano Zampini   }
1633927de2eSStefano Zampini   ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr);
1643927de2eSStefano Zampini   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr);
1653927de2eSStefano Zampini   /* All processes need to compute entire row ownership */
1663927de2eSStefano Zampini   ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr);
1673927de2eSStefano Zampini   ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr);
1683927de2eSStefano Zampini   for (i=0;i<nsubdomains;i++) {
1693927de2eSStefano Zampini     for (j=mat_ranges[i];j<mat_ranges[i+1];j++) {
1703927de2eSStefano Zampini       row_ownership[j]=i;
1713927de2eSStefano Zampini     }
1723927de2eSStefano Zampini   }
1733927de2eSStefano Zampini 
1743927de2eSStefano Zampini   /*
1753927de2eSStefano Zampini      my_dnz and my_onz contains exact contribution to preallocation from each local mat
1763927de2eSStefano Zampini      then, they will be summed up properly. This way, preallocation is always sufficient
1773927de2eSStefano Zampini   */
1783927de2eSStefano Zampini   ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr);
1793927de2eSStefano Zampini   /* preallocation as a MATAIJ */
1803927de2eSStefano Zampini   if (isdense) { /* special case for dense local matrices */
1813927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
1823927de2eSStefano Zampini       PetscInt index_row = global_indices[i];
1833927de2eSStefano Zampini       for (j=i;j<local_rows;j++) {
1843927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
1853927de2eSStefano Zampini         PetscInt index_col = global_indices[j];
1863927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
1873927de2eSStefano Zampini           my_dnz[i] += 1;
1883927de2eSStefano Zampini         } else { /* offdiag block */
1893927de2eSStefano Zampini           my_onz[i] += 1;
1903927de2eSStefano Zampini         }
1913927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
1923927de2eSStefano Zampini         if (i != j) {
1933927de2eSStefano Zampini           owner = row_ownership[index_col];
1943927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
1953927de2eSStefano Zampini             my_dnz[j] += 1;
1963927de2eSStefano Zampini           } else {
1973927de2eSStefano Zampini             my_onz[j] += 1;
1983927de2eSStefano Zampini           }
1993927de2eSStefano Zampini         }
2003927de2eSStefano Zampini       }
2013927de2eSStefano Zampini     }
2023927de2eSStefano Zampini   } else {
2033927de2eSStefano Zampini     for (i=0;i<local_rows;i++) {
2043927de2eSStefano Zampini       const PetscInt *cols;
2053927de2eSStefano Zampini       PetscInt       ncols,index_row = global_indices[i];
2063927de2eSStefano Zampini       ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2073927de2eSStefano Zampini       for (j=0;j<ncols;j++) {
2083927de2eSStefano Zampini         PetscInt owner = row_ownership[index_row];
2093927de2eSStefano Zampini         PetscInt index_col = global_indices[cols[j]];
2103927de2eSStefano Zampini         if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */
2113927de2eSStefano Zampini           my_dnz[i] += 1;
2123927de2eSStefano Zampini         } else { /* offdiag block */
2133927de2eSStefano Zampini           my_onz[i] += 1;
2143927de2eSStefano Zampini         }
2153927de2eSStefano Zampini         /* same as before, interchanging rows and cols */
216d9a9e74cSStefano Zampini         if (issbaij && index_col != index_row) {
2173927de2eSStefano Zampini           owner = row_ownership[index_col];
2183927de2eSStefano Zampini           if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) {
219d9a9e74cSStefano Zampini             my_dnz[cols[j]] += 1;
2203927de2eSStefano Zampini           } else {
221d9a9e74cSStefano Zampini             my_onz[cols[j]] += 1;
2223927de2eSStefano Zampini           }
2233927de2eSStefano Zampini         }
2243927de2eSStefano Zampini       }
2253927de2eSStefano Zampini       ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr);
2263927de2eSStefano Zampini     }
2273927de2eSStefano Zampini   }
2283927de2eSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&global_indices);CHKERRQ(ierr);
2293927de2eSStefano Zampini   ierr = PetscFree(row_ownership);CHKERRQ(ierr);
2303927de2eSStefano Zampini   if (maxreduce) {
2313927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2323927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr);
2333927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2343927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr);
2353927de2eSStefano Zampini   } else {
2363927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2373927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr);
2383927de2eSStefano Zampini     ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2393927de2eSStefano Zampini     ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr);
2403927de2eSStefano Zampini   }
2413927de2eSStefano Zampini   ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr);
2423927de2eSStefano Zampini 
2433927de2eSStefano Zampini   /* Resize preallocation if overestimated */
2443927de2eSStefano Zampini   for (i=0;i<lrows;i++) {
2453927de2eSStefano Zampini     dnz[i] = PetscMin(dnz[i],lcols);
2463927de2eSStefano Zampini     onz[i] = PetscMin(onz[i],cols-lcols);
2473927de2eSStefano Zampini   }
2483927de2eSStefano Zampini   /* set preallocation */
2493927de2eSStefano Zampini   ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr);
2503927de2eSStefano Zampini   for (i=0;i<lrows/bs;i++) {
2513927de2eSStefano Zampini     dnz[i] = dnz[i*bs]/bs;
2523927de2eSStefano Zampini     onz[i] = onz[i*bs]/bs;
2533927de2eSStefano Zampini   }
2543927de2eSStefano Zampini   ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2553927de2eSStefano Zampini   ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr);
2563927de2eSStefano Zampini   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
2573927de2eSStefano Zampini   if (issbaij) {
2583927de2eSStefano Zampini     ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr);
2593927de2eSStefano Zampini   }
2603927de2eSStefano Zampini   PetscFunctionReturn(0);
2613927de2eSStefano Zampini }
2623927de2eSStefano Zampini 
2633927de2eSStefano Zampini #undef __FUNCT__
264b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS"
265b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M)
266b7ce53b6SStefano Zampini {
267b7ce53b6SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
268d9a9e74cSStefano Zampini   Mat            local_mat;
269b7ce53b6SStefano Zampini   /* info on mat */
2703cfa4ea4SStefano Zampini   PetscInt       bs,rows,cols,lrows,lcols;
271b7ce53b6SStefano Zampini   PetscInt       local_rows,local_cols;
272686e3a49SStefano Zampini   PetscBool      isdense,issbaij,isseqaij;
273686e3a49SStefano Zampini   const PetscInt *global_indices_rows;
2747c03b4e8SStefano Zampini   PetscMPIInt    nsubdomains;
275b7ce53b6SStefano Zampini   /* values insertion */
276b7ce53b6SStefano Zampini   PetscScalar    *array;
277b7ce53b6SStefano Zampini   /* work */
278b7ce53b6SStefano Zampini   PetscErrorCode ierr;
279b7ce53b6SStefano Zampini 
280b7ce53b6SStefano Zampini   PetscFunctionBegin;
281b7ce53b6SStefano Zampini   /* MISSING CHECKS
282b7ce53b6SStefano Zampini     - rectangular case not covered (it is not allowed by MATIS)
283b7ce53b6SStefano Zampini   */
284b7ce53b6SStefano Zampini   /* get info from mat */
2857c03b4e8SStefano Zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr);
2867c03b4e8SStefano Zampini   if (nsubdomains == 1) {
2877c03b4e8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
2887c03b4e8SStefano Zampini       ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr);
2897c03b4e8SStefano Zampini     } else {
2907c03b4e8SStefano Zampini       ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
2917c03b4e8SStefano Zampini     }
2927c03b4e8SStefano Zampini     PetscFunctionReturn(0);
2937c03b4e8SStefano Zampini   }
294b7ce53b6SStefano Zampini   ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr);
295b7ce53b6SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
2963cfa4ea4SStefano Zampini   ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr);
297b7ce53b6SStefano Zampini   ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr);
298b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr);
299686e3a49SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
300b7ce53b6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr);
301b7ce53b6SStefano Zampini 
302686e3a49SStefano Zampini   /* map indices of local mat rows to global values */
303686e3a49SStefano Zampini   ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&global_indices_rows);CHKERRQ(ierr);
304b7ce53b6SStefano Zampini 
305b7ce53b6SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
306b7ce53b6SStefano Zampini     MatType     new_mat_type;
3073927de2eSStefano Zampini     PetscBool   issbaij_red;
308b7ce53b6SStefano Zampini 
309b7ce53b6SStefano Zampini     /* determining new matrix type */
310b7ce53b6SStefano Zampini     ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
311b7ce53b6SStefano Zampini     if (issbaij_red) {
312b7ce53b6SStefano Zampini       new_mat_type = MATSBAIJ;
313b7ce53b6SStefano Zampini     } else {
314b7ce53b6SStefano Zampini       if (bs>1) {
315b7ce53b6SStefano Zampini         new_mat_type = MATBAIJ;
316b7ce53b6SStefano Zampini       } else {
317b7ce53b6SStefano Zampini         new_mat_type = MATAIJ;
318b7ce53b6SStefano Zampini       }
319b7ce53b6SStefano Zampini     }
320b7ce53b6SStefano Zampini 
3213927de2eSStefano Zampini     ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr);
3223cfa4ea4SStefano Zampini     ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr);
3233927de2eSStefano Zampini     ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr);
3243927de2eSStefano Zampini     ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr);
3253927de2eSStefano Zampini     ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr);
326b7ce53b6SStefano Zampini   } else {
3273cfa4ea4SStefano Zampini     PetscInt mbs,mrows,mcols,mlrows,mlcols;
328b7ce53b6SStefano Zampini     /* some checks */
329b7ce53b6SStefano Zampini     ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr);
330b7ce53b6SStefano Zampini     ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr);
3313cfa4ea4SStefano Zampini     ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr);
332b7ce53b6SStefano Zampini     if (mrows != rows) {
333b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows);
334b7ce53b6SStefano Zampini     }
3353cfa4ea4SStefano Zampini     if (mcols != cols) {
336b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols);
337b7ce53b6SStefano Zampini     }
3383cfa4ea4SStefano Zampini     if (mlrows != lrows) {
3393cfa4ea4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows);
3403cfa4ea4SStefano Zampini     }
3413cfa4ea4SStefano Zampini     if (mlcols != lcols) {
3423cfa4ea4SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols);
3433cfa4ea4SStefano Zampini     }
344b7ce53b6SStefano Zampini     if (mbs != bs) {
345b7ce53b6SStefano Zampini       SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs);
346b7ce53b6SStefano Zampini     }
347b7ce53b6SStefano Zampini     ierr = MatZeroEntries(*M);CHKERRQ(ierr);
348b7ce53b6SStefano Zampini   }
349d9a9e74cSStefano Zampini 
350d9a9e74cSStefano Zampini   if (issbaij) {
351d9a9e74cSStefano Zampini     ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr);
352d9a9e74cSStefano Zampini   } else {
353d9a9e74cSStefano Zampini     ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr);
354d9a9e74cSStefano Zampini     local_mat = matis->A;
355d9a9e74cSStefano Zampini   }
356686e3a49SStefano Zampini 
357b7ce53b6SStefano Zampini   /* Set values */
358b7ce53b6SStefano Zampini   if (isdense) { /* special case for dense local matrices */
359b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
360d9a9e74cSStefano Zampini     ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr);
361686e3a49SStefano Zampini     ierr = MatSetValues(*M,local_rows,global_indices_rows,local_cols,global_indices_rows,array,ADD_VALUES);CHKERRQ(ierr);
362d9a9e74cSStefano Zampini     ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr);
363686e3a49SStefano Zampini   } else if (isseqaij) {
364686e3a49SStefano Zampini     PetscInt  i,nvtxs,*xadj,*adjncy,*global_indices_cols;
365686e3a49SStefano Zampini     PetscBool done;
366686e3a49SStefano Zampini 
367d9a9e74cSStefano Zampini     ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
368686e3a49SStefano Zampini     if (!done) {
369d9a9e74cSStefano Zampini       SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
370b7ce53b6SStefano Zampini     }
371d9a9e74cSStefano Zampini     ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr);
372686e3a49SStefano Zampini     ierr = PetscMalloc1(xadj[nvtxs],&global_indices_cols);CHKERRQ(ierr);
373686e3a49SStefano Zampini     ierr = ISLocalToGlobalMappingApply(matis->mapping,xadj[nvtxs],adjncy,global_indices_cols);CHKERRQ(ierr);
374686e3a49SStefano Zampini     for (i=0;i<nvtxs;i++) {
375686e3a49SStefano Zampini       PetscInt    row,*cols,ncols;
376686e3a49SStefano Zampini       PetscScalar *mat_vals;
377686e3a49SStefano Zampini 
378686e3a49SStefano Zampini       row = global_indices_rows[i];
379686e3a49SStefano Zampini       ncols = xadj[i+1]-xadj[i];
380686e3a49SStefano Zampini       cols = global_indices_cols+xadj[i];
381686e3a49SStefano Zampini       mat_vals = array+xadj[i];
382686e3a49SStefano Zampini       ierr = MatSetValues(*M,1,&row,ncols,cols,mat_vals,ADD_VALUES);CHKERRQ(ierr);
383686e3a49SStefano Zampini     }
384d9a9e74cSStefano Zampini     ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr);
385686e3a49SStefano Zampini     if (!done) {
386d9a9e74cSStefano Zampini       SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__);
387686e3a49SStefano Zampini     }
388d9a9e74cSStefano Zampini     ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr);
389686e3a49SStefano Zampini     ierr = PetscFree(global_indices_cols);CHKERRQ(ierr);
390686e3a49SStefano Zampini   } else { /* very basic values insertion for all other matrix types */
391c0962df8SStefano Zampini     PetscInt i,*global_indices_cols;
392c0962df8SStefano Zampini 
393c0962df8SStefano Zampini     ierr = PetscMalloc1(local_cols,&global_indices_cols);CHKERRQ(ierr);
394686e3a49SStefano Zampini     for (i=0;i<local_rows;i++) {
395686e3a49SStefano Zampini       PetscInt       j;
396686e3a49SStefano Zampini       const PetscInt *local_indices;
397686e3a49SStefano Zampini 
398d9a9e74cSStefano Zampini       ierr = MatGetRow(local_mat,i,&j,&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
399c0962df8SStefano Zampini       ierr = ISLocalToGlobalMappingApply(matis->mapping,j,local_indices,global_indices_cols);CHKERRQ(ierr);
400c0962df8SStefano Zampini       ierr = MatSetValues(*M,1,&global_indices_rows[i],j,global_indices_cols,array,ADD_VALUES);CHKERRQ(ierr);
401d9a9e74cSStefano Zampini       ierr = MatRestoreRow(local_mat,i,&j,&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr);
402686e3a49SStefano Zampini     }
403c0962df8SStefano Zampini     ierr = PetscFree(global_indices_cols);CHKERRQ(ierr);
404b7ce53b6SStefano Zampini   }
405b7ce53b6SStefano Zampini   ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
406d9a9e74cSStefano Zampini   ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&global_indices_rows);CHKERRQ(ierr);
407d9a9e74cSStefano Zampini   ierr = MatDestroy(&local_mat);CHKERRQ(ierr);
408b7ce53b6SStefano Zampini   ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
409b7ce53b6SStefano Zampini   if (isdense) {
410b7ce53b6SStefano Zampini     ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
411b7ce53b6SStefano Zampini   }
412b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
413b7ce53b6SStefano Zampini }
414b7ce53b6SStefano Zampini 
415b7ce53b6SStefano Zampini #undef __FUNCT__
416b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ"
417b7ce53b6SStefano Zampini /*@
418b7ce53b6SStefano Zampini     MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format
419b7ce53b6SStefano Zampini 
420b7ce53b6SStefano Zampini   Input Parameter:
421b7ce53b6SStefano Zampini .  mat - the matrix (should be of type MATIS)
422b7ce53b6SStefano Zampini .  reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
423b7ce53b6SStefano Zampini 
424b7ce53b6SStefano Zampini   Output Parameter:
425b7ce53b6SStefano Zampini .  newmat - the matrix in AIJ format
426b7ce53b6SStefano Zampini 
427b7ce53b6SStefano Zampini   Level: developer
428b7ce53b6SStefano Zampini 
429eb82efa4SStefano Zampini   Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested.
430b7ce53b6SStefano Zampini 
431b7ce53b6SStefano Zampini .seealso: MATIS
432b7ce53b6SStefano Zampini @*/
433b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat)
434b7ce53b6SStefano Zampini {
435b7ce53b6SStefano Zampini   PetscErrorCode ierr;
436b7ce53b6SStefano Zampini 
437b7ce53b6SStefano Zampini   PetscFunctionBegin;
438b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
439b7ce53b6SStefano Zampini   PetscValidLogicalCollectiveEnum(mat,reuse,2);
440b7ce53b6SStefano Zampini   PetscValidPointer(newmat,3);
441b7ce53b6SStefano Zampini   if (reuse != MAT_INITIAL_MATRIX) {
442b7ce53b6SStefano Zampini     PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3);
443b7ce53b6SStefano Zampini     PetscCheckSameComm(mat,1,*newmat,3);
444eb82efa4SStefano Zampini     if (mat == *newmat) {
445eb82efa4SStefano Zampini       SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix");
446eb82efa4SStefano Zampini     }
447b7ce53b6SStefano Zampini   }
448b7ce53b6SStefano Zampini   ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr);
449b7ce53b6SStefano Zampini   PetscFunctionReturn(0);
450b7ce53b6SStefano Zampini }
451b7ce53b6SStefano Zampini 
452b7ce53b6SStefano Zampini #undef __FUNCT__
453ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS"
454ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat)
455ad6194a2SStefano Zampini {
456ad6194a2SStefano Zampini   PetscErrorCode ierr;
457ad6194a2SStefano Zampini   Mat_IS         *matis = (Mat_IS*)(mat->data);
458ad6194a2SStefano Zampini   PetscInt       bs,m,n,M,N;
459ad6194a2SStefano Zampini   Mat            B,localmat;
460ad6194a2SStefano Zampini 
461ad6194a2SStefano Zampini   PetscFunctionBegin;
462ad6194a2SStefano Zampini   ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
463ad6194a2SStefano Zampini   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
464ad6194a2SStefano Zampini   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
465ad6194a2SStefano Zampini   ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr);
466ad6194a2SStefano Zampini   ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr);
467ad6194a2SStefano Zampini   ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr);
468b3317aa8SStefano Zampini   ierr = MatDestroy(&localmat);CHKERRQ(ierr);
469ad6194a2SStefano Zampini   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
470ad6194a2SStefano Zampini   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
471ad6194a2SStefano Zampini   *newmat = B;
472ad6194a2SStefano Zampini   PetscFunctionReturn(0);
473ad6194a2SStefano Zampini }
474ad6194a2SStefano Zampini 
475ad6194a2SStefano Zampini #undef __FUNCT__
47669796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS"
47769796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool  *flg)
47869796d55SStefano Zampini {
47969796d55SStefano Zampini   PetscErrorCode ierr;
48069796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
48169796d55SStefano Zampini   PetscBool      local_sym;
48269796d55SStefano Zampini 
48369796d55SStefano Zampini   PetscFunctionBegin;
48469796d55SStefano Zampini   ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr);
48569796d55SStefano Zampini   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
48669796d55SStefano Zampini   PetscFunctionReturn(0);
48769796d55SStefano Zampini }
48869796d55SStefano Zampini 
48969796d55SStefano Zampini #undef __FUNCT__
49069796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS"
49169796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool  *flg)
49269796d55SStefano Zampini {
49369796d55SStefano Zampini   PetscErrorCode ierr;
49469796d55SStefano Zampini   Mat_IS         *matis = (Mat_IS*)A->data;
49569796d55SStefano Zampini   PetscBool      local_sym;
49669796d55SStefano Zampini 
49769796d55SStefano Zampini   PetscFunctionBegin;
49869796d55SStefano Zampini   ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr);
49969796d55SStefano Zampini   ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
50069796d55SStefano Zampini   PetscFunctionReturn(0);
50169796d55SStefano Zampini }
50269796d55SStefano Zampini 
50369796d55SStefano Zampini #undef __FUNCT__
504b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS"
505dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A)
506b4319ba4SBarry Smith {
507dfbe8321SBarry Smith   PetscErrorCode ierr;
508b4319ba4SBarry Smith   Mat_IS         *b = (Mat_IS*)A->data;
509b4319ba4SBarry Smith 
510b4319ba4SBarry Smith   PetscFunctionBegin;
5116bf464f9SBarry Smith   ierr = MatDestroy(&b->A);CHKERRQ(ierr);
5126bf464f9SBarry Smith   ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr);
5136bf464f9SBarry Smith   ierr = VecDestroy(&b->x);CHKERRQ(ierr);
5146bf464f9SBarry Smith   ierr = VecDestroy(&b->y);CHKERRQ(ierr);
5156bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr);
51628f4e0baSStefano Zampini   ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr);
51728f4e0baSStefano Zampini   ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr);
518bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
519dbd8c25aSHong Zhang   ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr);
520bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr);
521b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr);
522b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr);
5232e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr);
524b4319ba4SBarry Smith   PetscFunctionReturn(0);
525b4319ba4SBarry Smith }
526b4319ba4SBarry Smith 
527b4319ba4SBarry Smith #undef __FUNCT__
528b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS"
529dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y)
530b4319ba4SBarry Smith {
531dfbe8321SBarry Smith   PetscErrorCode ierr;
532b4319ba4SBarry Smith   Mat_IS         *is  = (Mat_IS*)A->data;
533b4319ba4SBarry Smith   PetscScalar    zero = 0.0;
534b4319ba4SBarry Smith 
535b4319ba4SBarry Smith   PetscFunctionBegin;
536b4319ba4SBarry Smith   /*  scatter the global vector x into the local work vector */
537ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
538ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
539b4319ba4SBarry Smith 
540b4319ba4SBarry Smith   /* multiply the local matrix */
541b4319ba4SBarry Smith   ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr);
542b4319ba4SBarry Smith 
543b4319ba4SBarry Smith   /* scatter product back into global memory */
5442dcb1b2aSMatthew Knepley   ierr = VecSet(y,zero);CHKERRQ(ierr);
545ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
546ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
547b4319ba4SBarry Smith   PetscFunctionReturn(0);
548b4319ba4SBarry Smith }
549b4319ba4SBarry Smith 
550b4319ba4SBarry Smith #undef __FUNCT__
5512e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS"
552b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5532e74eeadSLisandro Dalcin {
554650997f4SStefano Zampini   Vec            temp_vec;
5552e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5562e74eeadSLisandro Dalcin 
5572e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A * v1.*/
558650997f4SStefano Zampini   if (v3 != v2) {
559650997f4SStefano Zampini     ierr = MatMult(A,v1,v3);CHKERRQ(ierr);
560650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
561650997f4SStefano Zampini   } else {
562650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
563650997f4SStefano Zampini     ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr);
564650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
565650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
566650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
567650997f4SStefano Zampini   }
5682e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5692e74eeadSLisandro Dalcin }
5702e74eeadSLisandro Dalcin 
5712e74eeadSLisandro Dalcin #undef __FUNCT__
5722e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS"
5732e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y)
5742e74eeadSLisandro Dalcin {
5752e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
5762e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5772e74eeadSLisandro Dalcin 
5782e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  y = A' * x */
5792e74eeadSLisandro Dalcin   /*  scatter the global vector x into the local work vector */
580ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
581ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5822e74eeadSLisandro Dalcin 
5832e74eeadSLisandro Dalcin   /* multiply the local matrix */
5842e74eeadSLisandro Dalcin   ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr);
5852e74eeadSLisandro Dalcin 
5862e74eeadSLisandro Dalcin   /* scatter product back into global vector */
5872e74eeadSLisandro Dalcin   ierr = VecSet(y,0);CHKERRQ(ierr);
588ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
589ca9f406cSSatish Balay   ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5902e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
5912e74eeadSLisandro Dalcin }
5922e74eeadSLisandro Dalcin 
5932e74eeadSLisandro Dalcin #undef __FUNCT__
5942e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS"
5952e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3)
5962e74eeadSLisandro Dalcin {
597650997f4SStefano Zampini   Vec            temp_vec;
5982e74eeadSLisandro Dalcin   PetscErrorCode ierr;
5992e74eeadSLisandro Dalcin 
6002e74eeadSLisandro Dalcin   PetscFunctionBegin; /*  v3 = v2 + A' * v1.*/
601650997f4SStefano Zampini   if (v3 != v2) {
602650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr);
603650997f4SStefano Zampini     ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr);
604650997f4SStefano Zampini   } else {
605650997f4SStefano Zampini     ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr);
606650997f4SStefano Zampini     ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr);
607650997f4SStefano Zampini     ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr);
608650997f4SStefano Zampini     ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr);
609650997f4SStefano Zampini     ierr = VecDestroy(&temp_vec);CHKERRQ(ierr);
610650997f4SStefano Zampini   }
6112e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
6122e74eeadSLisandro Dalcin }
6132e74eeadSLisandro Dalcin 
6142e74eeadSLisandro Dalcin #undef __FUNCT__
615b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS"
616dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer)
617b4319ba4SBarry Smith {
618b4319ba4SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
619dfbe8321SBarry Smith   PetscErrorCode ierr;
620b4319ba4SBarry Smith   PetscViewer    sviewer;
621b4319ba4SBarry Smith 
622b4319ba4SBarry Smith   PetscFunctionBegin;
623dd2fa690SBarry Smith   ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
624b4319ba4SBarry Smith   ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
625b4319ba4SBarry Smith   ierr = MatView(a->A,sviewer);CHKERRQ(ierr);
626b4319ba4SBarry Smith   ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
627b4319ba4SBarry Smith   PetscFunctionReturn(0);
628b4319ba4SBarry Smith }
629b4319ba4SBarry Smith 
630b4319ba4SBarry Smith #undef __FUNCT__
631b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS"
632784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping)
633b4319ba4SBarry Smith {
634dfbe8321SBarry Smith   PetscErrorCode ierr;
6354e4c7dbeSStefano Zampini   PetscInt       n,bs;
636b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
637b4319ba4SBarry Smith   IS             from,to;
638b4319ba4SBarry Smith   Vec            global;
639b4319ba4SBarry Smith 
640b4319ba4SBarry Smith   PetscFunctionBegin;
641784ac674SJed Brown   PetscCheckSameComm(A,1,rmapping,2);
642ce94432eSBarry Smith   if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical");
64370cf5478SStefano Zampini   if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */
64470cf5478SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
64570cf5478SStefano Zampini     ierr = VecDestroy(&is->x);CHKERRQ(ierr);
64670cf5478SStefano Zampini     ierr = VecDestroy(&is->y);CHKERRQ(ierr);
64770cf5478SStefano Zampini     ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr);
6481c47cb0fSStefano Zampini     ierr = MatDestroy(&is->A);CHKERRQ(ierr);
64928f4e0baSStefano Zampini     ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr);
65028f4e0baSStefano Zampini     ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr);
65170cf5478SStefano Zampini   }
652784ac674SJed Brown   ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr);
6536bf464f9SBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr);
654784ac674SJed Brown   is->mapping = rmapping;
655*fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
656*fc27028aSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
657fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr);
658fa7f1dd8SStefano Zampini   ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr);
659b4319ba4SBarry Smith 
660b4319ba4SBarry Smith   /* Create the local matrix A */
661784ac674SJed Brown   ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr);
6622e1947a5SStefano Zampini   ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&bs);CHKERRQ(ierr);
663f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr);
6642e1947a5SStefano Zampini   if (bs > 1) {
6652e1947a5SStefano Zampini     ierr = MatSetType(is->A,MATSEQBAIJ);CHKERRQ(ierr);
6662e1947a5SStefano Zampini   } else {
6672e1947a5SStefano Zampini     ierr = MatSetType(is->A,MATSEQAIJ);CHKERRQ(ierr);
6682e1947a5SStefano Zampini   }
669f69a0ea3SMatthew Knepley   ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr);
6704e4c7dbeSStefano Zampini   ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr);
671ff130e51SJed Brown   ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr);
672ff130e51SJed Brown   ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr);
673b4319ba4SBarry Smith   ierr = MatSetFromOptions(is->A);CHKERRQ(ierr);
674b4319ba4SBarry Smith 
675b4319ba4SBarry Smith   /* Create the local work vectors */
6764e4c7dbeSStefano Zampini   ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr);
6774e4c7dbeSStefano Zampini   ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr);
6784e4c7dbeSStefano Zampini   ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr);
679ff130e51SJed Brown   ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr);
680ff130e51SJed Brown   ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr);
6814e4c7dbeSStefano Zampini   ierr = VecSetFromOptions(is->x);CHKERRQ(ierr);
682b4319ba4SBarry Smith   ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr);
683b4319ba4SBarry Smith 
684b4319ba4SBarry Smith   /* setup the global to local scatter */
685b4319ba4SBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr);
686784ac674SJed Brown   ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr);
6872a7a6963SBarry Smith   ierr = MatCreateVecs(A,&global,NULL);CHKERRQ(ierr);
688b4319ba4SBarry Smith   ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr);
6896bf464f9SBarry Smith   ierr = VecDestroy(&global);CHKERRQ(ierr);
6906bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
6916bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
692b4319ba4SBarry Smith   PetscFunctionReturn(0);
693b4319ba4SBarry Smith }
694b4319ba4SBarry Smith 
6952e74eeadSLisandro Dalcin #undef __FUNCT__
6962e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS"
6972e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv)
6982e74eeadSLisandro Dalcin {
6992e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)mat->data;
7002e74eeadSLisandro Dalcin   PetscInt       rows_l[2048],cols_l[2048];
7012e74eeadSLisandro Dalcin   PetscErrorCode ierr;
7022e74eeadSLisandro Dalcin 
7032e74eeadSLisandro Dalcin   PetscFunctionBegin;
7042e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
705f23aa3ddSBarry 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);
7062e74eeadSLisandro Dalcin #endif
7072e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr);
7082e74eeadSLisandro Dalcin   ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr);
7092e74eeadSLisandro Dalcin   ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr);
7102e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7112e74eeadSLisandro Dalcin }
7122e74eeadSLisandro Dalcin 
7132e74eeadSLisandro Dalcin #undef ISG2LMapSetUp
7142e74eeadSLisandro Dalcin #undef ISG2LMapApply
715b4319ba4SBarry Smith 
716b4319ba4SBarry Smith #undef __FUNCT__
717b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS"
71813f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
719b4319ba4SBarry Smith {
720dfbe8321SBarry Smith   PetscErrorCode ierr;
721b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
722b4319ba4SBarry Smith 
723b4319ba4SBarry Smith   PetscFunctionBegin;
724b4319ba4SBarry Smith   ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
725b4319ba4SBarry Smith   PetscFunctionReturn(0);
726b4319ba4SBarry Smith }
727b4319ba4SBarry Smith 
728b4319ba4SBarry Smith #undef __FUNCT__
729f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS"
730f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv)
731f0006bf2SLisandro Dalcin {
732f0006bf2SLisandro Dalcin   PetscErrorCode ierr;
733f0006bf2SLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
734f0006bf2SLisandro Dalcin 
735f0006bf2SLisandro Dalcin   PetscFunctionBegin;
736f0006bf2SLisandro Dalcin   ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr);
737f0006bf2SLisandro Dalcin   PetscFunctionReturn(0);
738f0006bf2SLisandro Dalcin }
739f0006bf2SLisandro Dalcin 
740f0006bf2SLisandro Dalcin #undef __FUNCT__
7412e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS"
7422b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
7432e74eeadSLisandro Dalcin {
7442e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
7450298fd71SBarry Smith   PetscInt       n_l = 0, *rows_l = NULL;
7462e74eeadSLisandro Dalcin   PetscErrorCode ierr;
7472e74eeadSLisandro Dalcin 
7482e74eeadSLisandro Dalcin   PetscFunctionBegin;
749ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
7502e74eeadSLisandro Dalcin   if (n) {
751785e854fSJed Brown     ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr);
7522e74eeadSLisandro Dalcin     ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr);
7532e74eeadSLisandro Dalcin   }
7542b40b63fSBarry Smith   ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr);
755c31cb41cSBarry Smith   ierr = PetscFree(rows_l);CHKERRQ(ierr);
7562e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
7572e74eeadSLisandro Dalcin }
7582e74eeadSLisandro Dalcin 
7592e74eeadSLisandro Dalcin #undef __FUNCT__
760b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS"
7612b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
762b4319ba4SBarry Smith {
763b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
764dfbe8321SBarry Smith   PetscErrorCode ierr;
765f4df32b1SMatthew Knepley   PetscInt       i;
766b4319ba4SBarry Smith   PetscScalar    *array;
767b4319ba4SBarry Smith 
768b4319ba4SBarry Smith   PetscFunctionBegin;
769ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support");
770b4319ba4SBarry Smith   {
771b4319ba4SBarry Smith     /*
772b4319ba4SBarry Smith        Set up is->x as a "counting vector". This is in order to MatMult_IS
773b4319ba4SBarry Smith        work properly in the interface nodes.
774b4319ba4SBarry Smith     */
775b4319ba4SBarry Smith     Vec         counter;
776b4319ba4SBarry Smith     PetscScalar one=1.0, zero=0.0;
7772a7a6963SBarry Smith     ierr = MatCreateVecs(A,&counter,NULL);CHKERRQ(ierr);
7782dcb1b2aSMatthew Knepley     ierr = VecSet(counter,zero);CHKERRQ(ierr);
7792dcb1b2aSMatthew Knepley     ierr = VecSet(is->x,one);CHKERRQ(ierr);
780ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
781ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
782ca9f406cSSatish Balay     ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
783ca9f406cSSatish Balay     ierr = VecScatterEnd  (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
7846bf464f9SBarry Smith     ierr = VecDestroy(&counter);CHKERRQ(ierr);
785b4319ba4SBarry Smith   }
786958c9bccSBarry Smith   if (!n) {
787b4319ba4SBarry Smith     is->pure_neumann = PETSC_TRUE;
788b4319ba4SBarry Smith   } else {
789b4319ba4SBarry Smith     is->pure_neumann = PETSC_FALSE;
7902205254eSKarl Rupp 
791b4319ba4SBarry Smith     ierr = VecGetArray(is->x,&array);CHKERRQ(ierr);
7922b40b63fSBarry Smith     ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr);
793b4319ba4SBarry Smith     for (i=0; i<n; i++) {
794f4df32b1SMatthew Knepley       ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr);
795b4319ba4SBarry Smith     }
796b4319ba4SBarry Smith     ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
797b4319ba4SBarry Smith     ierr = MatAssemblyEnd  (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
798b4319ba4SBarry Smith     ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr);
799b4319ba4SBarry Smith   }
800b4319ba4SBarry Smith   PetscFunctionReturn(0);
801b4319ba4SBarry Smith }
802b4319ba4SBarry Smith 
803b4319ba4SBarry Smith #undef __FUNCT__
804b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS"
805dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type)
806b4319ba4SBarry Smith {
807b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
808dfbe8321SBarry Smith   PetscErrorCode ierr;
809b4319ba4SBarry Smith 
810b4319ba4SBarry Smith   PetscFunctionBegin;
811b4319ba4SBarry Smith   ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr);
812b4319ba4SBarry Smith   PetscFunctionReturn(0);
813b4319ba4SBarry Smith }
814b4319ba4SBarry Smith 
815b4319ba4SBarry Smith #undef __FUNCT__
816b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS"
817dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type)
818b4319ba4SBarry Smith {
819b4319ba4SBarry Smith   Mat_IS         *is = (Mat_IS*)A->data;
820dfbe8321SBarry Smith   PetscErrorCode ierr;
821b4319ba4SBarry Smith 
822b4319ba4SBarry Smith   PetscFunctionBegin;
823b4319ba4SBarry Smith   ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr);
824b4319ba4SBarry Smith   PetscFunctionReturn(0);
825b4319ba4SBarry Smith }
826b4319ba4SBarry Smith 
827b4319ba4SBarry Smith #undef __FUNCT__
828b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS"
8297087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local)
830b4319ba4SBarry Smith {
831b4319ba4SBarry Smith   Mat_IS *is = (Mat_IS*)mat->data;
832b4319ba4SBarry Smith 
833b4319ba4SBarry Smith   PetscFunctionBegin;
834b4319ba4SBarry Smith   *local = is->A;
835b4319ba4SBarry Smith   PetscFunctionReturn(0);
836b4319ba4SBarry Smith }
837b4319ba4SBarry Smith 
838b4319ba4SBarry Smith #undef __FUNCT__
839b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat"
840b4319ba4SBarry Smith /*@
841b4319ba4SBarry Smith     MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix.
842b4319ba4SBarry Smith 
843b4319ba4SBarry Smith   Input Parameter:
844b4319ba4SBarry Smith .  mat - the matrix
845b4319ba4SBarry Smith 
846b4319ba4SBarry Smith   Output Parameter:
847eb82efa4SStefano Zampini .  local - the local matrix
848b4319ba4SBarry Smith 
849b4319ba4SBarry Smith   Level: advanced
850b4319ba4SBarry Smith 
851b4319ba4SBarry Smith   Notes:
852b4319ba4SBarry Smith     This can be called if you have precomputed the nonzero structure of the
853b4319ba4SBarry Smith   matrix and want to provide it to the inner matrix object to improve the performance
854b4319ba4SBarry Smith   of the MatSetValues() operation.
855b4319ba4SBarry Smith 
856b4319ba4SBarry Smith .seealso: MATIS
857b4319ba4SBarry Smith @*/
8587087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local)
859b4319ba4SBarry Smith {
8604ac538c5SBarry Smith   PetscErrorCode ierr;
861b4319ba4SBarry Smith 
862b4319ba4SBarry Smith   PetscFunctionBegin;
8630700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
864b4319ba4SBarry Smith   PetscValidPointer(local,2);
8654ac538c5SBarry Smith   ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr);
866b4319ba4SBarry Smith   PetscFunctionReturn(0);
867b4319ba4SBarry Smith }
868b4319ba4SBarry Smith 
8693b03a366Sstefano_zampini #undef __FUNCT__
8703b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS"
8713b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local)
8723b03a366Sstefano_zampini {
8733b03a366Sstefano_zampini   Mat_IS         *is = (Mat_IS*)mat->data;
8743b03a366Sstefano_zampini   PetscInt       nrows,ncols,orows,ocols;
8753b03a366Sstefano_zampini   PetscErrorCode ierr;
8763b03a366Sstefano_zampini 
8773b03a366Sstefano_zampini   PetscFunctionBegin;
8784e4c7dbeSStefano Zampini   if (is->A) {
8793b03a366Sstefano_zampini     ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr);
8803b03a366Sstefano_zampini     ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr);
881f23aa3ddSBarry 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);
8824e4c7dbeSStefano Zampini   }
8833b03a366Sstefano_zampini   ierr  = PetscObjectReference((PetscObject)local);CHKERRQ(ierr);
8843b03a366Sstefano_zampini   ierr  = MatDestroy(&is->A);CHKERRQ(ierr);
8853b03a366Sstefano_zampini   is->A = local;
8863b03a366Sstefano_zampini   PetscFunctionReturn(0);
8873b03a366Sstefano_zampini }
8883b03a366Sstefano_zampini 
8893b03a366Sstefano_zampini #undef __FUNCT__
8903b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat"
8913b03a366Sstefano_zampini /*@
892eb82efa4SStefano Zampini     MatISSetLocalMat - Replace the local matrix stored inside a MATIS object.
8933b03a366Sstefano_zampini 
8943b03a366Sstefano_zampini   Input Parameter:
8953b03a366Sstefano_zampini .  mat - the matrix
896eb82efa4SStefano Zampini .  local - the local matrix
8973b03a366Sstefano_zampini 
8983b03a366Sstefano_zampini   Output Parameter:
8993b03a366Sstefano_zampini 
9003b03a366Sstefano_zampini   Level: advanced
9013b03a366Sstefano_zampini 
9023b03a366Sstefano_zampini   Notes:
9033b03a366Sstefano_zampini     This can be called if you have precomputed the local matrix and
9043b03a366Sstefano_zampini   want to provide it to the matrix object MATIS.
9053b03a366Sstefano_zampini 
9063b03a366Sstefano_zampini .seealso: MATIS
9073b03a366Sstefano_zampini @*/
9083b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local)
9093b03a366Sstefano_zampini {
9103b03a366Sstefano_zampini   PetscErrorCode ierr;
9113b03a366Sstefano_zampini 
9123b03a366Sstefano_zampini   PetscFunctionBegin;
9133b03a366Sstefano_zampini   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
914b7ce53b6SStefano Zampini   PetscValidHeaderSpecific(local,MAT_CLASSID,2);
9153b03a366Sstefano_zampini   ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr);
9163b03a366Sstefano_zampini   PetscFunctionReturn(0);
9173b03a366Sstefano_zampini }
9183b03a366Sstefano_zampini 
9196726f965SBarry Smith #undef __FUNCT__
9206726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS"
9216726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A)
9226726f965SBarry Smith {
9236726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
9246726f965SBarry Smith   PetscErrorCode ierr;
9256726f965SBarry Smith 
9266726f965SBarry Smith   PetscFunctionBegin;
9276726f965SBarry Smith   ierr = MatZeroEntries(a->A);CHKERRQ(ierr);
9286726f965SBarry Smith   PetscFunctionReturn(0);
9296726f965SBarry Smith }
9306726f965SBarry Smith 
9316726f965SBarry Smith #undef __FUNCT__
9322e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS"
9332e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a)
9342e74eeadSLisandro Dalcin {
9352e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9362e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9372e74eeadSLisandro Dalcin 
9382e74eeadSLisandro Dalcin   PetscFunctionBegin;
9392e74eeadSLisandro Dalcin   ierr = MatScale(is->A,a);CHKERRQ(ierr);
9402e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9412e74eeadSLisandro Dalcin }
9422e74eeadSLisandro Dalcin 
9432e74eeadSLisandro Dalcin #undef __FUNCT__
9442e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS"
9452e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v)
9462e74eeadSLisandro Dalcin {
9472e74eeadSLisandro Dalcin   Mat_IS         *is = (Mat_IS*)A->data;
9482e74eeadSLisandro Dalcin   PetscErrorCode ierr;
9492e74eeadSLisandro Dalcin 
9502e74eeadSLisandro Dalcin   PetscFunctionBegin;
9512e74eeadSLisandro Dalcin   /* get diagonal of the local matrix */
9522e74eeadSLisandro Dalcin   ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr);
9532e74eeadSLisandro Dalcin 
9542e74eeadSLisandro Dalcin   /* scatter diagonal back into global vector */
9552e74eeadSLisandro Dalcin   ierr = VecSet(v,0);CHKERRQ(ierr);
956ca9f406cSSatish Balay   ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
957ca9f406cSSatish Balay   ierr = VecScatterEnd  (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
9582e74eeadSLisandro Dalcin   PetscFunctionReturn(0);
9592e74eeadSLisandro Dalcin }
9602e74eeadSLisandro Dalcin 
9612e74eeadSLisandro Dalcin #undef __FUNCT__
9626726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS"
963ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg)
9646726f965SBarry Smith {
9656726f965SBarry Smith   Mat_IS         *a = (Mat_IS*)A->data;
9666726f965SBarry Smith   PetscErrorCode ierr;
9676726f965SBarry Smith 
9686726f965SBarry Smith   PetscFunctionBegin;
9694e0d8c25SBarry Smith   ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr);
9706726f965SBarry Smith   PetscFunctionReturn(0);
9716726f965SBarry Smith }
9726726f965SBarry Smith 
973284134d9SBarry Smith #undef __FUNCT__
974284134d9SBarry Smith #define __FUNCT__ "MatCreateIS"
975284134d9SBarry Smith /*@
976284134d9SBarry Smith     MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each
977284134d9SBarry Smith        process but not across processes.
978284134d9SBarry Smith 
979284134d9SBarry Smith    Input Parameters:
980284134d9SBarry Smith +     comm - MPI communicator that will share the matrix
9814e4c7dbeSStefano Zampini .     bs - local and global block size of the matrix
982df3898eeSBarry Smith .     m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products
983284134d9SBarry Smith -     map - mapping that defines the global number for each local number
984284134d9SBarry Smith 
985284134d9SBarry Smith    Output Parameter:
986284134d9SBarry Smith .    A - the resulting matrix
987284134d9SBarry Smith 
9888e6c10adSSatish Balay    Level: advanced
9898e6c10adSSatish Balay 
990284134d9SBarry Smith    Notes: See MATIS for more details
9918cda6cd7SBarry Smith           m and n are NOT related to the size of the map, they are the size of the part of the vector owned
9928cda6cd7SBarry Smith           by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points
9938cda6cd7SBarry Smith           plus the ghost points to global indices.
994284134d9SBarry Smith 
995284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping()
996284134d9SBarry Smith @*/
9974e4c7dbeSStefano Zampini PetscErrorCode  MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A)
998284134d9SBarry Smith {
999284134d9SBarry Smith   PetscErrorCode ierr;
1000284134d9SBarry Smith 
1001284134d9SBarry Smith   PetscFunctionBegin;
1002284134d9SBarry Smith   ierr = MatCreate(comm,A);CHKERRQ(ierr);
1003d3ee2243SStefano Zampini   ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr);
1004284134d9SBarry Smith   ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr);
1005284134d9SBarry Smith   ierr = MatSetType(*A,MATIS);CHKERRQ(ierr);
10063b03a366Sstefano_zampini   ierr = MatSetUp(*A);CHKERRQ(ierr);
1007784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr);
1008284134d9SBarry Smith   PetscFunctionReturn(0);
1009284134d9SBarry Smith }
1010284134d9SBarry Smith 
1011b4319ba4SBarry Smith /*MC
1012eb82efa4SStefano Zampini    MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC).
1013b4319ba4SBarry Smith    This stores the matrices in globally unassembled form. Each processor
1014b4319ba4SBarry Smith    assembles only its local Neumann problem and the parallel matrix vector
1015b4319ba4SBarry Smith    product is handled "implicitly".
1016b4319ba4SBarry Smith 
1017b4319ba4SBarry Smith    Operations Provided:
10186726f965SBarry Smith +  MatMult()
10192e74eeadSLisandro Dalcin .  MatMultAdd()
10202e74eeadSLisandro Dalcin .  MatMultTranspose()
10212e74eeadSLisandro Dalcin .  MatMultTransposeAdd()
10226726f965SBarry Smith .  MatZeroEntries()
10236726f965SBarry Smith .  MatSetOption()
10242e74eeadSLisandro Dalcin .  MatZeroRows()
10256726f965SBarry Smith .  MatZeroRowsLocal()
10262e74eeadSLisandro Dalcin .  MatSetValues()
10276726f965SBarry Smith .  MatSetValuesLocal()
10282e74eeadSLisandro Dalcin .  MatScale()
10292e74eeadSLisandro Dalcin .  MatGetDiagonal()
10306726f965SBarry Smith -  MatSetLocalToGlobalMapping()
1031b4319ba4SBarry Smith 
1032b4319ba4SBarry Smith    Options Database Keys:
1033b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions()
1034b4319ba4SBarry Smith 
1035b4319ba4SBarry Smith    Notes: Options prefix for the inner matrix are given by -is_mat_xxx
1036b4319ba4SBarry Smith 
1037b4319ba4SBarry Smith           You must call MatSetLocalToGlobalMapping() before using this matrix type.
1038b4319ba4SBarry Smith 
1039b4319ba4SBarry Smith           You can do matrix preallocation on the local matrix after you obtain it with
1040eb82efa4SStefano Zampini           MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation()
1041b4319ba4SBarry Smith 
1042b4319ba4SBarry Smith   Level: advanced
1043b4319ba4SBarry Smith 
1044eb82efa4SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC
1045b4319ba4SBarry Smith 
1046b4319ba4SBarry Smith M*/
1047b4319ba4SBarry Smith 
1048b4319ba4SBarry Smith #undef __FUNCT__
1049b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS"
10508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A)
1051b4319ba4SBarry Smith {
1052dfbe8321SBarry Smith   PetscErrorCode ierr;
1053b4319ba4SBarry Smith   Mat_IS         *b;
1054b4319ba4SBarry Smith 
1055b4319ba4SBarry Smith   PetscFunctionBegin;
1056b00a9115SJed Brown   ierr    = PetscNewLog(A,&b);CHKERRQ(ierr);
1057b4319ba4SBarry Smith   A->data = (void*)b;
1058b4319ba4SBarry Smith   ierr    = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1059b4319ba4SBarry Smith 
1060b4319ba4SBarry Smith   A->ops->mult                    = MatMult_IS;
10612e74eeadSLisandro Dalcin   A->ops->multadd                 = MatMultAdd_IS;
10622e74eeadSLisandro Dalcin   A->ops->multtranspose           = MatMultTranspose_IS;
10632e74eeadSLisandro Dalcin   A->ops->multtransposeadd        = MatMultTransposeAdd_IS;
1064b4319ba4SBarry Smith   A->ops->destroy                 = MatDestroy_IS;
1065b4319ba4SBarry Smith   A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS;
10662e74eeadSLisandro Dalcin   A->ops->setvalues               = MatSetValues_IS;
1067b4319ba4SBarry Smith   A->ops->setvalueslocal          = MatSetValuesLocal_IS;
1068f0006bf2SLisandro Dalcin   A->ops->setvaluesblockedlocal   = MatSetValuesBlockedLocal_IS;
10692e74eeadSLisandro Dalcin   A->ops->zerorows                = MatZeroRows_IS;
1070b4319ba4SBarry Smith   A->ops->zerorowslocal           = MatZeroRowsLocal_IS;
1071b4319ba4SBarry Smith   A->ops->assemblybegin           = MatAssemblyBegin_IS;
1072b4319ba4SBarry Smith   A->ops->assemblyend             = MatAssemblyEnd_IS;
1073b4319ba4SBarry Smith   A->ops->view                    = MatView_IS;
10746726f965SBarry Smith   A->ops->zeroentries             = MatZeroEntries_IS;
10752e74eeadSLisandro Dalcin   A->ops->scale                   = MatScale_IS;
10762e74eeadSLisandro Dalcin   A->ops->getdiagonal             = MatGetDiagonal_IS;
10776726f965SBarry Smith   A->ops->setoption               = MatSetOption_IS;
107869796d55SStefano Zampini   A->ops->ishermitian             = MatIsHermitian_IS;
107969796d55SStefano Zampini   A->ops->issymmetric             = MatIsSymmetric_IS;
1080ad6194a2SStefano Zampini   A->ops->duplicate               = MatDuplicate_IS;
1081b4319ba4SBarry Smith 
1082b7ce53b6SStefano Zampini   /* zeros pointer */
1083b4319ba4SBarry Smith   b->A           = 0;
1084b4319ba4SBarry Smith   b->ctx         = 0;
1085b4319ba4SBarry Smith   b->x           = 0;
1086b4319ba4SBarry Smith   b->y           = 0;
1087b09366fdSStefano Zampini   b->mapping     = 0;
108828f4e0baSStefano Zampini   b->sf          = 0;
108928f4e0baSStefano Zampini   b->sf_rootdata = 0;
109028f4e0baSStefano Zampini   b->sf_leafdata = 0;
1091b7ce53b6SStefano Zampini 
1092b7ce53b6SStefano Zampini   /* special MATIS functions */
1093bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr);
1094bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr);
1095b7ce53b6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr);
10962e1947a5SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr);
109717667f90SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr);
1098b4319ba4SBarry Smith   PetscFunctionReturn(0);
1099b4319ba4SBarry Smith }
1100