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*/ 164f2d7cafSStefano Zampini #include <petsc/private/sfimpl.h> 1728f4e0baSStefano Zampini 18cf0a3239SStefano Zampini #undef __FUNCT__ 19cf0a3239SStefano Zampini #define __FUNCT__ "MatISSetUpSF" 20cf0a3239SStefano Zampini /*@ 21cf0a3239SStefano Zampini MatISSetUpSF - Setup star forest object used by MatIS during MatISSetPreallocation. 22cf0a3239SStefano Zampini 23cf0a3239SStefano Zampini Collective on MPI_Comm 24cf0a3239SStefano Zampini 25cf0a3239SStefano Zampini Input Parameters: 26cf0a3239SStefano Zampini + A - the matrix 27cf0a3239SStefano Zampini 28cf0a3239SStefano Zampini Level: advanced 29cf0a3239SStefano Zampini 30cf0a3239SStefano Zampini Notes: This function does not be to be called by the user. 31cf0a3239SStefano Zampini 32cf0a3239SStefano Zampini .keywords: matrix 33cf0a3239SStefano Zampini 34cf0a3239SStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatISSetPreallocation(), MatISGetLocalMat() 35cf0a3239SStefano Zampini @*/ 36cf0a3239SStefano Zampini PetscErrorCode MatISSetUpSF(Mat A) 37cf0a3239SStefano Zampini { 38cf0a3239SStefano Zampini PetscErrorCode ierr; 39cf0a3239SStefano Zampini 40cf0a3239SStefano Zampini PetscFunctionBegin; 41cf0a3239SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 42cf0a3239SStefano Zampini PetscValidType(A,1); 43cf0a3239SStefano Zampini ierr = PetscTryMethod(A,"MatISSetUpSF_C",(Mat),(A));CHKERRQ(ierr); 44cf0a3239SStefano Zampini PetscFunctionReturn(0); 45cf0a3239SStefano Zampini } 46a72627d2SStefano Zampini 4728f4e0baSStefano Zampini #undef __FUNCT__ 48cf0a3239SStefano Zampini #define __FUNCT__ "MatISSetUpSF_IS" 49cf0a3239SStefano Zampini static PetscErrorCode MatISSetUpSF_IS(Mat B) 5028f4e0baSStefano Zampini { 5128f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 5228f4e0baSStefano Zampini const PetscInt *gidxs; 534f2d7cafSStefano Zampini PetscInt nleaves; 5428f4e0baSStefano Zampini PetscErrorCode ierr; 5528f4e0baSStefano Zampini 5628f4e0baSStefano Zampini PetscFunctionBegin; 574f2d7cafSStefano Zampini if (matis->sf) PetscFunctionReturn(0); 5828f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 593bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 604f2d7cafSStefano Zampini ierr = ISLocalToGlobalMappingGetSize(B->rmap->mapping,&nleaves);CHKERRQ(ierr); 6128f4e0baSStefano Zampini /* PETSC_OWN_POINTER refers to ilocal which is NULL */ 624f2d7cafSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 633bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 644f2d7cafSStefano Zampini ierr = PetscMalloc2(matis->sf->nroots,&matis->sf_rootdata,matis->sf->nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 6528f4e0baSStefano Zampini PetscFunctionReturn(0); 6628f4e0baSStefano Zampini } 672e1947a5SStefano Zampini 682e1947a5SStefano Zampini #undef __FUNCT__ 692e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 70eb82efa4SStefano Zampini /*@ 71a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 72a88811baSStefano Zampini 73a88811baSStefano Zampini Collective on MPI_Comm 74a88811baSStefano Zampini 75a88811baSStefano Zampini Input Parameters: 76a88811baSStefano Zampini + B - the matrix 77a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 78a88811baSStefano Zampini (same value is used for all local rows) 79a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 80a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 81a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 82a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 83a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 84a88811baSStefano Zampini the diagonal entry even if it is zero. 85a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 86a88811baSStefano Zampini submatrix (same value is used for all local rows). 87a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 88a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 89a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 90a88811baSStefano Zampini structure. The size of this array is equal to the number 91a88811baSStefano Zampini of local rows, i.e 'm'. 92a88811baSStefano Zampini 93a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 94a88811baSStefano Zampini 95a88811baSStefano Zampini Level: intermediate 96a88811baSStefano Zampini 97a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 98a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 99a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 100a88811baSStefano Zampini 101a88811baSStefano Zampini .keywords: matrix 102a88811baSStefano Zampini 103a88811baSStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat() 104a88811baSStefano Zampini @*/ 1052e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 1062e1947a5SStefano Zampini { 1072e1947a5SStefano Zampini PetscErrorCode ierr; 1082e1947a5SStefano Zampini 1092e1947a5SStefano Zampini PetscFunctionBegin; 1102e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 1112e1947a5SStefano Zampini PetscValidType(B,1); 1122e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 1132e1947a5SStefano Zampini PetscFunctionReturn(0); 1142e1947a5SStefano Zampini } 1152e1947a5SStefano Zampini 1162e1947a5SStefano Zampini #undef __FUNCT__ 1172e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 1182e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 1192e1947a5SStefano Zampini { 1202e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 12128f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 1222e1947a5SStefano Zampini PetscErrorCode ierr; 1232e1947a5SStefano Zampini 1242e1947a5SStefano Zampini PetscFunctionBegin; 1256c4ed002SBarry Smith if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 126cf0a3239SStefano Zampini ierr = MatISSetUpSF(B);CHKERRQ(ierr); 1274f2d7cafSStefano Zampini 1284f2d7cafSStefano Zampini if (!d_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nz; 1294f2d7cafSStefano Zampini else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 1304f2d7cafSStefano Zampini 1314f2d7cafSStefano Zampini if (!o_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nz; 1324f2d7cafSStefano Zampini else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 1334f2d7cafSStefano Zampini 13428f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 13528f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 13628f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 13728f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 1384f2d7cafSStefano Zampini 1394f2d7cafSStefano Zampini for (i=0;i<matis->sf->nleaves;i++) matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 14028f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 1414f2d7cafSStefano Zampini 1424f2d7cafSStefano Zampini for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 14328f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 1444f2d7cafSStefano Zampini 1454f2d7cafSStefano Zampini for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 14628f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 1472e1947a5SStefano Zampini PetscFunctionReturn(0); 1482e1947a5SStefano Zampini } 149b4319ba4SBarry Smith 150b4319ba4SBarry Smith #undef __FUNCT__ 1513927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private" 1523927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 1533927de2eSStefano Zampini { 1543927de2eSStefano Zampini Mat_IS *matis = (Mat_IS*)(A->data); 1553927de2eSStefano Zampini PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 156ecf5a873SStefano Zampini const PetscInt *global_indices_r,*global_indices_c; 1573927de2eSStefano Zampini PetscInt i,j,bs,rows,cols; 1583927de2eSStefano Zampini PetscInt lrows,lcols; 1593927de2eSStefano Zampini PetscInt local_rows,local_cols; 1603927de2eSStefano Zampini PetscMPIInt nsubdomains; 1613927de2eSStefano Zampini PetscBool isdense,issbaij; 1623927de2eSStefano Zampini PetscErrorCode ierr; 1633927de2eSStefano Zampini 1643927de2eSStefano Zampini PetscFunctionBegin; 1653927de2eSStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr); 1663927de2eSStefano Zampini ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr); 1673927de2eSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1683927de2eSStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 1693927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 1703927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 171ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 172ecf5a873SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 173ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_c);CHKERRQ(ierr); 174ecf5a873SStefano Zampini } else { 175ecf5a873SStefano Zampini global_indices_c = global_indices_r; 176ecf5a873SStefano Zampini } 177ecf5a873SStefano Zampini 1783927de2eSStefano Zampini if (issbaij) { 1793927de2eSStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 1803927de2eSStefano Zampini } 1813927de2eSStefano Zampini /* 182ecf5a873SStefano Zampini An SF reduce is needed to sum up properly on shared rows. 1833927de2eSStefano Zampini Note that generally preallocation is not exact, since it overestimates nonzeros 1843927de2eSStefano Zampini */ 185cf0a3239SStefano Zampini ierr = MatISSetUpSF(A);CHKERRQ(ierr); 1863927de2eSStefano Zampini ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr); 1873927de2eSStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr); 1883927de2eSStefano Zampini /* All processes need to compute entire row ownership */ 1893927de2eSStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 1903927de2eSStefano Zampini ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 1913927de2eSStefano Zampini for (i=0;i<nsubdomains;i++) { 1923927de2eSStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 1933927de2eSStefano Zampini row_ownership[j] = i; 1943927de2eSStefano Zampini } 1953927de2eSStefano Zampini } 1963927de2eSStefano Zampini 1973927de2eSStefano Zampini /* 1983927de2eSStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 1993927de2eSStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 2003927de2eSStefano Zampini */ 2013927de2eSStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 2023927de2eSStefano Zampini /* preallocation as a MATAIJ */ 2033927de2eSStefano Zampini if (isdense) { /* special case for dense local matrices */ 2043927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 205ecf5a873SStefano Zampini PetscInt index_row = global_indices_r[i]; 2063927de2eSStefano Zampini for (j=i;j<local_rows;j++) { 2073927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 208ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[j]; 2093927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 2103927de2eSStefano Zampini my_dnz[i] += 1; 2113927de2eSStefano Zampini } else { /* offdiag block */ 2123927de2eSStefano Zampini my_onz[i] += 1; 2133927de2eSStefano Zampini } 2143927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 2153927de2eSStefano Zampini if (i != j) { 2163927de2eSStefano Zampini owner = row_ownership[index_col]; 2173927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 2183927de2eSStefano Zampini my_dnz[j] += 1; 2193927de2eSStefano Zampini } else { 2203927de2eSStefano Zampini my_onz[j] += 1; 2213927de2eSStefano Zampini } 2223927de2eSStefano Zampini } 2233927de2eSStefano Zampini } 2243927de2eSStefano Zampini } 225ecf5a873SStefano Zampini } else { /* TODO: this could be optimized using MatGetRowIJ */ 2263927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 2273927de2eSStefano Zampini const PetscInt *cols; 228ecf5a873SStefano Zampini PetscInt ncols,index_row = global_indices_r[i]; 2293927de2eSStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 2303927de2eSStefano Zampini for (j=0;j<ncols;j++) { 2313927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 232ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[cols[j]]; 2333927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 2343927de2eSStefano Zampini my_dnz[i] += 1; 2353927de2eSStefano Zampini } else { /* offdiag block */ 2363927de2eSStefano Zampini my_onz[i] += 1; 2373927de2eSStefano Zampini } 2383927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 239d9a9e74cSStefano Zampini if (issbaij && index_col != index_row) { 2403927de2eSStefano Zampini owner = row_ownership[index_col]; 2413927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 242d9a9e74cSStefano Zampini my_dnz[cols[j]] += 1; 2433927de2eSStefano Zampini } else { 244d9a9e74cSStefano Zampini my_onz[cols[j]] += 1; 2453927de2eSStefano Zampini } 2463927de2eSStefano Zampini } 2473927de2eSStefano Zampini } 2483927de2eSStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 2493927de2eSStefano Zampini } 2503927de2eSStefano Zampini } 251ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 252ecf5a873SStefano Zampini if (global_indices_c != global_indices_r) { 253ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_c);CHKERRQ(ierr); 254ecf5a873SStefano Zampini } 2553927de2eSStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 256ecf5a873SStefano Zampini 257ecf5a873SStefano Zampini /* Reduce my_dnz and my_onz */ 2583927de2eSStefano Zampini if (maxreduce) { 2593927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 2603927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 2613927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 2623927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 2633927de2eSStefano Zampini } else { 2643927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 2653927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 2663927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 2673927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 2683927de2eSStefano Zampini } 2693927de2eSStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 2703927de2eSStefano Zampini 2713927de2eSStefano Zampini /* Resize preallocation if overestimated */ 2723927de2eSStefano Zampini for (i=0;i<lrows;i++) { 2733927de2eSStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 2743927de2eSStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 2753927de2eSStefano Zampini } 2763927de2eSStefano Zampini /* set preallocation */ 2773927de2eSStefano Zampini ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 2783927de2eSStefano Zampini for (i=0;i<lrows/bs;i++) { 2793927de2eSStefano Zampini dnz[i] = dnz[i*bs]/bs; 2803927de2eSStefano Zampini onz[i] = onz[i*bs]/bs; 2813927de2eSStefano Zampini } 2823927de2eSStefano Zampini ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 2833927de2eSStefano Zampini ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 2843927de2eSStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 2853927de2eSStefano Zampini if (issbaij) { 2863927de2eSStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 2873927de2eSStefano Zampini } 2883927de2eSStefano Zampini PetscFunctionReturn(0); 2893927de2eSStefano Zampini } 2903927de2eSStefano Zampini 2913927de2eSStefano Zampini #undef __FUNCT__ 292b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 293b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 294b7ce53b6SStefano Zampini { 295b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 296d9a9e74cSStefano Zampini Mat local_mat; 297b7ce53b6SStefano Zampini /* info on mat */ 2983cfa4ea4SStefano Zampini PetscInt bs,rows,cols,lrows,lcols; 299b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 300686e3a49SStefano Zampini PetscBool isdense,issbaij,isseqaij; 3017c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 302b7ce53b6SStefano Zampini /* values insertion */ 303b7ce53b6SStefano Zampini PetscScalar *array; 304b7ce53b6SStefano Zampini /* work */ 305b7ce53b6SStefano Zampini PetscErrorCode ierr; 306b7ce53b6SStefano Zampini 307b7ce53b6SStefano Zampini PetscFunctionBegin; 308b7ce53b6SStefano Zampini /* get info from mat */ 3097c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 3107c03b4e8SStefano Zampini if (nsubdomains == 1) { 3117c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 3127c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 3137c03b4e8SStefano Zampini } else { 3147c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 3157c03b4e8SStefano Zampini } 3167c03b4e8SStefano Zampini PetscFunctionReturn(0); 3177c03b4e8SStefano Zampini } 318b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 319b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 3203cfa4ea4SStefano Zampini ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 321b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 322b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 323686e3a49SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 324b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 325b7ce53b6SStefano Zampini 326b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 327b7ce53b6SStefano Zampini MatType new_mat_type; 3283927de2eSStefano Zampini PetscBool issbaij_red; 329b7ce53b6SStefano Zampini 330b7ce53b6SStefano Zampini /* determining new matrix type */ 331b2566f29SBarry Smith ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 332b7ce53b6SStefano Zampini if (issbaij_red) { 333b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 334b7ce53b6SStefano Zampini } else { 335b7ce53b6SStefano Zampini if (bs>1) { 336b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 337b7ce53b6SStefano Zampini } else { 338b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 339b7ce53b6SStefano Zampini } 340b7ce53b6SStefano Zampini } 341b7ce53b6SStefano Zampini 3423927de2eSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 3433cfa4ea4SStefano Zampini ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 3443927de2eSStefano Zampini ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 3453927de2eSStefano Zampini ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr); 3463927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 347b7ce53b6SStefano Zampini } else { 3483cfa4ea4SStefano Zampini PetscInt mbs,mrows,mcols,mlrows,mlcols; 349b7ce53b6SStefano Zampini /* some checks */ 350b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 351b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 3523cfa4ea4SStefano Zampini ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 3536c4ed002SBarry Smith if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 3546c4ed002SBarry Smith if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 3556c4ed002SBarry Smith if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 3566c4ed002SBarry Smith if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 3576c4ed002SBarry Smith if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 358b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 359b7ce53b6SStefano Zampini } 360d9a9e74cSStefano Zampini 361d9a9e74cSStefano Zampini if (issbaij) { 362d9a9e74cSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 363d9a9e74cSStefano Zampini } else { 364d9a9e74cSStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 365d9a9e74cSStefano Zampini local_mat = matis->A; 366d9a9e74cSStefano Zampini } 367686e3a49SStefano Zampini 368b7ce53b6SStefano Zampini /* Set values */ 369ecf5a873SStefano Zampini ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 370b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 371ecf5a873SStefano Zampini PetscInt i,*dummy_rows,*dummy_cols; 372ecf5a873SStefano Zampini 373ecf5a873SStefano Zampini if (local_rows != local_cols) { 374ecf5a873SStefano Zampini ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 375ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 376ecf5a873SStefano Zampini for (i=0;i<local_cols;i++) dummy_cols[i] = i; 377ecf5a873SStefano Zampini } else { 378ecf5a873SStefano Zampini ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 379ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 380ecf5a873SStefano Zampini dummy_cols = dummy_rows; 381ecf5a873SStefano Zampini } 382b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 383d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 384ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 385d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 386ecf5a873SStefano Zampini if (dummy_rows != dummy_cols) { 387ecf5a873SStefano Zampini ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 388ecf5a873SStefano Zampini } else { 389ecf5a873SStefano Zampini ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 390ecf5a873SStefano Zampini } 391686e3a49SStefano Zampini } else if (isseqaij) { 392ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 393686e3a49SStefano Zampini PetscBool done; 394686e3a49SStefano Zampini 395d9a9e74cSStefano Zampini ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 3966c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 397d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 398686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 399ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 400686e3a49SStefano Zampini } 401d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 4026c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 403d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 404686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 405ecf5a873SStefano Zampini PetscInt i; 406c0962df8SStefano Zampini 407686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 408686e3a49SStefano Zampini PetscInt j; 409ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 410686e3a49SStefano Zampini 411ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 412ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 413ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 414686e3a49SStefano Zampini } 415b7ce53b6SStefano Zampini } 416b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 417d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 418b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 419b7ce53b6SStefano Zampini if (isdense) { 420b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 421b7ce53b6SStefano Zampini } 422b7ce53b6SStefano Zampini PetscFunctionReturn(0); 423b7ce53b6SStefano Zampini } 424b7ce53b6SStefano Zampini 425b7ce53b6SStefano Zampini #undef __FUNCT__ 426b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 427b7ce53b6SStefano Zampini /*@ 428b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 429b7ce53b6SStefano Zampini 430b7ce53b6SStefano Zampini Input Parameter: 431b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 432b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 433b7ce53b6SStefano Zampini 434b7ce53b6SStefano Zampini Output Parameter: 435b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 436b7ce53b6SStefano Zampini 437b7ce53b6SStefano Zampini Level: developer 438b7ce53b6SStefano Zampini 439eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 440b7ce53b6SStefano Zampini 441b7ce53b6SStefano Zampini .seealso: MATIS 442b7ce53b6SStefano Zampini @*/ 443b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 444b7ce53b6SStefano Zampini { 445b7ce53b6SStefano Zampini PetscErrorCode ierr; 446b7ce53b6SStefano Zampini 447b7ce53b6SStefano Zampini PetscFunctionBegin; 448b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 449b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 450b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 451b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 452b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 453b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 4546c4ed002SBarry Smith if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 455b7ce53b6SStefano Zampini } 456b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 457b7ce53b6SStefano Zampini PetscFunctionReturn(0); 458b7ce53b6SStefano Zampini } 459b7ce53b6SStefano Zampini 460b7ce53b6SStefano Zampini #undef __FUNCT__ 461ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 462ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 463ad6194a2SStefano Zampini { 464ad6194a2SStefano Zampini PetscErrorCode ierr; 465ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 466ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 467ad6194a2SStefano Zampini Mat B,localmat; 468ad6194a2SStefano Zampini 469ad6194a2SStefano Zampini PetscFunctionBegin; 470ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 471ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 472ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 473e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 474ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 475ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 476b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 477ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 478ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 479ad6194a2SStefano Zampini *newmat = B; 480ad6194a2SStefano Zampini PetscFunctionReturn(0); 481ad6194a2SStefano Zampini } 482ad6194a2SStefano Zampini 483ad6194a2SStefano Zampini #undef __FUNCT__ 48469796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 48569796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 48669796d55SStefano Zampini { 48769796d55SStefano Zampini PetscErrorCode ierr; 48869796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 48969796d55SStefano Zampini PetscBool local_sym; 49069796d55SStefano Zampini 49169796d55SStefano Zampini PetscFunctionBegin; 49269796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 493b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 49469796d55SStefano Zampini PetscFunctionReturn(0); 49569796d55SStefano Zampini } 49669796d55SStefano Zampini 49769796d55SStefano Zampini #undef __FUNCT__ 49869796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 49969796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 50069796d55SStefano Zampini { 50169796d55SStefano Zampini PetscErrorCode ierr; 50269796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 50369796d55SStefano Zampini PetscBool local_sym; 50469796d55SStefano Zampini 50569796d55SStefano Zampini PetscFunctionBegin; 50669796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 507b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 50869796d55SStefano Zampini PetscFunctionReturn(0); 50969796d55SStefano Zampini } 51069796d55SStefano Zampini 51169796d55SStefano Zampini #undef __FUNCT__ 512b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 513dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 514b4319ba4SBarry Smith { 515dfbe8321SBarry Smith PetscErrorCode ierr; 516b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 517b4319ba4SBarry Smith 518b4319ba4SBarry Smith PetscFunctionBegin; 5196bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 520e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 521e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 5226bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 5236bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 52428f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 52528f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 526bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 527dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 528bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 529b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 530b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 5312e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 532cf0a3239SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",NULL);CHKERRQ(ierr); 533b4319ba4SBarry Smith PetscFunctionReturn(0); 534b4319ba4SBarry Smith } 535b4319ba4SBarry Smith 536b4319ba4SBarry Smith #undef __FUNCT__ 537b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 538dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 539b4319ba4SBarry Smith { 540dfbe8321SBarry Smith PetscErrorCode ierr; 541b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 542b4319ba4SBarry Smith PetscScalar zero = 0.0; 543b4319ba4SBarry Smith 544b4319ba4SBarry Smith PetscFunctionBegin; 545b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 546e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 547e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 548b4319ba4SBarry Smith 549b4319ba4SBarry Smith /* multiply the local matrix */ 550b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 551b4319ba4SBarry Smith 552b4319ba4SBarry Smith /* scatter product back into global memory */ 5532dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 554e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 555e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 556b4319ba4SBarry Smith PetscFunctionReturn(0); 557b4319ba4SBarry Smith } 558b4319ba4SBarry Smith 559b4319ba4SBarry Smith #undef __FUNCT__ 5602e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 561b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5622e74eeadSLisandro Dalcin { 563650997f4SStefano Zampini Vec temp_vec; 5642e74eeadSLisandro Dalcin PetscErrorCode ierr; 5652e74eeadSLisandro Dalcin 5662e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 567650997f4SStefano Zampini if (v3 != v2) { 568650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 569650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 570650997f4SStefano Zampini } else { 571650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 572650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 573650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 574650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 575650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 576650997f4SStefano Zampini } 5772e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5782e74eeadSLisandro Dalcin } 5792e74eeadSLisandro Dalcin 5802e74eeadSLisandro Dalcin #undef __FUNCT__ 5812e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 582e176bc59SStefano Zampini PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 5832e74eeadSLisandro Dalcin { 5842e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 5852e74eeadSLisandro Dalcin PetscErrorCode ierr; 5862e74eeadSLisandro Dalcin 587e176bc59SStefano Zampini PetscFunctionBegin; 5882e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 589e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 590e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5912e74eeadSLisandro Dalcin 5922e74eeadSLisandro Dalcin /* multiply the local matrix */ 593e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 5942e74eeadSLisandro Dalcin 5952e74eeadSLisandro Dalcin /* scatter product back into global vector */ 596e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 597e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 598e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5992e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6002e74eeadSLisandro Dalcin } 6012e74eeadSLisandro Dalcin 6022e74eeadSLisandro Dalcin #undef __FUNCT__ 6032e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 6042e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 6052e74eeadSLisandro Dalcin { 606650997f4SStefano Zampini Vec temp_vec; 6072e74eeadSLisandro Dalcin PetscErrorCode ierr; 6082e74eeadSLisandro Dalcin 6092e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 610650997f4SStefano Zampini if (v3 != v2) { 611650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 612650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 613650997f4SStefano Zampini } else { 614650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 615650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 616650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 617650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 618650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 619650997f4SStefano Zampini } 6202e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6212e74eeadSLisandro Dalcin } 6222e74eeadSLisandro Dalcin 6232e74eeadSLisandro Dalcin #undef __FUNCT__ 624b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 625dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 626b4319ba4SBarry Smith { 627b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 628dfbe8321SBarry Smith PetscErrorCode ierr; 629b4319ba4SBarry Smith PetscViewer sviewer; 630ee2491ecSStefano Zampini PetscBool isascii,view = PETSC_TRUE; 631b4319ba4SBarry Smith 632b4319ba4SBarry Smith PetscFunctionBegin; 633ee2491ecSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 634ee2491ecSStefano Zampini if (isascii) { 635ee2491ecSStefano Zampini PetscViewerFormat format; 636ee2491ecSStefano Zampini 637ee2491ecSStefano Zampini ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 638ee2491ecSStefano Zampini if (format == PETSC_VIEWER_ASCII_INFO) view = PETSC_FALSE; 639ee2491ecSStefano Zampini } 640ee2491ecSStefano Zampini if (!view) PetscFunctionReturn(0); 6413f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 642b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 6433f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 644b4319ba4SBarry Smith PetscFunctionReturn(0); 645b4319ba4SBarry Smith } 646b4319ba4SBarry Smith 647b4319ba4SBarry Smith #undef __FUNCT__ 648b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 649784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 650b4319ba4SBarry Smith { 651dfbe8321SBarry Smith PetscErrorCode ierr; 652e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 653b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 654b4319ba4SBarry Smith IS from,to; 655e176bc59SStefano Zampini Vec cglobal,rglobal; 656b4319ba4SBarry Smith 657b4319ba4SBarry Smith PetscFunctionBegin; 658784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 659e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 6603bbff08aSStefano Zampini /* Destroy any previous data */ 66170cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 66270cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 663e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 664e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 6651c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 66628f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 66728f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 6683bbff08aSStefano Zampini 6693bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 670fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 671fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 672fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 673fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 674b4319ba4SBarry Smith 675b4319ba4SBarry Smith /* Create the local matrix A */ 676e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 677e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 678e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 679e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 680f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 681e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 682e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 683e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 684ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 685ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 686b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 687b4319ba4SBarry Smith 688b4319ba4SBarry Smith /* Create the local work vectors */ 6893bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 690b4319ba4SBarry Smith 691e176bc59SStefano Zampini /* setup the global to local scatters */ 692e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 693e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 694784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 695e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 696e176bc59SStefano Zampini if (rmapping != cmapping) { 697e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 698e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 699e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 700e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 701e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 702e176bc59SStefano Zampini } else { 703e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 704e176bc59SStefano Zampini is->cctx = is->rctx; 705e176bc59SStefano Zampini } 706e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 707e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 7086bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 7096bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 710*9c0b00dbSStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 711b4319ba4SBarry Smith PetscFunctionReturn(0); 712b4319ba4SBarry Smith } 713b4319ba4SBarry Smith 7142e74eeadSLisandro Dalcin #undef __FUNCT__ 7152e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 7162e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 7172e74eeadSLisandro Dalcin { 7182e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 7192e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 7202e74eeadSLisandro Dalcin PetscErrorCode ierr; 7212e74eeadSLisandro Dalcin 7222e74eeadSLisandro Dalcin PetscFunctionBegin; 7232e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 724f23aa3ddSBarry 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); 7252e74eeadSLisandro Dalcin #endif 7263bbff08aSStefano Zampini ierr = ISG2LMapApply(mat->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 7273bbff08aSStefano Zampini ierr = ISG2LMapApply(mat->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 7282e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 7292e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7302e74eeadSLisandro Dalcin } 7312e74eeadSLisandro Dalcin 732b4319ba4SBarry Smith #undef __FUNCT__ 733b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 73413f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 735b4319ba4SBarry Smith { 736dfbe8321SBarry Smith PetscErrorCode ierr; 737b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 738b4319ba4SBarry Smith 739b4319ba4SBarry Smith PetscFunctionBegin; 740b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 741b4319ba4SBarry Smith PetscFunctionReturn(0); 742b4319ba4SBarry Smith } 743b4319ba4SBarry Smith 744b4319ba4SBarry Smith #undef __FUNCT__ 745f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 746f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 747f0006bf2SLisandro Dalcin { 748f0006bf2SLisandro Dalcin PetscErrorCode ierr; 749f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 750f0006bf2SLisandro Dalcin 751f0006bf2SLisandro Dalcin PetscFunctionBegin; 752f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 753f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 754f0006bf2SLisandro Dalcin } 755f0006bf2SLisandro Dalcin 756f0006bf2SLisandro Dalcin #undef __FUNCT__ 7572e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 7582b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 7592e74eeadSLisandro Dalcin { 7600298fd71SBarry Smith PetscInt n_l = 0, *rows_l = NULL; 7612e74eeadSLisandro Dalcin PetscErrorCode ierr; 7622e74eeadSLisandro Dalcin 7632e74eeadSLisandro Dalcin PetscFunctionBegin; 764ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 7652e74eeadSLisandro Dalcin if (n) { 766785e854fSJed Brown ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr); 7673bbff08aSStefano Zampini ierr = ISGlobalToLocalMappingApply(A->rmap->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 7682e74eeadSLisandro Dalcin } 7692b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 770c31cb41cSBarry Smith ierr = PetscFree(rows_l);CHKERRQ(ierr); 7712e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7722e74eeadSLisandro Dalcin } 7732e74eeadSLisandro Dalcin 7742e74eeadSLisandro Dalcin #undef __FUNCT__ 775b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 7762b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 777b4319ba4SBarry Smith { 778b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 779dfbe8321SBarry Smith PetscErrorCode ierr; 780f4df32b1SMatthew Knepley PetscInt i; 781b4319ba4SBarry Smith PetscScalar *array; 782b4319ba4SBarry Smith 783b4319ba4SBarry Smith PetscFunctionBegin; 784ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 785b4319ba4SBarry Smith { 786b4319ba4SBarry Smith /* 787b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 788b4319ba4SBarry Smith work properly in the interface nodes. 789b4319ba4SBarry Smith */ 790b4319ba4SBarry Smith Vec counter; 791e176bc59SStefano Zampini ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr); 792e176bc59SStefano Zampini ierr = VecSet(counter,0.);CHKERRQ(ierr); 793e176bc59SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 794e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 795e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 796e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 797e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7986bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 799b4319ba4SBarry Smith } 800958c9bccSBarry Smith if (!n) { 801b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 802b4319ba4SBarry Smith } else { 803b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 8042205254eSKarl Rupp 805e176bc59SStefano Zampini ierr = VecGetArray(is->y,&array);CHKERRQ(ierr); 8062b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 807b4319ba4SBarry Smith for (i=0; i<n; i++) { 808f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 809b4319ba4SBarry Smith } 810b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 811b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 812e176bc59SStefano Zampini ierr = VecRestoreArray(is->y,&array);CHKERRQ(ierr); 813b4319ba4SBarry Smith } 814b4319ba4SBarry Smith PetscFunctionReturn(0); 815b4319ba4SBarry Smith } 816b4319ba4SBarry Smith 817b4319ba4SBarry Smith #undef __FUNCT__ 818b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 819dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 820b4319ba4SBarry Smith { 821b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 822dfbe8321SBarry Smith PetscErrorCode ierr; 823b4319ba4SBarry Smith 824b4319ba4SBarry Smith PetscFunctionBegin; 825b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 826b4319ba4SBarry Smith PetscFunctionReturn(0); 827b4319ba4SBarry Smith } 828b4319ba4SBarry Smith 829b4319ba4SBarry Smith #undef __FUNCT__ 830b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 831dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 832b4319ba4SBarry Smith { 833b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 834dfbe8321SBarry Smith PetscErrorCode ierr; 835b4319ba4SBarry Smith 836b4319ba4SBarry Smith PetscFunctionBegin; 837b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 838b4319ba4SBarry Smith PetscFunctionReturn(0); 839b4319ba4SBarry Smith } 840b4319ba4SBarry Smith 841b4319ba4SBarry Smith #undef __FUNCT__ 842b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 8437087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 844b4319ba4SBarry Smith { 845b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 846b4319ba4SBarry Smith 847b4319ba4SBarry Smith PetscFunctionBegin; 848b4319ba4SBarry Smith *local = is->A; 849b4319ba4SBarry Smith PetscFunctionReturn(0); 850b4319ba4SBarry Smith } 851b4319ba4SBarry Smith 852b4319ba4SBarry Smith #undef __FUNCT__ 853b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 854b4319ba4SBarry Smith /*@ 855b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 856b4319ba4SBarry Smith 857b4319ba4SBarry Smith Input Parameter: 858b4319ba4SBarry Smith . mat - the matrix 859b4319ba4SBarry Smith 860b4319ba4SBarry Smith Output Parameter: 861eb82efa4SStefano Zampini . local - the local matrix 862b4319ba4SBarry Smith 863b4319ba4SBarry Smith Level: advanced 864b4319ba4SBarry Smith 865b4319ba4SBarry Smith Notes: 866b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 867b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 868b4319ba4SBarry Smith of the MatSetValues() operation. 869b4319ba4SBarry Smith 870b4319ba4SBarry Smith .seealso: MATIS 871b4319ba4SBarry Smith @*/ 8727087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 873b4319ba4SBarry Smith { 8744ac538c5SBarry Smith PetscErrorCode ierr; 875b4319ba4SBarry Smith 876b4319ba4SBarry Smith PetscFunctionBegin; 8770700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 878b4319ba4SBarry Smith PetscValidPointer(local,2); 8794ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 880b4319ba4SBarry Smith PetscFunctionReturn(0); 881b4319ba4SBarry Smith } 882b4319ba4SBarry Smith 8833b03a366Sstefano_zampini #undef __FUNCT__ 8843b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 8853b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 8863b03a366Sstefano_zampini { 8873b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 8883b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 8893b03a366Sstefano_zampini PetscErrorCode ierr; 8903b03a366Sstefano_zampini 8913b03a366Sstefano_zampini PetscFunctionBegin; 8924e4c7dbeSStefano Zampini if (is->A) { 8933b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 8943b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 895f23aa3ddSBarry 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); 8964e4c7dbeSStefano Zampini } 8973b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 8983b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 8993b03a366Sstefano_zampini is->A = local; 9003b03a366Sstefano_zampini PetscFunctionReturn(0); 9013b03a366Sstefano_zampini } 9023b03a366Sstefano_zampini 9033b03a366Sstefano_zampini #undef __FUNCT__ 9043b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 9053b03a366Sstefano_zampini /*@ 906eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 9073b03a366Sstefano_zampini 9083b03a366Sstefano_zampini Input Parameter: 9093b03a366Sstefano_zampini . mat - the matrix 910eb82efa4SStefano Zampini . local - the local matrix 9113b03a366Sstefano_zampini 9123b03a366Sstefano_zampini Output Parameter: 9133b03a366Sstefano_zampini 9143b03a366Sstefano_zampini Level: advanced 9153b03a366Sstefano_zampini 9163b03a366Sstefano_zampini Notes: 9173b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 9183b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 9193b03a366Sstefano_zampini 9203b03a366Sstefano_zampini .seealso: MATIS 9213b03a366Sstefano_zampini @*/ 9223b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 9233b03a366Sstefano_zampini { 9243b03a366Sstefano_zampini PetscErrorCode ierr; 9253b03a366Sstefano_zampini 9263b03a366Sstefano_zampini PetscFunctionBegin; 9273b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 928b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 9293b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 9303b03a366Sstefano_zampini PetscFunctionReturn(0); 9313b03a366Sstefano_zampini } 9323b03a366Sstefano_zampini 9336726f965SBarry Smith #undef __FUNCT__ 9346726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 9356726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 9366726f965SBarry Smith { 9376726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 9386726f965SBarry Smith PetscErrorCode ierr; 9396726f965SBarry Smith 9406726f965SBarry Smith PetscFunctionBegin; 9416726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 9426726f965SBarry Smith PetscFunctionReturn(0); 9436726f965SBarry Smith } 9446726f965SBarry Smith 9456726f965SBarry Smith #undef __FUNCT__ 9462e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 9472e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 9482e74eeadSLisandro Dalcin { 9492e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9502e74eeadSLisandro Dalcin PetscErrorCode ierr; 9512e74eeadSLisandro Dalcin 9522e74eeadSLisandro Dalcin PetscFunctionBegin; 9532e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 9542e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9552e74eeadSLisandro Dalcin } 9562e74eeadSLisandro Dalcin 9572e74eeadSLisandro Dalcin #undef __FUNCT__ 9582e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 9592e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 9602e74eeadSLisandro Dalcin { 9612e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9622e74eeadSLisandro Dalcin PetscErrorCode ierr; 9632e74eeadSLisandro Dalcin 9642e74eeadSLisandro Dalcin PetscFunctionBegin; 9652e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 966e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 9672e74eeadSLisandro Dalcin 9682e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 9692e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 970e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 971e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9722e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9732e74eeadSLisandro Dalcin } 9742e74eeadSLisandro Dalcin 9752e74eeadSLisandro Dalcin #undef __FUNCT__ 9766726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 977ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 9786726f965SBarry Smith { 9796726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 9806726f965SBarry Smith PetscErrorCode ierr; 9816726f965SBarry Smith 9826726f965SBarry Smith PetscFunctionBegin; 9834e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 9846726f965SBarry Smith PetscFunctionReturn(0); 9856726f965SBarry Smith } 9866726f965SBarry Smith 987284134d9SBarry Smith #undef __FUNCT__ 988284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 989284134d9SBarry Smith /*@ 990284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 991284134d9SBarry Smith process but not across processes. 992284134d9SBarry Smith 993284134d9SBarry Smith Input Parameters: 994284134d9SBarry Smith + comm - MPI communicator that will share the matrix 995e176bc59SStefano Zampini . bs - block size of the matrix 996df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 997e176bc59SStefano Zampini . rmap - local to global map for rows 998e176bc59SStefano Zampini - cmap - local to global map for cols 999284134d9SBarry Smith 1000284134d9SBarry Smith Output Parameter: 1001284134d9SBarry Smith . A - the resulting matrix 1002284134d9SBarry Smith 10038e6c10adSSatish Balay Level: advanced 10048e6c10adSSatish Balay 1005284134d9SBarry Smith Notes: See MATIS for more details 10068cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 1007e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 1008e176bc59SStefano Zampini If either rmap or cmap are NULL, than the matrix is assumed to be square 1009284134d9SBarry Smith 1010284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 1011284134d9SBarry Smith @*/ 1012e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 1013284134d9SBarry Smith { 1014284134d9SBarry Smith PetscErrorCode ierr; 1015284134d9SBarry Smith 1016284134d9SBarry Smith PetscFunctionBegin; 1017e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 1018284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1019d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 1020284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1021284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 1022e176bc59SStefano Zampini if (rmap && cmap) { 1023e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 1024e176bc59SStefano Zampini } else if (!rmap) { 1025e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 1026e176bc59SStefano Zampini } else { 1027e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 1028e176bc59SStefano Zampini } 1029284134d9SBarry Smith PetscFunctionReturn(0); 1030284134d9SBarry Smith } 1031284134d9SBarry Smith 1032b4319ba4SBarry Smith /*MC 1033eb82efa4SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC). 1034b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 1035b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 1036b4319ba4SBarry Smith product is handled "implicitly". 1037b4319ba4SBarry Smith 1038b4319ba4SBarry Smith Operations Provided: 10396726f965SBarry Smith + MatMult() 10402e74eeadSLisandro Dalcin . MatMultAdd() 10412e74eeadSLisandro Dalcin . MatMultTranspose() 10422e74eeadSLisandro Dalcin . MatMultTransposeAdd() 10436726f965SBarry Smith . MatZeroEntries() 10446726f965SBarry Smith . MatSetOption() 10452e74eeadSLisandro Dalcin . MatZeroRows() 10466726f965SBarry Smith . MatZeroRowsLocal() 10472e74eeadSLisandro Dalcin . MatSetValues() 10486726f965SBarry Smith . MatSetValuesLocal() 10492e74eeadSLisandro Dalcin . MatScale() 10502e74eeadSLisandro Dalcin . MatGetDiagonal() 10516726f965SBarry Smith - MatSetLocalToGlobalMapping() 1052b4319ba4SBarry Smith 1053b4319ba4SBarry Smith Options Database Keys: 1054b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1055b4319ba4SBarry Smith 1056b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1057b4319ba4SBarry Smith 1058b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1059b4319ba4SBarry Smith 1060b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1061eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1062b4319ba4SBarry Smith 1063b4319ba4SBarry Smith Level: advanced 1064b4319ba4SBarry Smith 1065eb82efa4SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC 1066b4319ba4SBarry Smith 1067b4319ba4SBarry Smith M*/ 1068b4319ba4SBarry Smith 1069b4319ba4SBarry Smith #undef __FUNCT__ 1070b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 10718cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1072b4319ba4SBarry Smith { 1073dfbe8321SBarry Smith PetscErrorCode ierr; 1074b4319ba4SBarry Smith Mat_IS *b; 1075b4319ba4SBarry Smith 1076b4319ba4SBarry Smith PetscFunctionBegin; 1077b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1078b4319ba4SBarry Smith A->data = (void*)b; 1079b4319ba4SBarry Smith 1080e176bc59SStefano Zampini /* matrix ops */ 1081e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1082b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 10832e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 10842e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 10852e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1086b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1087b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 10882e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 1089b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1090f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 10912e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1092b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 1093b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1094b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1095b4319ba4SBarry Smith A->ops->view = MatView_IS; 10966726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 10972e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 10982e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 10996726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 110069796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 110169796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1102ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 1103b4319ba4SBarry Smith 1104b7ce53b6SStefano Zampini /* special MATIS functions */ 1105bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1106bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1107b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 11082e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 1109cf0a3239SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",MatISSetUpSF_IS);CHKERRQ(ierr); 111017667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1111b4319ba4SBarry Smith PetscFunctionReturn(0); 1112b4319ba4SBarry Smith } 1113