1be1d678aSKris Buschelman 2b4319ba4SBarry Smith /* 3b4319ba4SBarry Smith Creates a matrix class for using the Neumann-Neumann type preconditioners. 4b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 5b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 6b4319ba4SBarry Smith product is handled "implicitly". 7b4319ba4SBarry Smith 8b4319ba4SBarry Smith Currently this allows for only one subdomain per processor. 9b4319ba4SBarry Smith */ 10b4319ba4SBarry Smith 11c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/ 1228f4e0baSStefano Zampini 13a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat); 147230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat,PetscInt,const PetscInt[],PetscScalar); 15a72627d2SStefano Zampini 1628f4e0baSStefano Zampini #undef __FUNCT__ 17659959c5SStefano Zampini #define __FUNCT__ "MatGetLocalSubMatrix_IS" 18659959c5SStefano Zampini PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat) 19659959c5SStefano Zampini { 20659959c5SStefano Zampini Mat lA,lsubmat; 21659959c5SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 22659959c5SStefano Zampini IS is,isn; 23659959c5SStefano Zampini const PetscInt *rg,*rl; 24659959c5SStefano Zampini #if defined(PETSC_USE_DEBUG) 25659959c5SStefano Zampini PetscInt nrg; 26659959c5SStefano Zampini #endif 27659959c5SStefano Zampini PetscInt N,M,nrl,i,*idxs; 28659959c5SStefano Zampini PetscErrorCode ierr; 29659959c5SStefano Zampini 30659959c5SStefano Zampini PetscFunctionBegin; 31659959c5SStefano Zampini /* extract local submatrix (it will complain if row and col exceed the local sizes) */ 32659959c5SStefano Zampini ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr); 33659959c5SStefano Zampini ierr = MatGetSubMatrix(lA,row,col,MAT_INITIAL_MATRIX,&lsubmat);CHKERRQ(ierr); 34659959c5SStefano Zampini /* compute new l2g map for rows */ 35659959c5SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 36659959c5SStefano Zampini ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr); 37659959c5SStefano Zampini ierr = ISGetIndices(row,&rl);CHKERRQ(ierr); 38659959c5SStefano Zampini #if defined(PETSC_USE_DEBUG) 39659959c5SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr); 40659959c5SStefano Zampini for (i=0;i<nrl;i++) if (rl[i]>=nrg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local row index %d -> %d greater than maximum possible %d\n",i,rl[i],nrg); 41659959c5SStefano Zampini #endif 42659959c5SStefano Zampini ierr = PetscMalloc1(nrl,&idxs);CHKERRQ(ierr); 43659959c5SStefano Zampini for (i=0;i<nrl;i++) idxs[i] = rg[rl[i]]; 44659959c5SStefano Zampini ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr); 45659959c5SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 46659959c5SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrl,idxs,PETSC_COPY_VALUES,&is);CHKERRQ(ierr); 47659959c5SStefano Zampini ierr = PetscFree(idxs);CHKERRQ(ierr); 48659959c5SStefano Zampini ierr = ISRenumber(is,NULL,&M,&isn);CHKERRQ(ierr); 49659959c5SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 50659959c5SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(isn,&rl2g);CHKERRQ(ierr); 51659959c5SStefano Zampini ierr = ISDestroy(&isn);CHKERRQ(ierr); 52659959c5SStefano Zampini /* compute new l2g map for columns */ 53659959c5SStefano Zampini if (col != row || A->rmap->mapping != A->cmap->mapping) { 54659959c5SStefano Zampini const PetscInt *cg,*cl; 55659959c5SStefano Zampini #if defined(PETSC_USE_DEBUG) 56659959c5SStefano Zampini PetscInt ncg; 57659959c5SStefano Zampini #endif 58659959c5SStefano Zampini PetscInt ncl; 59659959c5SStefano Zampini 60659959c5SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 61659959c5SStefano Zampini ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr); 62659959c5SStefano Zampini ierr = ISGetIndices(col,&cl);CHKERRQ(ierr); 63659959c5SStefano Zampini #if defined(PETSC_USE_DEBUG) 64659959c5SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr); 65659959c5SStefano Zampini for (i=0;i<ncl;i++) if (cl[i]>=ncg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local column index %d -> %d greater than maximum possible %d\n",i,cl[i],ncg); 66659959c5SStefano Zampini #endif 67659959c5SStefano Zampini ierr = PetscMalloc1(ncl,&idxs);CHKERRQ(ierr); 68659959c5SStefano Zampini for (i=0;i<ncl;i++) idxs[i] = cg[cl[i]]; 69659959c5SStefano Zampini ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr); 70659959c5SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 71659959c5SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncl,idxs,PETSC_COPY_VALUES,&is);CHKERRQ(ierr); 72659959c5SStefano Zampini ierr = PetscFree(idxs);CHKERRQ(ierr); 73659959c5SStefano Zampini ierr = ISRenumber(is,NULL,&N,&isn);CHKERRQ(ierr); 74659959c5SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 75659959c5SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(isn,&cl2g);CHKERRQ(ierr); 76659959c5SStefano Zampini ierr = ISDestroy(&isn);CHKERRQ(ierr); 77659959c5SStefano Zampini } else { 78*0a40dfa3SStefano Zampini ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr); 79659959c5SStefano Zampini cl2g = rl2g; 80659959c5SStefano Zampini N = M; 81659959c5SStefano Zampini } 82659959c5SStefano Zampini /* create the MATIS submatrix */ 83659959c5SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr); 84659959c5SStefano Zampini ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 85659959c5SStefano Zampini ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr); 86659959c5SStefano Zampini ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr); 87659959c5SStefano Zampini ierr = MatISSetLocalMat(*submat,lsubmat);CHKERRQ(ierr); 88659959c5SStefano Zampini ierr = MatDestroy(&lsubmat);CHKERRQ(ierr); 89659959c5SStefano Zampini ierr = MatSetUp(*submat);CHKERRQ(ierr); 90659959c5SStefano Zampini ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 91659959c5SStefano Zampini ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 92659959c5SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 93659959c5SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 94659959c5SStefano Zampini PetscFunctionReturn(0); 95659959c5SStefano Zampini } 96659959c5SStefano Zampini 97659959c5SStefano Zampini #undef __FUNCT__ 982b404112SStefano Zampini #define __FUNCT__ "MatCopy_IS" 992b404112SStefano Zampini PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str) 1002b404112SStefano Zampini { 1012b404112SStefano Zampini Mat_IS *a = (Mat_IS*)A->data,*b; 1022b404112SStefano Zampini PetscBool ismatis; 1032b404112SStefano Zampini PetscErrorCode ierr; 1042b404112SStefano Zampini 1052b404112SStefano Zampini PetscFunctionBegin; 1062b404112SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr); 1072b404112SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented"); 1082b404112SStefano Zampini b = (Mat_IS*)B->data; 1092b404112SStefano Zampini ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr); 1102b404112SStefano Zampini PetscFunctionReturn(0); 1112b404112SStefano Zampini } 1122b404112SStefano Zampini 1132b404112SStefano Zampini #undef __FUNCT__ 1146bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS" 1156bd84002SStefano Zampini PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool *missing,PetscInt *d) 1166bd84002SStefano Zampini { 1176bd84002SStefano Zampini Mat_IS *a = (Mat_IS*)A->data; 1186bd84002SStefano Zampini PetscErrorCode ierr; 1196bd84002SStefano Zampini 1206bd84002SStefano Zampini PetscFunctionBegin; 1216bd84002SStefano Zampini if (d) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Need to be implemented"); 1226bd84002SStefano Zampini ierr = MatMissingDiagonal(a->A,missing,NULL);CHKERRQ(ierr); 1236bd84002SStefano Zampini PetscFunctionReturn(0); 1246bd84002SStefano Zampini } 1256bd84002SStefano Zampini 1266bd84002SStefano Zampini #undef __FUNCT__ 12728f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private" 128a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B) 12928f4e0baSStefano Zampini { 13028f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 13128f4e0baSStefano Zampini const PetscInt *gidxs; 13228f4e0baSStefano Zampini PetscErrorCode ierr; 13328f4e0baSStefano Zampini 13428f4e0baSStefano Zampini PetscFunctionBegin; 13528f4e0baSStefano Zampini ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr); 13628f4e0baSStefano Zampini ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr); 13728f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 1383bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 13928f4e0baSStefano Zampini /* PETSC_OWN_POINTER refers to ilocal which is NULL */ 14028f4e0baSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 1413bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 14228f4e0baSStefano Zampini ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 14328f4e0baSStefano Zampini PetscFunctionReturn(0); 14428f4e0baSStefano Zampini } 1452e1947a5SStefano Zampini 1462e1947a5SStefano Zampini #undef __FUNCT__ 1472e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 148eb82efa4SStefano Zampini /*@ 149a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 150a88811baSStefano Zampini 151a88811baSStefano Zampini Collective on MPI_Comm 152a88811baSStefano Zampini 153a88811baSStefano Zampini Input Parameters: 154a88811baSStefano Zampini + B - the matrix 155a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 156a88811baSStefano Zampini (same value is used for all local rows) 157a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 158a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 159a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 160a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 161a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 162a88811baSStefano Zampini the diagonal entry even if it is zero. 163a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 164a88811baSStefano Zampini submatrix (same value is used for all local rows). 165a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 166a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 167a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 168a88811baSStefano Zampini structure. The size of this array is equal to the number 169a88811baSStefano Zampini of local rows, i.e 'm'. 170a88811baSStefano Zampini 171a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 172a88811baSStefano Zampini 173a88811baSStefano Zampini Level: intermediate 174a88811baSStefano Zampini 175a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 176a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 177a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 178a88811baSStefano Zampini 179a88811baSStefano Zampini .keywords: matrix 180a88811baSStefano Zampini 1813c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS 182a88811baSStefano Zampini @*/ 1832e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 1842e1947a5SStefano Zampini { 1852e1947a5SStefano Zampini PetscErrorCode ierr; 1862e1947a5SStefano Zampini 1872e1947a5SStefano Zampini PetscFunctionBegin; 1882e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 1892e1947a5SStefano Zampini PetscValidType(B,1); 1902e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 1912e1947a5SStefano Zampini PetscFunctionReturn(0); 1922e1947a5SStefano Zampini } 1932e1947a5SStefano Zampini 1942e1947a5SStefano Zampini #undef __FUNCT__ 1952e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 1967230de76SStefano Zampini static PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 1972e1947a5SStefano Zampini { 1982e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 19928f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 2002e1947a5SStefano Zampini PetscErrorCode ierr; 2012e1947a5SStefano Zampini 2022e1947a5SStefano Zampini PetscFunctionBegin; 2036c4ed002SBarry Smith if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 20428f4e0baSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 20528f4e0baSStefano Zampini ierr = MatISComputeSF_Private(B);CHKERRQ(ierr); 20628f4e0baSStefano Zampini } 2072e1947a5SStefano Zampini if (!d_nnz) { 20828f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz; 2092e1947a5SStefano Zampini } else { 21028f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 2112e1947a5SStefano Zampini } 2122e1947a5SStefano Zampini if (!o_nnz) { 21328f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz; 2142e1947a5SStefano Zampini } else { 21528f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 2162e1947a5SStefano Zampini } 21728f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 21828f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 21928f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 22028f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 22128f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves;i++) { 22228f4e0baSStefano Zampini matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 2232e1947a5SStefano Zampini } 22428f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 22528f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 22628f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 2272e1947a5SStefano Zampini } 22828f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 22928f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 23028f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 2312e1947a5SStefano Zampini } 23228f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 2332e1947a5SStefano Zampini PetscFunctionReturn(0); 2342e1947a5SStefano Zampini } 235b4319ba4SBarry Smith 236b4319ba4SBarry Smith #undef __FUNCT__ 2373927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private" 2383927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 2393927de2eSStefano Zampini { 2403927de2eSStefano Zampini Mat_IS *matis = (Mat_IS*)(A->data); 2413927de2eSStefano Zampini PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 242ecf5a873SStefano Zampini const PetscInt *global_indices_r,*global_indices_c; 2433927de2eSStefano Zampini PetscInt i,j,bs,rows,cols; 2443927de2eSStefano Zampini PetscInt lrows,lcols; 2453927de2eSStefano Zampini PetscInt local_rows,local_cols; 2463927de2eSStefano Zampini PetscMPIInt nsubdomains; 2473927de2eSStefano Zampini PetscBool isdense,issbaij; 2483927de2eSStefano Zampini PetscErrorCode ierr; 2493927de2eSStefano Zampini 2503927de2eSStefano Zampini PetscFunctionBegin; 2513927de2eSStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr); 2523927de2eSStefano Zampini ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr); 2533927de2eSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 2543927de2eSStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 2553927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 2563927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 257ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 258ecf5a873SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 2597230de76SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 260ecf5a873SStefano Zampini } else { 261ecf5a873SStefano Zampini global_indices_c = global_indices_r; 262ecf5a873SStefano Zampini } 263ecf5a873SStefano Zampini 2643927de2eSStefano Zampini if (issbaij) { 2653927de2eSStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 2663927de2eSStefano Zampini } 2673927de2eSStefano Zampini /* 268ecf5a873SStefano Zampini An SF reduce is needed to sum up properly on shared rows. 2693927de2eSStefano Zampini Note that generally preallocation is not exact, since it overestimates nonzeros 2703927de2eSStefano Zampini */ 2713927de2eSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 2723927de2eSStefano Zampini ierr = MatISComputeSF_Private(A);CHKERRQ(ierr); 2733927de2eSStefano Zampini } 2743927de2eSStefano Zampini ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr); 2753927de2eSStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr); 2763927de2eSStefano Zampini /* All processes need to compute entire row ownership */ 2773927de2eSStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 2783927de2eSStefano Zampini ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 2793927de2eSStefano Zampini for (i=0;i<nsubdomains;i++) { 2803927de2eSStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 2813927de2eSStefano Zampini row_ownership[j] = i; 2823927de2eSStefano Zampini } 2833927de2eSStefano Zampini } 2847230de76SStefano Zampini ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 2853927de2eSStefano Zampini 2863927de2eSStefano Zampini /* 2873927de2eSStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 2883927de2eSStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 2893927de2eSStefano Zampini */ 2903927de2eSStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 2913927de2eSStefano Zampini /* preallocation as a MATAIJ */ 2923927de2eSStefano Zampini if (isdense) { /* special case for dense local matrices */ 2933927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 294ecf5a873SStefano Zampini PetscInt index_row = global_indices_r[i]; 2953927de2eSStefano Zampini for (j=i;j<local_rows;j++) { 2963927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 297ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[j]; 2983927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 2993927de2eSStefano Zampini my_dnz[i] += 1; 3003927de2eSStefano Zampini } else { /* offdiag block */ 3013927de2eSStefano Zampini my_onz[i] += 1; 3023927de2eSStefano Zampini } 3033927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 3043927de2eSStefano Zampini if (i != j) { 3053927de2eSStefano Zampini owner = row_ownership[index_col]; 3063927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 3073927de2eSStefano Zampini my_dnz[j] += 1; 3083927de2eSStefano Zampini } else { 3093927de2eSStefano Zampini my_onz[j] += 1; 3103927de2eSStefano Zampini } 3113927de2eSStefano Zampini } 3123927de2eSStefano Zampini } 3133927de2eSStefano Zampini } 314ecf5a873SStefano Zampini } else { /* TODO: this could be optimized using MatGetRowIJ */ 3153927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 3163927de2eSStefano Zampini const PetscInt *cols; 317ecf5a873SStefano Zampini PetscInt ncols,index_row = global_indices_r[i]; 3183927de2eSStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 3193927de2eSStefano Zampini for (j=0;j<ncols;j++) { 3203927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 321ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[cols[j]]; 3223927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 3233927de2eSStefano Zampini my_dnz[i] += 1; 3243927de2eSStefano Zampini } else { /* offdiag block */ 3253927de2eSStefano Zampini my_onz[i] += 1; 3263927de2eSStefano Zampini } 3273927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 328d9a9e74cSStefano Zampini if (issbaij && index_col != index_row) { 3293927de2eSStefano Zampini owner = row_ownership[index_col]; 3303927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 331d9a9e74cSStefano Zampini my_dnz[cols[j]] += 1; 3323927de2eSStefano Zampini } else { 333d9a9e74cSStefano Zampini my_onz[cols[j]] += 1; 3343927de2eSStefano Zampini } 3353927de2eSStefano Zampini } 3363927de2eSStefano Zampini } 3373927de2eSStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 3383927de2eSStefano Zampini } 3393927de2eSStefano Zampini } 340ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 341ecf5a873SStefano Zampini if (global_indices_c != global_indices_r) { 3427230de76SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 343ecf5a873SStefano Zampini } 3443927de2eSStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 345ecf5a873SStefano Zampini 346ecf5a873SStefano Zampini /* Reduce my_dnz and my_onz */ 3473927de2eSStefano Zampini if (maxreduce) { 3483927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 3493927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 3503927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 3513927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 3523927de2eSStefano Zampini } else { 3533927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 3543927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 3553927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 3563927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 3573927de2eSStefano Zampini } 3583927de2eSStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 3593927de2eSStefano Zampini 3603927de2eSStefano Zampini /* Resize preallocation if overestimated */ 3613927de2eSStefano Zampini for (i=0;i<lrows;i++) { 3623927de2eSStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 3633927de2eSStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 3643927de2eSStefano Zampini } 3653927de2eSStefano Zampini /* set preallocation */ 3663927de2eSStefano Zampini ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 3673927de2eSStefano Zampini for (i=0;i<lrows/bs;i++) { 3683927de2eSStefano Zampini dnz[i] = dnz[i*bs]/bs; 3693927de2eSStefano Zampini onz[i] = onz[i*bs]/bs; 3703927de2eSStefano Zampini } 3713927de2eSStefano Zampini ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 3723927de2eSStefano Zampini ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 3733927de2eSStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 3743927de2eSStefano Zampini if (issbaij) { 3753927de2eSStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 3763927de2eSStefano Zampini } 3773927de2eSStefano Zampini PetscFunctionReturn(0); 3783927de2eSStefano Zampini } 3793927de2eSStefano Zampini 3803927de2eSStefano Zampini #undef __FUNCT__ 381b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 3827230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 383b7ce53b6SStefano Zampini { 384b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 385d9a9e74cSStefano Zampini Mat local_mat; 386b7ce53b6SStefano Zampini /* info on mat */ 3873cfa4ea4SStefano Zampini PetscInt bs,rows,cols,lrows,lcols; 388b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 389686e3a49SStefano Zampini PetscBool isdense,issbaij,isseqaij; 3907c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 391b7ce53b6SStefano Zampini /* values insertion */ 392b7ce53b6SStefano Zampini PetscScalar *array; 393b7ce53b6SStefano Zampini /* work */ 394b7ce53b6SStefano Zampini PetscErrorCode ierr; 395b7ce53b6SStefano Zampini 396b7ce53b6SStefano Zampini PetscFunctionBegin; 397b7ce53b6SStefano Zampini /* get info from mat */ 3987c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 3997c03b4e8SStefano Zampini if (nsubdomains == 1) { 4007c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 4017c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 4027c03b4e8SStefano Zampini } else { 4037c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 4047c03b4e8SStefano Zampini } 4057c03b4e8SStefano Zampini PetscFunctionReturn(0); 4067c03b4e8SStefano Zampini } 407b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 408b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 4093cfa4ea4SStefano Zampini ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 410b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 411b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 412686e3a49SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 413b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 414b7ce53b6SStefano Zampini 415b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 416b7ce53b6SStefano Zampini MatType new_mat_type; 4173927de2eSStefano Zampini PetscBool issbaij_red; 418b7ce53b6SStefano Zampini 419b7ce53b6SStefano Zampini /* determining new matrix type */ 420b2566f29SBarry Smith ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 421b7ce53b6SStefano Zampini if (issbaij_red) { 422b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 423b7ce53b6SStefano Zampini } else { 424b7ce53b6SStefano Zampini if (bs>1) { 425b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 426b7ce53b6SStefano Zampini } else { 427b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 428b7ce53b6SStefano Zampini } 429b7ce53b6SStefano Zampini } 430b7ce53b6SStefano Zampini 4313927de2eSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 4323cfa4ea4SStefano Zampini ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 4333927de2eSStefano Zampini ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 4343927de2eSStefano Zampini ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr); 4353927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 436b7ce53b6SStefano Zampini } else { 4373cfa4ea4SStefano Zampini PetscInt mbs,mrows,mcols,mlrows,mlcols; 438b7ce53b6SStefano Zampini /* some checks */ 439b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 440b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 4413cfa4ea4SStefano Zampini ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 4426c4ed002SBarry Smith if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 4436c4ed002SBarry Smith if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 4446c4ed002SBarry Smith if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 4456c4ed002SBarry Smith if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 4466c4ed002SBarry Smith if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 447b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 448b7ce53b6SStefano Zampini } 449d9a9e74cSStefano Zampini 450d9a9e74cSStefano Zampini if (issbaij) { 451d9a9e74cSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 452d9a9e74cSStefano Zampini } else { 453d9a9e74cSStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 454d9a9e74cSStefano Zampini local_mat = matis->A; 455d9a9e74cSStefano Zampini } 456686e3a49SStefano Zampini 457b7ce53b6SStefano Zampini /* Set values */ 458ecf5a873SStefano Zampini ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 459b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 460ecf5a873SStefano Zampini PetscInt i,*dummy_rows,*dummy_cols; 461ecf5a873SStefano Zampini 462ecf5a873SStefano Zampini if (local_rows != local_cols) { 463ecf5a873SStefano Zampini ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 464ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 465ecf5a873SStefano Zampini for (i=0;i<local_cols;i++) dummy_cols[i] = i; 466ecf5a873SStefano Zampini } else { 467ecf5a873SStefano Zampini ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 468ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 469ecf5a873SStefano Zampini dummy_cols = dummy_rows; 470ecf5a873SStefano Zampini } 471b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 472d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 473ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 474d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 475ecf5a873SStefano Zampini if (dummy_rows != dummy_cols) { 476ecf5a873SStefano Zampini ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 477ecf5a873SStefano Zampini } else { 478ecf5a873SStefano Zampini ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 479ecf5a873SStefano Zampini } 480686e3a49SStefano Zampini } else if (isseqaij) { 481ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 482686e3a49SStefano Zampini PetscBool done; 483686e3a49SStefano Zampini 484d9a9e74cSStefano Zampini ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 4856c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 486d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 487686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 488ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 489686e3a49SStefano Zampini } 490d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 4916c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 492d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 493686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 494ecf5a873SStefano Zampini PetscInt i; 495c0962df8SStefano Zampini 496686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 497686e3a49SStefano Zampini PetscInt j; 498ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 499686e3a49SStefano Zampini 500ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 501ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 502ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 503686e3a49SStefano Zampini } 504b7ce53b6SStefano Zampini } 505b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 506d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 507b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 508b7ce53b6SStefano Zampini if (isdense) { 509b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 510b7ce53b6SStefano Zampini } 511b7ce53b6SStefano Zampini PetscFunctionReturn(0); 512b7ce53b6SStefano Zampini } 513b7ce53b6SStefano Zampini 514b7ce53b6SStefano Zampini #undef __FUNCT__ 515b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 516b7ce53b6SStefano Zampini /*@ 517b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 518b7ce53b6SStefano Zampini 519b7ce53b6SStefano Zampini Input Parameter: 520b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 521b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 522b7ce53b6SStefano Zampini 523b7ce53b6SStefano Zampini Output Parameter: 524b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 525b7ce53b6SStefano Zampini 526b7ce53b6SStefano Zampini Level: developer 527b7ce53b6SStefano Zampini 528eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 529b7ce53b6SStefano Zampini 530b7ce53b6SStefano Zampini .seealso: MATIS 531b7ce53b6SStefano Zampini @*/ 532b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 533b7ce53b6SStefano Zampini { 534b7ce53b6SStefano Zampini PetscErrorCode ierr; 535b7ce53b6SStefano Zampini 536b7ce53b6SStefano Zampini PetscFunctionBegin; 537b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 538b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 539b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 540b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 541b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 542b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 5436c4ed002SBarry Smith if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 544b7ce53b6SStefano Zampini } 545b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 546b7ce53b6SStefano Zampini PetscFunctionReturn(0); 547b7ce53b6SStefano Zampini } 548b7ce53b6SStefano Zampini 549b7ce53b6SStefano Zampini #undef __FUNCT__ 550ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 551ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 552ad6194a2SStefano Zampini { 553ad6194a2SStefano Zampini PetscErrorCode ierr; 554ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 555ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 556ad6194a2SStefano Zampini Mat B,localmat; 557ad6194a2SStefano Zampini 558ad6194a2SStefano Zampini PetscFunctionBegin; 559ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 560ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 561ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 562e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 563ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 564ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 565b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 566ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 567ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 568ad6194a2SStefano Zampini *newmat = B; 569ad6194a2SStefano Zampini PetscFunctionReturn(0); 570ad6194a2SStefano Zampini } 571ad6194a2SStefano Zampini 572ad6194a2SStefano Zampini #undef __FUNCT__ 57369796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 57469796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 57569796d55SStefano Zampini { 57669796d55SStefano Zampini PetscErrorCode ierr; 57769796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 57869796d55SStefano Zampini PetscBool local_sym; 57969796d55SStefano Zampini 58069796d55SStefano Zampini PetscFunctionBegin; 58169796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 582b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 58369796d55SStefano Zampini PetscFunctionReturn(0); 58469796d55SStefano Zampini } 58569796d55SStefano Zampini 58669796d55SStefano Zampini #undef __FUNCT__ 58769796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 58869796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 58969796d55SStefano Zampini { 59069796d55SStefano Zampini PetscErrorCode ierr; 59169796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 59269796d55SStefano Zampini PetscBool local_sym; 59369796d55SStefano Zampini 59469796d55SStefano Zampini PetscFunctionBegin; 59569796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 596b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 59769796d55SStefano Zampini PetscFunctionReturn(0); 59869796d55SStefano Zampini } 59969796d55SStefano Zampini 60069796d55SStefano Zampini #undef __FUNCT__ 601b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 602dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 603b4319ba4SBarry Smith { 604dfbe8321SBarry Smith PetscErrorCode ierr; 605b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 606b4319ba4SBarry Smith 607b4319ba4SBarry Smith PetscFunctionBegin; 6086bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 609e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 610e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 6116bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 6126bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 61328f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 61428f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 615bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 616dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 617bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 618b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 619b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 6202e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 621b4319ba4SBarry Smith PetscFunctionReturn(0); 622b4319ba4SBarry Smith } 623b4319ba4SBarry Smith 624b4319ba4SBarry Smith #undef __FUNCT__ 625b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 626dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 627b4319ba4SBarry Smith { 628dfbe8321SBarry Smith PetscErrorCode ierr; 629b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 630b4319ba4SBarry Smith PetscScalar zero = 0.0; 631b4319ba4SBarry Smith 632b4319ba4SBarry Smith PetscFunctionBegin; 633b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 634e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 635e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 636b4319ba4SBarry Smith 637b4319ba4SBarry Smith /* multiply the local matrix */ 638b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 639b4319ba4SBarry Smith 640b4319ba4SBarry Smith /* scatter product back into global memory */ 6412dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 642e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 643e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 644b4319ba4SBarry Smith PetscFunctionReturn(0); 645b4319ba4SBarry Smith } 646b4319ba4SBarry Smith 647b4319ba4SBarry Smith #undef __FUNCT__ 6482e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 649b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 6502e74eeadSLisandro Dalcin { 651650997f4SStefano Zampini Vec temp_vec; 6522e74eeadSLisandro Dalcin PetscErrorCode ierr; 6532e74eeadSLisandro Dalcin 6542e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 655650997f4SStefano Zampini if (v3 != v2) { 656650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 657650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 658650997f4SStefano Zampini } else { 659650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 660650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 661650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 662650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 663650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 664650997f4SStefano Zampini } 6652e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6662e74eeadSLisandro Dalcin } 6672e74eeadSLisandro Dalcin 6682e74eeadSLisandro Dalcin #undef __FUNCT__ 6692e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 670e176bc59SStefano Zampini PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 6712e74eeadSLisandro Dalcin { 6722e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 6732e74eeadSLisandro Dalcin PetscErrorCode ierr; 6742e74eeadSLisandro Dalcin 675e176bc59SStefano Zampini PetscFunctionBegin; 6762e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 677e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 678e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6792e74eeadSLisandro Dalcin 6802e74eeadSLisandro Dalcin /* multiply the local matrix */ 681e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 6822e74eeadSLisandro Dalcin 6832e74eeadSLisandro Dalcin /* scatter product back into global vector */ 684e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 685e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 686e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 6872e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6882e74eeadSLisandro Dalcin } 6892e74eeadSLisandro Dalcin 6902e74eeadSLisandro Dalcin #undef __FUNCT__ 6912e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 6922e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 6932e74eeadSLisandro Dalcin { 694650997f4SStefano Zampini Vec temp_vec; 6952e74eeadSLisandro Dalcin PetscErrorCode ierr; 6962e74eeadSLisandro Dalcin 6972e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 698650997f4SStefano Zampini if (v3 != v2) { 699650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 700650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 701650997f4SStefano Zampini } else { 702650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 703650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 704650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 705650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 706650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 707650997f4SStefano Zampini } 7082e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7092e74eeadSLisandro Dalcin } 7102e74eeadSLisandro Dalcin 7112e74eeadSLisandro Dalcin #undef __FUNCT__ 712b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 713dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 714b4319ba4SBarry Smith { 715b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 716dfbe8321SBarry Smith PetscErrorCode ierr; 717b4319ba4SBarry Smith PetscViewer sviewer; 718b4319ba4SBarry Smith 719b4319ba4SBarry Smith PetscFunctionBegin; 7203f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 721b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 7223f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 7236e520ac8SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 724b4319ba4SBarry Smith PetscFunctionReturn(0); 725b4319ba4SBarry Smith } 726b4319ba4SBarry Smith 727b4319ba4SBarry Smith #undef __FUNCT__ 728b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 729784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 730b4319ba4SBarry Smith { 731dfbe8321SBarry Smith PetscErrorCode ierr; 732e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 733b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 734b4319ba4SBarry Smith IS from,to; 735e176bc59SStefano Zampini Vec cglobal,rglobal; 736b4319ba4SBarry Smith 737b4319ba4SBarry Smith PetscFunctionBegin; 738784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 739e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 7403bbff08aSStefano Zampini /* Destroy any previous data */ 74170cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 74270cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 743e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 744e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 7451c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 74628f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 74728f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 7483bbff08aSStefano Zampini 7493bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 750fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 751fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 752fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 753fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 754b4319ba4SBarry Smith 755b4319ba4SBarry Smith /* Create the local matrix A */ 756e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 757e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 758e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 759e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 760f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 761e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 762e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 763e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 764ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 765ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 766b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 767b4319ba4SBarry Smith 768b4319ba4SBarry Smith /* Create the local work vectors */ 7693bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 770b4319ba4SBarry Smith 771e176bc59SStefano Zampini /* setup the global to local scatters */ 772e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 773e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 774784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 775e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 776e176bc59SStefano Zampini if (rmapping != cmapping) { 777e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 778e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 779e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 780e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 781e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 782e176bc59SStefano Zampini } else { 783e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 784e176bc59SStefano Zampini is->cctx = is->rctx; 785e176bc59SStefano Zampini } 786e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 787e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 7886bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 7896bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 79048ff6bf3SStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 791b4319ba4SBarry Smith PetscFunctionReturn(0); 792b4319ba4SBarry Smith } 793b4319ba4SBarry Smith 7942e74eeadSLisandro Dalcin #undef __FUNCT__ 7952e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 7962e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 7972e74eeadSLisandro Dalcin { 7982e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 7992e74eeadSLisandro Dalcin PetscErrorCode ierr; 80097563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 80197563a80SStefano Zampini PetscInt i,zm,zn; 80297563a80SStefano Zampini #endif 80397563a80SStefano Zampini PetscInt rows_l[2048],cols_l[2048]; 8042e74eeadSLisandro Dalcin 8052e74eeadSLisandro Dalcin PetscFunctionBegin; 8062e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 807f23aa3ddSBarry 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); 80897563a80SStefano Zampini /* count negative indices */ 80997563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 81097563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 8112e74eeadSLisandro Dalcin #endif 81297563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 81397563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 81497563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 81597563a80SStefano Zampini /* count negative indices (should be the same as before) */ 81697563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 81797563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 81897563a80SStefano Zampini if (zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the row indices can not be mapped! Maybe you should not use MATIS"); 81997563a80SStefano Zampini if (zn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the column indices can not be mapped! Maybe you should not use MATIS"); 82097563a80SStefano Zampini #endif 8212e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 8222e74eeadSLisandro Dalcin PetscFunctionReturn(0); 8232e74eeadSLisandro Dalcin } 8242e74eeadSLisandro Dalcin 825b4319ba4SBarry Smith #undef __FUNCT__ 82697563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS" 82797563a80SStefano Zampini PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 82897563a80SStefano Zampini { 82997563a80SStefano Zampini Mat_IS *is = (Mat_IS*)mat->data; 83097563a80SStefano Zampini PetscErrorCode ierr; 83197563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 83297563a80SStefano Zampini PetscInt i,zm,zn; 83397563a80SStefano Zampini #endif 83497563a80SStefano Zampini PetscInt rows_l[2048],cols_l[2048]; 83597563a80SStefano Zampini 83697563a80SStefano Zampini PetscFunctionBegin; 83797563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 83897563a80SStefano Zampini if (m > 2048 || n > 2048) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n); 83997563a80SStefano Zampini /* count negative indices */ 84097563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 84197563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 84297563a80SStefano Zampini #endif 84397563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 84497563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 84597563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 84697563a80SStefano Zampini /* count negative indices (should be the same as before) */ 84797563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 84897563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 84997563a80SStefano Zampini if (zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the row indices can not be mapped! Maybe you should not use MATIS"); 85097563a80SStefano Zampini if (zn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the column indices can not be mapped! Maybe you should not use MATIS"); 85197563a80SStefano Zampini #endif 85297563a80SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 85397563a80SStefano Zampini PetscFunctionReturn(0); 85497563a80SStefano Zampini } 85597563a80SStefano Zampini 85697563a80SStefano Zampini #undef __FUNCT__ 857b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 85813f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 859b4319ba4SBarry Smith { 860dfbe8321SBarry Smith PetscErrorCode ierr; 861b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 862b4319ba4SBarry Smith 863b4319ba4SBarry Smith PetscFunctionBegin; 864b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 865b4319ba4SBarry Smith PetscFunctionReturn(0); 866b4319ba4SBarry Smith } 867b4319ba4SBarry Smith 868b4319ba4SBarry Smith #undef __FUNCT__ 869f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 870f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 871f0006bf2SLisandro Dalcin { 872f0006bf2SLisandro Dalcin PetscErrorCode ierr; 873f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 874f0006bf2SLisandro Dalcin 875f0006bf2SLisandro Dalcin PetscFunctionBegin; 876f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 877f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 878f0006bf2SLisandro Dalcin } 879f0006bf2SLisandro Dalcin 880f0006bf2SLisandro Dalcin #undef __FUNCT__ 8812e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 8822b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 8832e74eeadSLisandro Dalcin { 8846e520ac8SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 8856e520ac8SStefano Zampini PetscInt nr,nl,len,i; 8866e520ac8SStefano Zampini PetscInt *lrows; 8872e74eeadSLisandro Dalcin PetscErrorCode ierr; 8882e74eeadSLisandro Dalcin 8892e74eeadSLisandro Dalcin PetscFunctionBegin; 8906e520ac8SStefano Zampini /* get locally owned rows */ 8916e520ac8SStefano Zampini ierr = MatZeroRowsMapLocal_Private(A,n,rows,&len,&lrows);CHKERRQ(ierr); 8926e520ac8SStefano Zampini /* fix right hand side if needed */ 8936e520ac8SStefano Zampini if (x && b) { 8946e520ac8SStefano Zampini const PetscScalar *xx; 8956e520ac8SStefano Zampini PetscScalar *bb; 8966e520ac8SStefano Zampini 8976e520ac8SStefano Zampini ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr); 8986e520ac8SStefano Zampini ierr = VecGetArray(b, &bb);CHKERRQ(ierr); 8996e520ac8SStefano Zampini for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]]; 9006e520ac8SStefano Zampini ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr); 9016e520ac8SStefano Zampini ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr); 9022e74eeadSLisandro Dalcin } 9036e520ac8SStefano Zampini /* get rows associated to the local matrices */ 9046e520ac8SStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 9056e520ac8SStefano Zampini ierr = MatISComputeSF_Private(A);CHKERRQ(ierr); 9066e520ac8SStefano Zampini } 9076e520ac8SStefano Zampini ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr); 9086e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr); 9096e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 9106e520ac8SStefano Zampini for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1; 9116e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 9126e520ac8SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 9136e520ac8SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 9146e520ac8SStefano Zampini ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr); 9156e520ac8SStefano Zampini for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i; 9167230de76SStefano Zampini ierr = MatISZeroRowsLocal_Private(A,nr,lrows,diag);CHKERRQ(ierr); 9176e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 9182e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9192e74eeadSLisandro Dalcin } 9202e74eeadSLisandro Dalcin 9212e74eeadSLisandro Dalcin #undef __FUNCT__ 9227230de76SStefano Zampini #define __FUNCT__ "MatISZeroRowsLocal_Private" 9237230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag) 924b4319ba4SBarry Smith { 925b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 926dfbe8321SBarry Smith PetscErrorCode ierr; 927b4319ba4SBarry Smith 928b4319ba4SBarry Smith PetscFunctionBegin; 9296e520ac8SStefano Zampini if (diag != 0.) { 930b4319ba4SBarry Smith /* 931b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 932b4319ba4SBarry Smith work properly in the interface nodes. 933b4319ba4SBarry Smith */ 934b4319ba4SBarry Smith Vec counter; 935e176bc59SStefano Zampini ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr); 936e176bc59SStefano Zampini ierr = VecSet(counter,0.);CHKERRQ(ierr); 937e176bc59SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 938e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 939e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 940e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 941e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9426bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 943b4319ba4SBarry Smith } 9442b404112SStefano Zampini 945958c9bccSBarry Smith if (!n) { 946b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 947b4319ba4SBarry Smith } else { 9487230de76SStefano Zampini PetscInt i; 949b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 9502205254eSKarl Rupp 9512b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 9526e520ac8SStefano Zampini if (diag != 0.) { 9536e520ac8SStefano Zampini const PetscScalar *array; 9546e520ac8SStefano Zampini ierr = VecGetArrayRead(is->y,&array);CHKERRQ(ierr); 955b4319ba4SBarry Smith for (i=0; i<n; i++) { 956f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 957b4319ba4SBarry Smith } 9586e520ac8SStefano Zampini ierr = VecRestoreArrayRead(is->y,&array);CHKERRQ(ierr); 9596e520ac8SStefano Zampini } 960b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 961b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 962b4319ba4SBarry Smith } 963b4319ba4SBarry Smith PetscFunctionReturn(0); 964b4319ba4SBarry Smith } 965b4319ba4SBarry Smith 966b4319ba4SBarry Smith #undef __FUNCT__ 967b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 968dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 969b4319ba4SBarry Smith { 970b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 971dfbe8321SBarry Smith PetscErrorCode ierr; 972b4319ba4SBarry Smith 973b4319ba4SBarry Smith PetscFunctionBegin; 974b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 975b4319ba4SBarry Smith PetscFunctionReturn(0); 976b4319ba4SBarry Smith } 977b4319ba4SBarry Smith 978b4319ba4SBarry Smith #undef __FUNCT__ 979b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 980dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 981b4319ba4SBarry Smith { 982b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 983dfbe8321SBarry Smith PetscErrorCode ierr; 984b4319ba4SBarry Smith 985b4319ba4SBarry Smith PetscFunctionBegin; 986b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 987b4319ba4SBarry Smith PetscFunctionReturn(0); 988b4319ba4SBarry Smith } 989b4319ba4SBarry Smith 990b4319ba4SBarry Smith #undef __FUNCT__ 991b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 9927087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 993b4319ba4SBarry Smith { 994b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 995b4319ba4SBarry Smith 996b4319ba4SBarry Smith PetscFunctionBegin; 997b4319ba4SBarry Smith *local = is->A; 998b4319ba4SBarry Smith PetscFunctionReturn(0); 999b4319ba4SBarry Smith } 1000b4319ba4SBarry Smith 1001b4319ba4SBarry Smith #undef __FUNCT__ 1002b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 1003b4319ba4SBarry Smith /*@ 1004b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 1005b4319ba4SBarry Smith 1006b4319ba4SBarry Smith Input Parameter: 1007b4319ba4SBarry Smith . mat - the matrix 1008b4319ba4SBarry Smith 1009b4319ba4SBarry Smith Output Parameter: 1010eb82efa4SStefano Zampini . local - the local matrix 1011b4319ba4SBarry Smith 1012b4319ba4SBarry Smith Level: advanced 1013b4319ba4SBarry Smith 1014b4319ba4SBarry Smith Notes: 1015b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 1016b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 1017b4319ba4SBarry Smith of the MatSetValues() operation. 1018b4319ba4SBarry Smith 1019b4319ba4SBarry Smith .seealso: MATIS 1020b4319ba4SBarry Smith @*/ 10217087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 1022b4319ba4SBarry Smith { 10234ac538c5SBarry Smith PetscErrorCode ierr; 1024b4319ba4SBarry Smith 1025b4319ba4SBarry Smith PetscFunctionBegin; 10260700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1027b4319ba4SBarry Smith PetscValidPointer(local,2); 10284ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 1029b4319ba4SBarry Smith PetscFunctionReturn(0); 1030b4319ba4SBarry Smith } 1031b4319ba4SBarry Smith 10323b03a366Sstefano_zampini #undef __FUNCT__ 10333b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 10343b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 10353b03a366Sstefano_zampini { 10363b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 10373b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 10383b03a366Sstefano_zampini PetscErrorCode ierr; 10393b03a366Sstefano_zampini 10403b03a366Sstefano_zampini PetscFunctionBegin; 10414e4c7dbeSStefano Zampini if (is->A) { 10423b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 10433b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 1044f23aa3ddSBarry 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); 10454e4c7dbeSStefano Zampini } 10463b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 10473b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 10483b03a366Sstefano_zampini is->A = local; 10493b03a366Sstefano_zampini PetscFunctionReturn(0); 10503b03a366Sstefano_zampini } 10513b03a366Sstefano_zampini 10523b03a366Sstefano_zampini #undef __FUNCT__ 10533b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 10543b03a366Sstefano_zampini /*@ 1055eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 10563b03a366Sstefano_zampini 10573b03a366Sstefano_zampini Input Parameter: 10583b03a366Sstefano_zampini . mat - the matrix 1059eb82efa4SStefano Zampini . local - the local matrix 10603b03a366Sstefano_zampini 10613b03a366Sstefano_zampini Output Parameter: 10623b03a366Sstefano_zampini 10633b03a366Sstefano_zampini Level: advanced 10643b03a366Sstefano_zampini 10653b03a366Sstefano_zampini Notes: 10663b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 10673b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 10683b03a366Sstefano_zampini 10693b03a366Sstefano_zampini .seealso: MATIS 10703b03a366Sstefano_zampini @*/ 10713b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 10723b03a366Sstefano_zampini { 10733b03a366Sstefano_zampini PetscErrorCode ierr; 10743b03a366Sstefano_zampini 10753b03a366Sstefano_zampini PetscFunctionBegin; 10763b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1077b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 10783b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 10793b03a366Sstefano_zampini PetscFunctionReturn(0); 10803b03a366Sstefano_zampini } 10813b03a366Sstefano_zampini 10826726f965SBarry Smith #undef __FUNCT__ 10836726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 10846726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 10856726f965SBarry Smith { 10866726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 10876726f965SBarry Smith PetscErrorCode ierr; 10886726f965SBarry Smith 10896726f965SBarry Smith PetscFunctionBegin; 10906726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 10916726f965SBarry Smith PetscFunctionReturn(0); 10926726f965SBarry Smith } 10936726f965SBarry Smith 10946726f965SBarry Smith #undef __FUNCT__ 10952e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 10962e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 10972e74eeadSLisandro Dalcin { 10982e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 10992e74eeadSLisandro Dalcin PetscErrorCode ierr; 11002e74eeadSLisandro Dalcin 11012e74eeadSLisandro Dalcin PetscFunctionBegin; 11022e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 11032e74eeadSLisandro Dalcin PetscFunctionReturn(0); 11042e74eeadSLisandro Dalcin } 11052e74eeadSLisandro Dalcin 11062e74eeadSLisandro Dalcin #undef __FUNCT__ 11072e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 11082e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 11092e74eeadSLisandro Dalcin { 11102e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 11112e74eeadSLisandro Dalcin PetscErrorCode ierr; 11122e74eeadSLisandro Dalcin 11132e74eeadSLisandro Dalcin PetscFunctionBegin; 11142e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 1115e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 11162e74eeadSLisandro Dalcin 11172e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 11182e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 1119e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1120e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11212e74eeadSLisandro Dalcin PetscFunctionReturn(0); 11222e74eeadSLisandro Dalcin } 11232e74eeadSLisandro Dalcin 11242e74eeadSLisandro Dalcin #undef __FUNCT__ 11256726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 1126ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 11276726f965SBarry Smith { 11286726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 11296726f965SBarry Smith PetscErrorCode ierr; 11306726f965SBarry Smith 11316726f965SBarry Smith PetscFunctionBegin; 11324e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 11336726f965SBarry Smith PetscFunctionReturn(0); 11346726f965SBarry Smith } 11356726f965SBarry Smith 1136284134d9SBarry Smith #undef __FUNCT__ 1137284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 1138284134d9SBarry Smith /*@ 11393c212e90SHong Zhang MatCreateIS - Creates a "process" unassembled matrix, assembled on each 1140284134d9SBarry Smith process but not across processes. 1141284134d9SBarry Smith 1142284134d9SBarry Smith Input Parameters: 1143284134d9SBarry Smith + comm - MPI communicator that will share the matrix 1144e176bc59SStefano Zampini . bs - block size of the matrix 1145df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 1146e176bc59SStefano Zampini . rmap - local to global map for rows 1147e176bc59SStefano Zampini - cmap - local to global map for cols 1148284134d9SBarry Smith 1149284134d9SBarry Smith Output Parameter: 1150284134d9SBarry Smith . A - the resulting matrix 1151284134d9SBarry Smith 11528e6c10adSSatish Balay Level: advanced 11538e6c10adSSatish Balay 11543c212e90SHong Zhang Notes: See MATIS for more details. 11553c212e90SHong Zhang m and n are NOT related to the size of the map; they are the size of the part of the vector owned 1156e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 11573c212e90SHong Zhang If either rmap or cmap are NULL, then the matrix is assumed to be square. 1158284134d9SBarry Smith 1159284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 1160284134d9SBarry Smith @*/ 1161e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 1162284134d9SBarry Smith { 1163284134d9SBarry Smith PetscErrorCode ierr; 1164284134d9SBarry Smith 1165284134d9SBarry Smith PetscFunctionBegin; 1166e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 1167284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1168d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 1169284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1170284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 1171e176bc59SStefano Zampini if (rmap && cmap) { 1172e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 1173e176bc59SStefano Zampini } else if (!rmap) { 1174e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 1175e176bc59SStefano Zampini } else { 1176e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 1177e176bc59SStefano Zampini } 1178284134d9SBarry Smith PetscFunctionReturn(0); 1179284134d9SBarry Smith } 1180284134d9SBarry Smith 1181b4319ba4SBarry Smith /*MC 1182eb82efa4SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC). 1183b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 1184b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 1185b4319ba4SBarry Smith product is handled "implicitly". 1186b4319ba4SBarry Smith 1187b4319ba4SBarry Smith Operations Provided: 11886726f965SBarry Smith + MatMult() 11892e74eeadSLisandro Dalcin . MatMultAdd() 11902e74eeadSLisandro Dalcin . MatMultTranspose() 11912e74eeadSLisandro Dalcin . MatMultTransposeAdd() 11926726f965SBarry Smith . MatZeroEntries() 11936726f965SBarry Smith . MatSetOption() 11942e74eeadSLisandro Dalcin . MatZeroRows() 11952e74eeadSLisandro Dalcin . MatSetValues() 119697563a80SStefano Zampini . MatSetValuesBlocked() 11976726f965SBarry Smith . MatSetValuesLocal() 119897563a80SStefano Zampini . MatSetValuesBlockedLocal() 11992e74eeadSLisandro Dalcin . MatScale() 12002e74eeadSLisandro Dalcin . MatGetDiagonal() 12012b404112SStefano Zampini . MatMissingDiagonal() 12022b404112SStefano Zampini . MatDuplicate() 12032b404112SStefano Zampini . MatCopy() 12046726f965SBarry Smith - MatSetLocalToGlobalMapping() 1205b4319ba4SBarry Smith 1206b4319ba4SBarry Smith Options Database Keys: 1207b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1208b4319ba4SBarry Smith 1209b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1210b4319ba4SBarry Smith 1211b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1212b4319ba4SBarry Smith 1213b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1214eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1215b4319ba4SBarry Smith 1216b4319ba4SBarry Smith Level: advanced 1217b4319ba4SBarry Smith 12183c212e90SHong Zhang .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC, MatCreateIS() 1219b4319ba4SBarry Smith 1220b4319ba4SBarry Smith M*/ 1221b4319ba4SBarry Smith 1222b4319ba4SBarry Smith #undef __FUNCT__ 1223b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 12248cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1225b4319ba4SBarry Smith { 1226dfbe8321SBarry Smith PetscErrorCode ierr; 1227b4319ba4SBarry Smith Mat_IS *b; 1228b4319ba4SBarry Smith 1229b4319ba4SBarry Smith PetscFunctionBegin; 1230b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1231b4319ba4SBarry Smith A->data = (void*)b; 1232b4319ba4SBarry Smith 1233e176bc59SStefano Zampini /* matrix ops */ 1234e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1235b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 12362e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 12372e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 12382e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1239b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1240b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 12412e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 124298921651SStefano Zampini A->ops->setvaluesblocked = MatSetValuesBlocked_IS; 1243b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1244f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 12452e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1246b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1247b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1248b4319ba4SBarry Smith A->ops->view = MatView_IS; 12496726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 12502e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 12512e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 12526726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 125369796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 125469796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1255ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 12566bd84002SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_IS; 12572b404112SStefano Zampini A->ops->copy = MatCopy_IS; 1258659959c5SStefano Zampini A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_IS; 1259b4319ba4SBarry Smith 1260b7ce53b6SStefano Zampini /* special MATIS functions */ 1261bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1262bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1263b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 12642e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 126517667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1266b4319ba4SBarry Smith PetscFunctionReturn(0); 1267b4319ba4SBarry Smith } 1268