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 17*cf0a3239SStefano Zampini #undef __FUNCT__ 18*cf0a3239SStefano Zampini #define __FUNCT__ "MatISSetUpSF" 19*cf0a3239SStefano Zampini /*@ 20*cf0a3239SStefano Zampini MatISSetUpSF - Setup star forest object used by MatIS during MatISSetPreallocation. 21*cf0a3239SStefano Zampini 22*cf0a3239SStefano Zampini Collective on MPI_Comm 23*cf0a3239SStefano Zampini 24*cf0a3239SStefano Zampini Input Parameters: 25*cf0a3239SStefano Zampini + A - the matrix 26*cf0a3239SStefano Zampini 27*cf0a3239SStefano Zampini Level: advanced 28*cf0a3239SStefano Zampini 29*cf0a3239SStefano Zampini Notes: This function does not be to be called by the user. 30*cf0a3239SStefano Zampini 31*cf0a3239SStefano Zampini .keywords: matrix 32*cf0a3239SStefano Zampini 33*cf0a3239SStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatISSetPreallocation(), MatISGetLocalMat() 34*cf0a3239SStefano Zampini @*/ 35*cf0a3239SStefano Zampini PetscErrorCode MatISSetUpSF(Mat A) 36*cf0a3239SStefano Zampini { 37*cf0a3239SStefano Zampini PetscErrorCode ierr; 38*cf0a3239SStefano Zampini 39*cf0a3239SStefano Zampini PetscFunctionBegin; 40*cf0a3239SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 41*cf0a3239SStefano Zampini PetscValidType(A,1); 42*cf0a3239SStefano Zampini ierr = PetscTryMethod(A,"MatISSetUpSF_C",(Mat),(A));CHKERRQ(ierr); 43*cf0a3239SStefano Zampini PetscFunctionReturn(0); 44*cf0a3239SStefano Zampini } 45a72627d2SStefano Zampini 4628f4e0baSStefano Zampini #undef __FUNCT__ 47*cf0a3239SStefano Zampini #define __FUNCT__ "MatISSetUpSF_IS" 48*cf0a3239SStefano Zampini static PetscErrorCode MatISSetUpSF_IS(Mat B) 4928f4e0baSStefano Zampini { 5028f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 5128f4e0baSStefano Zampini const PetscInt *gidxs; 5228f4e0baSStefano Zampini PetscErrorCode ierr; 5328f4e0baSStefano Zampini 5428f4e0baSStefano Zampini PetscFunctionBegin; 55*cf0a3239SStefano Zampini if (matis->sf) { 56*cf0a3239SStefano Zampini PetscFunctionReturn(0); 57*cf0a3239SStefano Zampini } 5828f4e0baSStefano Zampini ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr); 5928f4e0baSStefano Zampini ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr); 6028f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 613bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 6228f4e0baSStefano Zampini /* PETSC_OWN_POINTER refers to ilocal which is NULL */ 6328f4e0baSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 643bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 6528f4e0baSStefano Zampini ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 6628f4e0baSStefano Zampini PetscFunctionReturn(0); 6728f4e0baSStefano Zampini } 682e1947a5SStefano Zampini 692e1947a5SStefano Zampini #undef __FUNCT__ 702e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 71eb82efa4SStefano Zampini /*@ 72a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 73a88811baSStefano Zampini 74a88811baSStefano Zampini Collective on MPI_Comm 75a88811baSStefano Zampini 76a88811baSStefano Zampini Input Parameters: 77a88811baSStefano Zampini + B - the matrix 78a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 79a88811baSStefano Zampini (same value is used for all local rows) 80a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 81a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 82a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 83a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 84a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 85a88811baSStefano Zampini the diagonal entry even if it is zero. 86a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 87a88811baSStefano Zampini submatrix (same value is used for all local rows). 88a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 89a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 90a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 91a88811baSStefano Zampini structure. The size of this array is equal to the number 92a88811baSStefano Zampini of local rows, i.e 'm'. 93a88811baSStefano Zampini 94a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 95a88811baSStefano Zampini 96a88811baSStefano Zampini Level: intermediate 97a88811baSStefano Zampini 98a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 99a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 100a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 101a88811baSStefano Zampini 102a88811baSStefano Zampini .keywords: matrix 103a88811baSStefano Zampini 104a88811baSStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat() 105a88811baSStefano Zampini @*/ 1062e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 1072e1947a5SStefano Zampini { 1082e1947a5SStefano Zampini PetscErrorCode ierr; 1092e1947a5SStefano Zampini 1102e1947a5SStefano Zampini PetscFunctionBegin; 1112e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 1122e1947a5SStefano Zampini PetscValidType(B,1); 1132e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 1142e1947a5SStefano Zampini PetscFunctionReturn(0); 1152e1947a5SStefano Zampini } 1162e1947a5SStefano Zampini 1172e1947a5SStefano Zampini #undef __FUNCT__ 1182e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 1192e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 1202e1947a5SStefano Zampini { 1212e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 12228f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 1232e1947a5SStefano Zampini PetscErrorCode ierr; 1242e1947a5SStefano Zampini 1252e1947a5SStefano Zampini PetscFunctionBegin; 1262e1947a5SStefano Zampini if (!matis->A) { 1272e1947a5SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 1282e1947a5SStefano Zampini } 129*cf0a3239SStefano Zampini ierr = MatISSetUpSF(B);CHKERRQ(ierr); 1302e1947a5SStefano Zampini if (!d_nnz) { 13128f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz; 1322e1947a5SStefano Zampini } else { 13328f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 1342e1947a5SStefano Zampini } 1352e1947a5SStefano Zampini if (!o_nnz) { 13628f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz; 1372e1947a5SStefano Zampini } else { 13828f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 1392e1947a5SStefano Zampini } 14028f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 14128f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 14228f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 14328f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 14428f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves;i++) { 14528f4e0baSStefano Zampini matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 1462e1947a5SStefano Zampini } 14728f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 14828f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 14928f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 1502e1947a5SStefano Zampini } 15128f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 15228f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 15328f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 1542e1947a5SStefano Zampini } 15528f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 1562e1947a5SStefano Zampini PetscFunctionReturn(0); 1572e1947a5SStefano Zampini } 158b4319ba4SBarry Smith 159b4319ba4SBarry Smith #undef __FUNCT__ 1603927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private" 1613927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 1623927de2eSStefano Zampini { 1633927de2eSStefano Zampini Mat_IS *matis = (Mat_IS*)(A->data); 1643927de2eSStefano Zampini PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 165ecf5a873SStefano Zampini const PetscInt *global_indices_r,*global_indices_c; 1663927de2eSStefano Zampini PetscInt i,j,bs,rows,cols; 1673927de2eSStefano Zampini PetscInt lrows,lcols; 1683927de2eSStefano Zampini PetscInt local_rows,local_cols; 1693927de2eSStefano Zampini PetscMPIInt nsubdomains; 1703927de2eSStefano Zampini PetscBool isdense,issbaij; 1713927de2eSStefano Zampini PetscErrorCode ierr; 1723927de2eSStefano Zampini 1733927de2eSStefano Zampini PetscFunctionBegin; 1743927de2eSStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr); 1753927de2eSStefano Zampini ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr); 1763927de2eSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1773927de2eSStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 1783927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 1793927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 180ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 181ecf5a873SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 182ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_c);CHKERRQ(ierr); 183ecf5a873SStefano Zampini } else { 184ecf5a873SStefano Zampini global_indices_c = global_indices_r; 185ecf5a873SStefano Zampini } 186ecf5a873SStefano Zampini 1873927de2eSStefano Zampini if (issbaij) { 1883927de2eSStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 1893927de2eSStefano Zampini } 1903927de2eSStefano Zampini /* 191ecf5a873SStefano Zampini An SF reduce is needed to sum up properly on shared rows. 1923927de2eSStefano Zampini Note that generally preallocation is not exact, since it overestimates nonzeros 1933927de2eSStefano Zampini */ 194*cf0a3239SStefano Zampini ierr = MatISSetUpSF(A);CHKERRQ(ierr); 1953927de2eSStefano Zampini ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr); 1963927de2eSStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr); 1973927de2eSStefano Zampini /* All processes need to compute entire row ownership */ 1983927de2eSStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 1993927de2eSStefano Zampini ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 2003927de2eSStefano Zampini for (i=0;i<nsubdomains;i++) { 2013927de2eSStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 2023927de2eSStefano Zampini row_ownership[j] = i; 2033927de2eSStefano Zampini } 2043927de2eSStefano Zampini } 2053927de2eSStefano Zampini 2063927de2eSStefano Zampini /* 2073927de2eSStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 2083927de2eSStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 2093927de2eSStefano Zampini */ 2103927de2eSStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 2113927de2eSStefano Zampini /* preallocation as a MATAIJ */ 2123927de2eSStefano Zampini if (isdense) { /* special case for dense local matrices */ 2133927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 214ecf5a873SStefano Zampini PetscInt index_row = global_indices_r[i]; 2153927de2eSStefano Zampini for (j=i;j<local_rows;j++) { 2163927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 217ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[j]; 2183927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 2193927de2eSStefano Zampini my_dnz[i] += 1; 2203927de2eSStefano Zampini } else { /* offdiag block */ 2213927de2eSStefano Zampini my_onz[i] += 1; 2223927de2eSStefano Zampini } 2233927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 2243927de2eSStefano Zampini if (i != j) { 2253927de2eSStefano Zampini owner = row_ownership[index_col]; 2263927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 2273927de2eSStefano Zampini my_dnz[j] += 1; 2283927de2eSStefano Zampini } else { 2293927de2eSStefano Zampini my_onz[j] += 1; 2303927de2eSStefano Zampini } 2313927de2eSStefano Zampini } 2323927de2eSStefano Zampini } 2333927de2eSStefano Zampini } 234ecf5a873SStefano Zampini } else { /* TODO: this could be optimized using MatGetRowIJ */ 2353927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 2363927de2eSStefano Zampini const PetscInt *cols; 237ecf5a873SStefano Zampini PetscInt ncols,index_row = global_indices_r[i]; 2383927de2eSStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 2393927de2eSStefano Zampini for (j=0;j<ncols;j++) { 2403927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 241ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[cols[j]]; 2423927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 2433927de2eSStefano Zampini my_dnz[i] += 1; 2443927de2eSStefano Zampini } else { /* offdiag block */ 2453927de2eSStefano Zampini my_onz[i] += 1; 2463927de2eSStefano Zampini } 2473927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 248d9a9e74cSStefano Zampini if (issbaij && index_col != index_row) { 2493927de2eSStefano Zampini owner = row_ownership[index_col]; 2503927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 251d9a9e74cSStefano Zampini my_dnz[cols[j]] += 1; 2523927de2eSStefano Zampini } else { 253d9a9e74cSStefano Zampini my_onz[cols[j]] += 1; 2543927de2eSStefano Zampini } 2553927de2eSStefano Zampini } 2563927de2eSStefano Zampini } 2573927de2eSStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 2583927de2eSStefano Zampini } 2593927de2eSStefano Zampini } 260ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 261ecf5a873SStefano Zampini if (global_indices_c != global_indices_r) { 262ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_c);CHKERRQ(ierr); 263ecf5a873SStefano Zampini } 2643927de2eSStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 265ecf5a873SStefano Zampini 266ecf5a873SStefano Zampini /* Reduce my_dnz and my_onz */ 2673927de2eSStefano Zampini if (maxreduce) { 2683927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 2693927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 2703927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 2713927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 2723927de2eSStefano Zampini } else { 2733927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 2743927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 2753927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 2763927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 2773927de2eSStefano Zampini } 2783927de2eSStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 2793927de2eSStefano Zampini 2803927de2eSStefano Zampini /* Resize preallocation if overestimated */ 2813927de2eSStefano Zampini for (i=0;i<lrows;i++) { 2823927de2eSStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 2833927de2eSStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 2843927de2eSStefano Zampini } 2853927de2eSStefano Zampini /* set preallocation */ 2863927de2eSStefano Zampini ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 2873927de2eSStefano Zampini for (i=0;i<lrows/bs;i++) { 2883927de2eSStefano Zampini dnz[i] = dnz[i*bs]/bs; 2893927de2eSStefano Zampini onz[i] = onz[i*bs]/bs; 2903927de2eSStefano Zampini } 2913927de2eSStefano Zampini ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 2923927de2eSStefano Zampini ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 2933927de2eSStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 2943927de2eSStefano Zampini if (issbaij) { 2953927de2eSStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 2963927de2eSStefano Zampini } 2973927de2eSStefano Zampini PetscFunctionReturn(0); 2983927de2eSStefano Zampini } 2993927de2eSStefano Zampini 3003927de2eSStefano Zampini #undef __FUNCT__ 301b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 302b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 303b7ce53b6SStefano Zampini { 304b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 305d9a9e74cSStefano Zampini Mat local_mat; 306b7ce53b6SStefano Zampini /* info on mat */ 3073cfa4ea4SStefano Zampini PetscInt bs,rows,cols,lrows,lcols; 308b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 309686e3a49SStefano Zampini PetscBool isdense,issbaij,isseqaij; 3107c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 311b7ce53b6SStefano Zampini /* values insertion */ 312b7ce53b6SStefano Zampini PetscScalar *array; 313b7ce53b6SStefano Zampini /* work */ 314b7ce53b6SStefano Zampini PetscErrorCode ierr; 315b7ce53b6SStefano Zampini 316b7ce53b6SStefano Zampini PetscFunctionBegin; 317b7ce53b6SStefano Zampini /* get info from mat */ 3187c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 3197c03b4e8SStefano Zampini if (nsubdomains == 1) { 3207c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 3217c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 3227c03b4e8SStefano Zampini } else { 3237c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 3247c03b4e8SStefano Zampini } 3257c03b4e8SStefano Zampini PetscFunctionReturn(0); 3267c03b4e8SStefano Zampini } 327b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 328b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 3293cfa4ea4SStefano Zampini ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 330b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 331b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 332686e3a49SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 333b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 334b7ce53b6SStefano Zampini 335b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 336b7ce53b6SStefano Zampini MatType new_mat_type; 3373927de2eSStefano Zampini PetscBool issbaij_red; 338b7ce53b6SStefano Zampini 339b7ce53b6SStefano Zampini /* determining new matrix type */ 340b7ce53b6SStefano Zampini ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 341b7ce53b6SStefano Zampini if (issbaij_red) { 342b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 343b7ce53b6SStefano Zampini } else { 344b7ce53b6SStefano Zampini if (bs>1) { 345b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 346b7ce53b6SStefano Zampini } else { 347b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 348b7ce53b6SStefano Zampini } 349b7ce53b6SStefano Zampini } 350b7ce53b6SStefano Zampini 3513927de2eSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 3523cfa4ea4SStefano Zampini ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 3533927de2eSStefano Zampini ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 3543927de2eSStefano Zampini ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr); 3553927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 356b7ce53b6SStefano Zampini } else { 3573cfa4ea4SStefano Zampini PetscInt mbs,mrows,mcols,mlrows,mlcols; 358b7ce53b6SStefano Zampini /* some checks */ 359b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 360b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 3613cfa4ea4SStefano Zampini ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 362b7ce53b6SStefano Zampini if (mrows != rows) { 363b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 364b7ce53b6SStefano Zampini } 3653cfa4ea4SStefano Zampini if (mcols != cols) { 366b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 367b7ce53b6SStefano Zampini } 3683cfa4ea4SStefano Zampini if (mlrows != lrows) { 3693cfa4ea4SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 3703cfa4ea4SStefano Zampini } 3713cfa4ea4SStefano Zampini if (mlcols != lcols) { 3723cfa4ea4SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 3733cfa4ea4SStefano Zampini } 374b7ce53b6SStefano Zampini if (mbs != bs) { 375b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 376b7ce53b6SStefano Zampini } 377b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 378b7ce53b6SStefano Zampini } 379d9a9e74cSStefano Zampini 380d9a9e74cSStefano Zampini if (issbaij) { 381d9a9e74cSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 382d9a9e74cSStefano Zampini } else { 383d9a9e74cSStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 384d9a9e74cSStefano Zampini local_mat = matis->A; 385d9a9e74cSStefano Zampini } 386686e3a49SStefano Zampini 387b7ce53b6SStefano Zampini /* Set values */ 388ecf5a873SStefano Zampini ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 389b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 390ecf5a873SStefano Zampini PetscInt i,*dummy_rows,*dummy_cols; 391ecf5a873SStefano Zampini 392ecf5a873SStefano Zampini if (local_rows != local_cols) { 393ecf5a873SStefano Zampini ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 394ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 395ecf5a873SStefano Zampini for (i=0;i<local_cols;i++) dummy_cols[i] = i; 396ecf5a873SStefano Zampini } else { 397ecf5a873SStefano Zampini ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 398ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 399ecf5a873SStefano Zampini dummy_cols = dummy_rows; 400ecf5a873SStefano Zampini } 401b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 402d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 403ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 404d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 405ecf5a873SStefano Zampini if (dummy_rows != dummy_cols) { 406ecf5a873SStefano Zampini ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 407ecf5a873SStefano Zampini } else { 408ecf5a873SStefano Zampini ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 409ecf5a873SStefano Zampini } 410686e3a49SStefano Zampini } else if (isseqaij) { 411ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 412686e3a49SStefano Zampini PetscBool done; 413686e3a49SStefano Zampini 414d9a9e74cSStefano Zampini ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 415686e3a49SStefano Zampini if (!done) { 416d9a9e74cSStefano Zampini SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 417b7ce53b6SStefano Zampini } 418d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 419686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 420ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 421686e3a49SStefano Zampini } 422d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 423686e3a49SStefano Zampini if (!done) { 424d9a9e74cSStefano Zampini SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 425686e3a49SStefano Zampini } 426d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 427686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 428ecf5a873SStefano Zampini PetscInt i; 429c0962df8SStefano Zampini 430686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 431686e3a49SStefano Zampini PetscInt j; 432ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 433686e3a49SStefano Zampini 434ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 435ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 436ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 437686e3a49SStefano Zampini } 438b7ce53b6SStefano Zampini } 439b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 440d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 441b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 442b7ce53b6SStefano Zampini if (isdense) { 443b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 444b7ce53b6SStefano Zampini } 445b7ce53b6SStefano Zampini PetscFunctionReturn(0); 446b7ce53b6SStefano Zampini } 447b7ce53b6SStefano Zampini 448b7ce53b6SStefano Zampini #undef __FUNCT__ 449b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 450b7ce53b6SStefano Zampini /*@ 451b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 452b7ce53b6SStefano Zampini 453b7ce53b6SStefano Zampini Input Parameter: 454b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 455b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 456b7ce53b6SStefano Zampini 457b7ce53b6SStefano Zampini Output Parameter: 458b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 459b7ce53b6SStefano Zampini 460b7ce53b6SStefano Zampini Level: developer 461b7ce53b6SStefano Zampini 462eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 463b7ce53b6SStefano Zampini 464b7ce53b6SStefano Zampini .seealso: MATIS 465b7ce53b6SStefano Zampini @*/ 466b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 467b7ce53b6SStefano Zampini { 468b7ce53b6SStefano Zampini PetscErrorCode ierr; 469b7ce53b6SStefano Zampini 470b7ce53b6SStefano Zampini PetscFunctionBegin; 471b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 472b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 473b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 474b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 475b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 476b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 477eb82efa4SStefano Zampini if (mat == *newmat) { 478eb82efa4SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 479eb82efa4SStefano Zampini } 480b7ce53b6SStefano Zampini } 481b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 482b7ce53b6SStefano Zampini PetscFunctionReturn(0); 483b7ce53b6SStefano Zampini } 484b7ce53b6SStefano Zampini 485b7ce53b6SStefano Zampini #undef __FUNCT__ 486ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 487ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 488ad6194a2SStefano Zampini { 489ad6194a2SStefano Zampini PetscErrorCode ierr; 490ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 491ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 492ad6194a2SStefano Zampini Mat B,localmat; 493ad6194a2SStefano Zampini 494ad6194a2SStefano Zampini PetscFunctionBegin; 495ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 496ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 497ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 498e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 499ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 500ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 501b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 502ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 503ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 504ad6194a2SStefano Zampini *newmat = B; 505ad6194a2SStefano Zampini PetscFunctionReturn(0); 506ad6194a2SStefano Zampini } 507ad6194a2SStefano Zampini 508ad6194a2SStefano Zampini #undef __FUNCT__ 50969796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 51069796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 51169796d55SStefano Zampini { 51269796d55SStefano Zampini PetscErrorCode ierr; 51369796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 51469796d55SStefano Zampini PetscBool local_sym; 51569796d55SStefano Zampini 51669796d55SStefano Zampini PetscFunctionBegin; 51769796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 51869796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 51969796d55SStefano Zampini PetscFunctionReturn(0); 52069796d55SStefano Zampini } 52169796d55SStefano Zampini 52269796d55SStefano Zampini #undef __FUNCT__ 52369796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 52469796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 52569796d55SStefano Zampini { 52669796d55SStefano Zampini PetscErrorCode ierr; 52769796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 52869796d55SStefano Zampini PetscBool local_sym; 52969796d55SStefano Zampini 53069796d55SStefano Zampini PetscFunctionBegin; 53169796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 53269796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 53369796d55SStefano Zampini PetscFunctionReturn(0); 53469796d55SStefano Zampini } 53569796d55SStefano Zampini 53669796d55SStefano Zampini #undef __FUNCT__ 537b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 538dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 539b4319ba4SBarry Smith { 540dfbe8321SBarry Smith PetscErrorCode ierr; 541b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 542b4319ba4SBarry Smith 543b4319ba4SBarry Smith PetscFunctionBegin; 5446bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 545e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 546e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 5476bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 5486bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 54928f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 55028f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 551bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 552dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 553bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 554b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 555b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 5562e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 557*cf0a3239SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",NULL);CHKERRQ(ierr); 558b4319ba4SBarry Smith PetscFunctionReturn(0); 559b4319ba4SBarry Smith } 560b4319ba4SBarry Smith 561b4319ba4SBarry Smith #undef __FUNCT__ 562b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 563dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 564b4319ba4SBarry Smith { 565dfbe8321SBarry Smith PetscErrorCode ierr; 566b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 567b4319ba4SBarry Smith PetscScalar zero = 0.0; 568b4319ba4SBarry Smith 569b4319ba4SBarry Smith PetscFunctionBegin; 570b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 571e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 572e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 573b4319ba4SBarry Smith 574b4319ba4SBarry Smith /* multiply the local matrix */ 575b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 576b4319ba4SBarry Smith 577b4319ba4SBarry Smith /* scatter product back into global memory */ 5782dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 579e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 580e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 581b4319ba4SBarry Smith PetscFunctionReturn(0); 582b4319ba4SBarry Smith } 583b4319ba4SBarry Smith 584b4319ba4SBarry Smith #undef __FUNCT__ 5852e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 586b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5872e74eeadSLisandro Dalcin { 588650997f4SStefano Zampini Vec temp_vec; 5892e74eeadSLisandro Dalcin PetscErrorCode ierr; 5902e74eeadSLisandro Dalcin 5912e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 592650997f4SStefano Zampini if (v3 != v2) { 593650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 594650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 595650997f4SStefano Zampini } else { 596650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 597650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 598650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 599650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 600650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 601650997f4SStefano Zampini } 6022e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6032e74eeadSLisandro Dalcin } 6042e74eeadSLisandro Dalcin 6052e74eeadSLisandro Dalcin #undef __FUNCT__ 6062e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 607e176bc59SStefano Zampini PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 6082e74eeadSLisandro Dalcin { 6092e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 6102e74eeadSLisandro Dalcin PetscErrorCode ierr; 6112e74eeadSLisandro Dalcin 612e176bc59SStefano Zampini PetscFunctionBegin; 6132e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 614e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 615e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6162e74eeadSLisandro Dalcin 6172e74eeadSLisandro Dalcin /* multiply the local matrix */ 618e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 6192e74eeadSLisandro Dalcin 6202e74eeadSLisandro Dalcin /* scatter product back into global vector */ 621e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 622e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 623e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6242e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6252e74eeadSLisandro Dalcin } 6262e74eeadSLisandro Dalcin 6272e74eeadSLisandro Dalcin #undef __FUNCT__ 6282e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 6292e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 6302e74eeadSLisandro Dalcin { 631650997f4SStefano Zampini Vec temp_vec; 6322e74eeadSLisandro Dalcin PetscErrorCode ierr; 6332e74eeadSLisandro Dalcin 6342e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 635650997f4SStefano Zampini if (v3 != v2) { 636650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 637650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 638650997f4SStefano Zampini } else { 639650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 640650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 641650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 642650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 643650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 644650997f4SStefano Zampini } 6452e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6462e74eeadSLisandro Dalcin } 6472e74eeadSLisandro Dalcin 6482e74eeadSLisandro Dalcin #undef __FUNCT__ 649b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 650dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 651b4319ba4SBarry Smith { 652b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 653dfbe8321SBarry Smith PetscErrorCode ierr; 654b4319ba4SBarry Smith PetscViewer sviewer; 655b4319ba4SBarry Smith 656b4319ba4SBarry Smith PetscFunctionBegin; 6573f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 658b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 6593f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 660b4319ba4SBarry Smith PetscFunctionReturn(0); 661b4319ba4SBarry Smith } 662b4319ba4SBarry Smith 663b4319ba4SBarry Smith #undef __FUNCT__ 664b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 665784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 666b4319ba4SBarry Smith { 667dfbe8321SBarry Smith PetscErrorCode ierr; 668e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 669b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 670b4319ba4SBarry Smith IS from,to; 671e176bc59SStefano Zampini Vec cglobal,rglobal; 672b4319ba4SBarry Smith 673b4319ba4SBarry Smith PetscFunctionBegin; 674784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 675e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 6763bbff08aSStefano Zampini /* Destroy any previous data */ 67770cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 67870cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 679e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 680e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 6811c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 68228f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 68328f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 6843bbff08aSStefano Zampini 6853bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 686fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 687fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 688fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 689fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 690b4319ba4SBarry Smith 691b4319ba4SBarry Smith /* Create the local matrix A */ 692e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 693e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 694e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 695e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 696f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 697e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 698e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 699e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 700ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 701ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 702b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 703b4319ba4SBarry Smith 704b4319ba4SBarry Smith /* Create the local work vectors */ 7053bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 706b4319ba4SBarry Smith 707e176bc59SStefano Zampini /* setup the global to local scatters */ 708e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 709e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 710784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 711e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 712e176bc59SStefano Zampini if (rmapping != cmapping) { 713e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 714e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 715e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 716e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 717e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 718e176bc59SStefano Zampini } else { 719e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 720e176bc59SStefano Zampini is->cctx = is->rctx; 721e176bc59SStefano Zampini } 722e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 723e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 7246bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 7256bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 726b4319ba4SBarry Smith PetscFunctionReturn(0); 727b4319ba4SBarry Smith } 728b4319ba4SBarry Smith 7292e74eeadSLisandro Dalcin #undef __FUNCT__ 7302e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 7312e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 7322e74eeadSLisandro Dalcin { 7332e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 7342e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 7352e74eeadSLisandro Dalcin PetscErrorCode ierr; 7362e74eeadSLisandro Dalcin 7372e74eeadSLisandro Dalcin PetscFunctionBegin; 7382e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 739f23aa3ddSBarry 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); 7402e74eeadSLisandro Dalcin #endif 7413bbff08aSStefano Zampini ierr = ISG2LMapApply(mat->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 7423bbff08aSStefano Zampini ierr = ISG2LMapApply(mat->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 7432e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 7442e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7452e74eeadSLisandro Dalcin } 7462e74eeadSLisandro Dalcin 747b4319ba4SBarry Smith #undef __FUNCT__ 748b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 74913f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 750b4319ba4SBarry Smith { 751dfbe8321SBarry Smith PetscErrorCode ierr; 752b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 753b4319ba4SBarry Smith 754b4319ba4SBarry Smith PetscFunctionBegin; 755b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 756b4319ba4SBarry Smith PetscFunctionReturn(0); 757b4319ba4SBarry Smith } 758b4319ba4SBarry Smith 759b4319ba4SBarry Smith #undef __FUNCT__ 760f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 761f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 762f0006bf2SLisandro Dalcin { 763f0006bf2SLisandro Dalcin PetscErrorCode ierr; 764f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 765f0006bf2SLisandro Dalcin 766f0006bf2SLisandro Dalcin PetscFunctionBegin; 767f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 768f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 769f0006bf2SLisandro Dalcin } 770f0006bf2SLisandro Dalcin 771f0006bf2SLisandro Dalcin #undef __FUNCT__ 7722e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 7732b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 7742e74eeadSLisandro Dalcin { 7750298fd71SBarry Smith PetscInt n_l = 0, *rows_l = NULL; 7762e74eeadSLisandro Dalcin PetscErrorCode ierr; 7772e74eeadSLisandro Dalcin 7782e74eeadSLisandro Dalcin PetscFunctionBegin; 779ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 7802e74eeadSLisandro Dalcin if (n) { 781785e854fSJed Brown ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr); 7823bbff08aSStefano Zampini ierr = ISGlobalToLocalMappingApply(A->rmap->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 7832e74eeadSLisandro Dalcin } 7842b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 785c31cb41cSBarry Smith ierr = PetscFree(rows_l);CHKERRQ(ierr); 7862e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7872e74eeadSLisandro Dalcin } 7882e74eeadSLisandro Dalcin 7892e74eeadSLisandro Dalcin #undef __FUNCT__ 790b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 7912b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 792b4319ba4SBarry Smith { 793b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 794dfbe8321SBarry Smith PetscErrorCode ierr; 795f4df32b1SMatthew Knepley PetscInt i; 796b4319ba4SBarry Smith PetscScalar *array; 797b4319ba4SBarry Smith 798b4319ba4SBarry Smith PetscFunctionBegin; 799ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 800b4319ba4SBarry Smith { 801b4319ba4SBarry Smith /* 802b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 803b4319ba4SBarry Smith work properly in the interface nodes. 804b4319ba4SBarry Smith */ 805b4319ba4SBarry Smith Vec counter; 806e176bc59SStefano Zampini ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr); 807e176bc59SStefano Zampini ierr = VecSet(counter,0.);CHKERRQ(ierr); 808e176bc59SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 809e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 810e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 811e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 812e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 8136bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 814b4319ba4SBarry Smith } 815958c9bccSBarry Smith if (!n) { 816b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 817b4319ba4SBarry Smith } else { 818b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 8192205254eSKarl Rupp 820e176bc59SStefano Zampini ierr = VecGetArray(is->y,&array);CHKERRQ(ierr); 8212b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 822b4319ba4SBarry Smith for (i=0; i<n; i++) { 823f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 824b4319ba4SBarry Smith } 825b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 826b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 827e176bc59SStefano Zampini ierr = VecRestoreArray(is->y,&array);CHKERRQ(ierr); 828b4319ba4SBarry Smith } 829b4319ba4SBarry Smith PetscFunctionReturn(0); 830b4319ba4SBarry Smith } 831b4319ba4SBarry Smith 832b4319ba4SBarry Smith #undef __FUNCT__ 833b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 834dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 835b4319ba4SBarry Smith { 836b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 837dfbe8321SBarry Smith PetscErrorCode ierr; 838b4319ba4SBarry Smith 839b4319ba4SBarry Smith PetscFunctionBegin; 840b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 841b4319ba4SBarry Smith PetscFunctionReturn(0); 842b4319ba4SBarry Smith } 843b4319ba4SBarry Smith 844b4319ba4SBarry Smith #undef __FUNCT__ 845b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 846dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 847b4319ba4SBarry Smith { 848b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 849dfbe8321SBarry Smith PetscErrorCode ierr; 850b4319ba4SBarry Smith 851b4319ba4SBarry Smith PetscFunctionBegin; 852b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 853b4319ba4SBarry Smith PetscFunctionReturn(0); 854b4319ba4SBarry Smith } 855b4319ba4SBarry Smith 856b4319ba4SBarry Smith #undef __FUNCT__ 857b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 8587087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 859b4319ba4SBarry Smith { 860b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 861b4319ba4SBarry Smith 862b4319ba4SBarry Smith PetscFunctionBegin; 863b4319ba4SBarry Smith *local = is->A; 864b4319ba4SBarry Smith PetscFunctionReturn(0); 865b4319ba4SBarry Smith } 866b4319ba4SBarry Smith 867b4319ba4SBarry Smith #undef __FUNCT__ 868b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 869b4319ba4SBarry Smith /*@ 870b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 871b4319ba4SBarry Smith 872b4319ba4SBarry Smith Input Parameter: 873b4319ba4SBarry Smith . mat - the matrix 874b4319ba4SBarry Smith 875b4319ba4SBarry Smith Output Parameter: 876eb82efa4SStefano Zampini . local - the local matrix 877b4319ba4SBarry Smith 878b4319ba4SBarry Smith Level: advanced 879b4319ba4SBarry Smith 880b4319ba4SBarry Smith Notes: 881b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 882b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 883b4319ba4SBarry Smith of the MatSetValues() operation. 884b4319ba4SBarry Smith 885b4319ba4SBarry Smith .seealso: MATIS 886b4319ba4SBarry Smith @*/ 8877087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 888b4319ba4SBarry Smith { 8894ac538c5SBarry Smith PetscErrorCode ierr; 890b4319ba4SBarry Smith 891b4319ba4SBarry Smith PetscFunctionBegin; 8920700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 893b4319ba4SBarry Smith PetscValidPointer(local,2); 8944ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 895b4319ba4SBarry Smith PetscFunctionReturn(0); 896b4319ba4SBarry Smith } 897b4319ba4SBarry Smith 8983b03a366Sstefano_zampini #undef __FUNCT__ 8993b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 9003b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 9013b03a366Sstefano_zampini { 9023b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 9033b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 9043b03a366Sstefano_zampini PetscErrorCode ierr; 9053b03a366Sstefano_zampini 9063b03a366Sstefano_zampini PetscFunctionBegin; 9074e4c7dbeSStefano Zampini if (is->A) { 9083b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 9093b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 910f23aa3ddSBarry 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); 9114e4c7dbeSStefano Zampini } 9123b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 9133b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 9143b03a366Sstefano_zampini is->A = local; 9153b03a366Sstefano_zampini PetscFunctionReturn(0); 9163b03a366Sstefano_zampini } 9173b03a366Sstefano_zampini 9183b03a366Sstefano_zampini #undef __FUNCT__ 9193b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 9203b03a366Sstefano_zampini /*@ 921eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 9223b03a366Sstefano_zampini 9233b03a366Sstefano_zampini Input Parameter: 9243b03a366Sstefano_zampini . mat - the matrix 925eb82efa4SStefano Zampini . local - the local matrix 9263b03a366Sstefano_zampini 9273b03a366Sstefano_zampini Output Parameter: 9283b03a366Sstefano_zampini 9293b03a366Sstefano_zampini Level: advanced 9303b03a366Sstefano_zampini 9313b03a366Sstefano_zampini Notes: 9323b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 9333b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 9343b03a366Sstefano_zampini 9353b03a366Sstefano_zampini .seealso: MATIS 9363b03a366Sstefano_zampini @*/ 9373b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 9383b03a366Sstefano_zampini { 9393b03a366Sstefano_zampini PetscErrorCode ierr; 9403b03a366Sstefano_zampini 9413b03a366Sstefano_zampini PetscFunctionBegin; 9423b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 943b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 9443b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 9453b03a366Sstefano_zampini PetscFunctionReturn(0); 9463b03a366Sstefano_zampini } 9473b03a366Sstefano_zampini 9486726f965SBarry Smith #undef __FUNCT__ 9496726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 9506726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 9516726f965SBarry Smith { 9526726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 9536726f965SBarry Smith PetscErrorCode ierr; 9546726f965SBarry Smith 9556726f965SBarry Smith PetscFunctionBegin; 9566726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 9576726f965SBarry Smith PetscFunctionReturn(0); 9586726f965SBarry Smith } 9596726f965SBarry Smith 9606726f965SBarry Smith #undef __FUNCT__ 9612e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 9622e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 9632e74eeadSLisandro Dalcin { 9642e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9652e74eeadSLisandro Dalcin PetscErrorCode ierr; 9662e74eeadSLisandro Dalcin 9672e74eeadSLisandro Dalcin PetscFunctionBegin; 9682e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 9692e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9702e74eeadSLisandro Dalcin } 9712e74eeadSLisandro Dalcin 9722e74eeadSLisandro Dalcin #undef __FUNCT__ 9732e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 9742e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 9752e74eeadSLisandro Dalcin { 9762e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9772e74eeadSLisandro Dalcin PetscErrorCode ierr; 9782e74eeadSLisandro Dalcin 9792e74eeadSLisandro Dalcin PetscFunctionBegin; 9802e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 981e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 9822e74eeadSLisandro Dalcin 9832e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 9842e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 985e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 986e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9872e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9882e74eeadSLisandro Dalcin } 9892e74eeadSLisandro Dalcin 9902e74eeadSLisandro Dalcin #undef __FUNCT__ 9916726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 992ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 9936726f965SBarry Smith { 9946726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 9956726f965SBarry Smith PetscErrorCode ierr; 9966726f965SBarry Smith 9976726f965SBarry Smith PetscFunctionBegin; 9984e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 9996726f965SBarry Smith PetscFunctionReturn(0); 10006726f965SBarry Smith } 10016726f965SBarry Smith 1002284134d9SBarry Smith #undef __FUNCT__ 1003284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 1004284134d9SBarry Smith /*@ 1005284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 1006284134d9SBarry Smith process but not across processes. 1007284134d9SBarry Smith 1008284134d9SBarry Smith Input Parameters: 1009284134d9SBarry Smith + comm - MPI communicator that will share the matrix 1010e176bc59SStefano Zampini . bs - block size of the matrix 1011df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 1012e176bc59SStefano Zampini . rmap - local to global map for rows 1013e176bc59SStefano Zampini - cmap - local to global map for cols 1014284134d9SBarry Smith 1015284134d9SBarry Smith Output Parameter: 1016284134d9SBarry Smith . A - the resulting matrix 1017284134d9SBarry Smith 10188e6c10adSSatish Balay Level: advanced 10198e6c10adSSatish Balay 1020284134d9SBarry Smith Notes: See MATIS for more details 10218cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 1022e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 1023e176bc59SStefano Zampini If either rmap or cmap are NULL, than the matrix is assumed to be square 1024284134d9SBarry Smith 1025284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 1026284134d9SBarry Smith @*/ 1027e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 1028284134d9SBarry Smith { 1029284134d9SBarry Smith PetscErrorCode ierr; 1030284134d9SBarry Smith 1031284134d9SBarry Smith PetscFunctionBegin; 1032e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 1033284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1034d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 1035284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1036284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 10373b03a366Sstefano_zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 1038e176bc59SStefano Zampini if (rmap && cmap) { 1039e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 1040e176bc59SStefano Zampini } else if (!rmap) { 1041e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 1042e176bc59SStefano Zampini } else { 1043e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 1044e176bc59SStefano Zampini } 1045284134d9SBarry Smith PetscFunctionReturn(0); 1046284134d9SBarry Smith } 1047284134d9SBarry Smith 1048b4319ba4SBarry Smith /*MC 1049eb82efa4SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC). 1050b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 1051b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 1052b4319ba4SBarry Smith product is handled "implicitly". 1053b4319ba4SBarry Smith 1054b4319ba4SBarry Smith Operations Provided: 10556726f965SBarry Smith + MatMult() 10562e74eeadSLisandro Dalcin . MatMultAdd() 10572e74eeadSLisandro Dalcin . MatMultTranspose() 10582e74eeadSLisandro Dalcin . MatMultTransposeAdd() 10596726f965SBarry Smith . MatZeroEntries() 10606726f965SBarry Smith . MatSetOption() 10612e74eeadSLisandro Dalcin . MatZeroRows() 10626726f965SBarry Smith . MatZeroRowsLocal() 10632e74eeadSLisandro Dalcin . MatSetValues() 10646726f965SBarry Smith . MatSetValuesLocal() 10652e74eeadSLisandro Dalcin . MatScale() 10662e74eeadSLisandro Dalcin . MatGetDiagonal() 10676726f965SBarry Smith - MatSetLocalToGlobalMapping() 1068b4319ba4SBarry Smith 1069b4319ba4SBarry Smith Options Database Keys: 1070b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1071b4319ba4SBarry Smith 1072b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1073b4319ba4SBarry Smith 1074b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1075b4319ba4SBarry Smith 1076b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1077eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1078b4319ba4SBarry Smith 1079b4319ba4SBarry Smith Level: advanced 1080b4319ba4SBarry Smith 1081eb82efa4SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC 1082b4319ba4SBarry Smith 1083b4319ba4SBarry Smith M*/ 1084b4319ba4SBarry Smith 1085b4319ba4SBarry Smith #undef __FUNCT__ 1086b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 10878cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1088b4319ba4SBarry Smith { 1089dfbe8321SBarry Smith PetscErrorCode ierr; 1090b4319ba4SBarry Smith Mat_IS *b; 1091b4319ba4SBarry Smith 1092b4319ba4SBarry Smith PetscFunctionBegin; 1093b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1094b4319ba4SBarry Smith A->data = (void*)b; 1095b4319ba4SBarry Smith 1096e176bc59SStefano Zampini /* matrix ops */ 1097e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1098b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 10992e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 11002e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 11012e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1102b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1103b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 11042e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 1105b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1106f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 11072e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1108b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 1109b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1110b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1111b4319ba4SBarry Smith A->ops->view = MatView_IS; 11126726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 11132e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 11142e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 11156726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 111669796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 111769796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1118ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 1119b4319ba4SBarry Smith 1120b7ce53b6SStefano Zampini /* special MATIS functions */ 1121bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1122bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1123b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 11242e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 1125*cf0a3239SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",MatISSetUpSF_IS);CHKERRQ(ierr); 112617667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1127b4319ba4SBarry Smith PetscFunctionReturn(0); 1128b4319ba4SBarry Smith } 1129