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 { 780a40dfa3SStefano 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 { 117*527b2640SStefano Zampini Vec v; 118*527b2640SStefano Zampini const PetscScalar *array; 119*527b2640SStefano Zampini PetscInt i,n; 1206bd84002SStefano Zampini PetscErrorCode ierr; 1216bd84002SStefano Zampini 1226bd84002SStefano Zampini PetscFunctionBegin; 123*527b2640SStefano Zampini *missing = PETSC_FALSE; 124*527b2640SStefano Zampini ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr); 125*527b2640SStefano Zampini ierr = MatGetDiagonal(A,v);CHKERRQ(ierr); 126*527b2640SStefano Zampini ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 127*527b2640SStefano Zampini ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr); 128*527b2640SStefano Zampini for (i=0;i<n;i++) if (array[i] == 0.) break; 129*527b2640SStefano Zampini ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr); 130*527b2640SStefano Zampini ierr = VecDestroy(&v);CHKERRQ(ierr); 131*527b2640SStefano Zampini if (i != n) *missing = PETSC_TRUE; 132*527b2640SStefano Zampini if (d) { 133*527b2640SStefano Zampini *d = -1; 134*527b2640SStefano Zampini if (*missing) { 135*527b2640SStefano Zampini PetscInt rstart; 136*527b2640SStefano Zampini ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 137*527b2640SStefano Zampini *d = i+rstart; 138*527b2640SStefano Zampini } 139*527b2640SStefano Zampini } 1406bd84002SStefano Zampini PetscFunctionReturn(0); 1416bd84002SStefano Zampini } 1426bd84002SStefano Zampini 1436bd84002SStefano Zampini #undef __FUNCT__ 14428f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private" 145a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B) 14628f4e0baSStefano Zampini { 14728f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 14828f4e0baSStefano Zampini const PetscInt *gidxs; 14928f4e0baSStefano Zampini PetscErrorCode ierr; 15028f4e0baSStefano Zampini 15128f4e0baSStefano Zampini PetscFunctionBegin; 15228f4e0baSStefano Zampini ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr); 15328f4e0baSStefano Zampini ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr); 15428f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 1553bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 15628f4e0baSStefano Zampini /* PETSC_OWN_POINTER refers to ilocal which is NULL */ 15728f4e0baSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 1583bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 15928f4e0baSStefano Zampini ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 16028f4e0baSStefano Zampini PetscFunctionReturn(0); 16128f4e0baSStefano Zampini } 1622e1947a5SStefano Zampini 1632e1947a5SStefano Zampini #undef __FUNCT__ 1642e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 165eb82efa4SStefano Zampini /*@ 166a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 167a88811baSStefano Zampini 168a88811baSStefano Zampini Collective on MPI_Comm 169a88811baSStefano Zampini 170a88811baSStefano Zampini Input Parameters: 171a88811baSStefano Zampini + B - the matrix 172a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 173a88811baSStefano Zampini (same value is used for all local rows) 174a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 175a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 176a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 177a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 178a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 179a88811baSStefano Zampini the diagonal entry even if it is zero. 180a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 181a88811baSStefano Zampini submatrix (same value is used for all local rows). 182a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 183a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 184a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 185a88811baSStefano Zampini structure. The size of this array is equal to the number 186a88811baSStefano Zampini of local rows, i.e 'm'. 187a88811baSStefano Zampini 188a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 189a88811baSStefano Zampini 190a88811baSStefano Zampini Level: intermediate 191a88811baSStefano Zampini 192a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 193a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 194a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 195a88811baSStefano Zampini 196a88811baSStefano Zampini .keywords: matrix 197a88811baSStefano Zampini 1983c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS 199a88811baSStefano Zampini @*/ 2002e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 2012e1947a5SStefano Zampini { 2022e1947a5SStefano Zampini PetscErrorCode ierr; 2032e1947a5SStefano Zampini 2042e1947a5SStefano Zampini PetscFunctionBegin; 2052e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 2062e1947a5SStefano Zampini PetscValidType(B,1); 2072e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 2082e1947a5SStefano Zampini PetscFunctionReturn(0); 2092e1947a5SStefano Zampini } 2102e1947a5SStefano Zampini 2112e1947a5SStefano Zampini #undef __FUNCT__ 2122e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 2137230de76SStefano Zampini static PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 2142e1947a5SStefano Zampini { 2152e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 21628f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 2172e1947a5SStefano Zampini PetscErrorCode ierr; 2182e1947a5SStefano Zampini 2192e1947a5SStefano Zampini PetscFunctionBegin; 2206c4ed002SBarry Smith if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 22128f4e0baSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 22228f4e0baSStefano Zampini ierr = MatISComputeSF_Private(B);CHKERRQ(ierr); 22328f4e0baSStefano Zampini } 2242e1947a5SStefano Zampini if (!d_nnz) { 22528f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz; 2262e1947a5SStefano Zampini } else { 22728f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 2282e1947a5SStefano Zampini } 2292e1947a5SStefano Zampini if (!o_nnz) { 23028f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz; 2312e1947a5SStefano Zampini } else { 23228f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 2332e1947a5SStefano Zampini } 23428f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 23528f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 23628f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 23728f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 23828f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves;i++) { 23928f4e0baSStefano Zampini matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 2402e1947a5SStefano Zampini } 24128f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 24228f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 24328f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 2442e1947a5SStefano Zampini } 24528f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 24628f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 24728f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 2482e1947a5SStefano Zampini } 24928f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 2502e1947a5SStefano Zampini PetscFunctionReturn(0); 2512e1947a5SStefano Zampini } 252b4319ba4SBarry Smith 253b4319ba4SBarry Smith #undef __FUNCT__ 2543927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private" 2553927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 2563927de2eSStefano Zampini { 2573927de2eSStefano Zampini Mat_IS *matis = (Mat_IS*)(A->data); 2583927de2eSStefano Zampini PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 259ecf5a873SStefano Zampini const PetscInt *global_indices_r,*global_indices_c; 2603927de2eSStefano Zampini PetscInt i,j,bs,rows,cols; 2613927de2eSStefano Zampini PetscInt lrows,lcols; 2623927de2eSStefano Zampini PetscInt local_rows,local_cols; 2633927de2eSStefano Zampini PetscMPIInt nsubdomains; 2643927de2eSStefano Zampini PetscBool isdense,issbaij; 2653927de2eSStefano Zampini PetscErrorCode ierr; 2663927de2eSStefano Zampini 2673927de2eSStefano Zampini PetscFunctionBegin; 2683927de2eSStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr); 2693927de2eSStefano Zampini ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr); 2703927de2eSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 2713927de2eSStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 2723927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 2733927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 274ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 275ecf5a873SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 2767230de76SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 277ecf5a873SStefano Zampini } else { 278ecf5a873SStefano Zampini global_indices_c = global_indices_r; 279ecf5a873SStefano Zampini } 280ecf5a873SStefano Zampini 2813927de2eSStefano Zampini if (issbaij) { 2823927de2eSStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 2833927de2eSStefano Zampini } 2843927de2eSStefano Zampini /* 285ecf5a873SStefano Zampini An SF reduce is needed to sum up properly on shared rows. 2863927de2eSStefano Zampini Note that generally preallocation is not exact, since it overestimates nonzeros 2873927de2eSStefano Zampini */ 2883927de2eSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 2893927de2eSStefano Zampini ierr = MatISComputeSF_Private(A);CHKERRQ(ierr); 2903927de2eSStefano Zampini } 2913927de2eSStefano Zampini ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr); 2923927de2eSStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr); 2933927de2eSStefano Zampini /* All processes need to compute entire row ownership */ 2943927de2eSStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 2953927de2eSStefano Zampini ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 2963927de2eSStefano Zampini for (i=0;i<nsubdomains;i++) { 2973927de2eSStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 2983927de2eSStefano Zampini row_ownership[j] = i; 2993927de2eSStefano Zampini } 3003927de2eSStefano Zampini } 3017230de76SStefano Zampini ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 3023927de2eSStefano Zampini 3033927de2eSStefano Zampini /* 3043927de2eSStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 3053927de2eSStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 3063927de2eSStefano Zampini */ 3073927de2eSStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 3083927de2eSStefano Zampini /* preallocation as a MATAIJ */ 3093927de2eSStefano Zampini if (isdense) { /* special case for dense local matrices */ 3103927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 311ecf5a873SStefano Zampini PetscInt index_row = global_indices_r[i]; 3123927de2eSStefano Zampini for (j=i;j<local_rows;j++) { 3133927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 314ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[j]; 3153927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 3163927de2eSStefano Zampini my_dnz[i] += 1; 3173927de2eSStefano Zampini } else { /* offdiag block */ 3183927de2eSStefano Zampini my_onz[i] += 1; 3193927de2eSStefano Zampini } 3203927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 3213927de2eSStefano Zampini if (i != j) { 3223927de2eSStefano Zampini owner = row_ownership[index_col]; 3233927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 3243927de2eSStefano Zampini my_dnz[j] += 1; 3253927de2eSStefano Zampini } else { 3263927de2eSStefano Zampini my_onz[j] += 1; 3273927de2eSStefano Zampini } 3283927de2eSStefano Zampini } 3293927de2eSStefano Zampini } 3303927de2eSStefano Zampini } 331ecf5a873SStefano Zampini } else { /* TODO: this could be optimized using MatGetRowIJ */ 3323927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 3333927de2eSStefano Zampini const PetscInt *cols; 334ecf5a873SStefano Zampini PetscInt ncols,index_row = global_indices_r[i]; 3353927de2eSStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 3363927de2eSStefano Zampini for (j=0;j<ncols;j++) { 3373927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 338ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[cols[j]]; 3393927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 3403927de2eSStefano Zampini my_dnz[i] += 1; 3413927de2eSStefano Zampini } else { /* offdiag block */ 3423927de2eSStefano Zampini my_onz[i] += 1; 3433927de2eSStefano Zampini } 3443927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 345d9a9e74cSStefano Zampini if (issbaij && index_col != index_row) { 3463927de2eSStefano Zampini owner = row_ownership[index_col]; 3473927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 348d9a9e74cSStefano Zampini my_dnz[cols[j]] += 1; 3493927de2eSStefano Zampini } else { 350d9a9e74cSStefano Zampini my_onz[cols[j]] += 1; 3513927de2eSStefano Zampini } 3523927de2eSStefano Zampini } 3533927de2eSStefano Zampini } 3543927de2eSStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 3553927de2eSStefano Zampini } 3563927de2eSStefano Zampini } 357ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 358ecf5a873SStefano Zampini if (global_indices_c != global_indices_r) { 3597230de76SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 360ecf5a873SStefano Zampini } 3613927de2eSStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 362ecf5a873SStefano Zampini 363ecf5a873SStefano Zampini /* Reduce my_dnz and my_onz */ 3643927de2eSStefano Zampini if (maxreduce) { 3653927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 3663927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 3673927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 3683927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 3693927de2eSStefano Zampini } else { 3703927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 3713927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 3723927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 3733927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 3743927de2eSStefano Zampini } 3753927de2eSStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 3763927de2eSStefano Zampini 3773927de2eSStefano Zampini /* Resize preallocation if overestimated */ 3783927de2eSStefano Zampini for (i=0;i<lrows;i++) { 3793927de2eSStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 3803927de2eSStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 3813927de2eSStefano Zampini } 3823927de2eSStefano Zampini /* set preallocation */ 3833927de2eSStefano Zampini ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 3843927de2eSStefano Zampini for (i=0;i<lrows/bs;i++) { 3853927de2eSStefano Zampini dnz[i] = dnz[i*bs]/bs; 3863927de2eSStefano Zampini onz[i] = onz[i*bs]/bs; 3873927de2eSStefano Zampini } 3883927de2eSStefano Zampini ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 3893927de2eSStefano Zampini ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 3903927de2eSStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 3913927de2eSStefano Zampini if (issbaij) { 3923927de2eSStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 3933927de2eSStefano Zampini } 3943927de2eSStefano Zampini PetscFunctionReturn(0); 3953927de2eSStefano Zampini } 3963927de2eSStefano Zampini 3973927de2eSStefano Zampini #undef __FUNCT__ 398b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 3997230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 400b7ce53b6SStefano Zampini { 401b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 402d9a9e74cSStefano Zampini Mat local_mat; 403b7ce53b6SStefano Zampini /* info on mat */ 4043cfa4ea4SStefano Zampini PetscInt bs,rows,cols,lrows,lcols; 405b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 406686e3a49SStefano Zampini PetscBool isdense,issbaij,isseqaij; 4077c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 408b7ce53b6SStefano Zampini /* values insertion */ 409b7ce53b6SStefano Zampini PetscScalar *array; 410b7ce53b6SStefano Zampini /* work */ 411b7ce53b6SStefano Zampini PetscErrorCode ierr; 412b7ce53b6SStefano Zampini 413b7ce53b6SStefano Zampini PetscFunctionBegin; 414b7ce53b6SStefano Zampini /* get info from mat */ 4157c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 4167c03b4e8SStefano Zampini if (nsubdomains == 1) { 4177c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 4187c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 4197c03b4e8SStefano Zampini } else { 4207c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 4217c03b4e8SStefano Zampini } 4227c03b4e8SStefano Zampini PetscFunctionReturn(0); 4237c03b4e8SStefano Zampini } 424b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 425b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 4263cfa4ea4SStefano Zampini ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 427b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 428b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 429686e3a49SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 430b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 431b7ce53b6SStefano Zampini 432b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 433b7ce53b6SStefano Zampini MatType new_mat_type; 4343927de2eSStefano Zampini PetscBool issbaij_red; 435b7ce53b6SStefano Zampini 436b7ce53b6SStefano Zampini /* determining new matrix type */ 437b2566f29SBarry Smith ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 438b7ce53b6SStefano Zampini if (issbaij_red) { 439b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 440b7ce53b6SStefano Zampini } else { 441b7ce53b6SStefano Zampini if (bs>1) { 442b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 443b7ce53b6SStefano Zampini } else { 444b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 445b7ce53b6SStefano Zampini } 446b7ce53b6SStefano Zampini } 447b7ce53b6SStefano Zampini 4483927de2eSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 4493cfa4ea4SStefano Zampini ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 4503927de2eSStefano Zampini ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 4513927de2eSStefano Zampini ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr); 4523927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 453b7ce53b6SStefano Zampini } else { 4543cfa4ea4SStefano Zampini PetscInt mbs,mrows,mcols,mlrows,mlcols; 455b7ce53b6SStefano Zampini /* some checks */ 456b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 457b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 4583cfa4ea4SStefano Zampini ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 4596c4ed002SBarry Smith if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 4606c4ed002SBarry Smith if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 4616c4ed002SBarry Smith if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 4626c4ed002SBarry Smith if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 4636c4ed002SBarry Smith if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 464b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 465b7ce53b6SStefano Zampini } 466d9a9e74cSStefano Zampini 467d9a9e74cSStefano Zampini if (issbaij) { 468d9a9e74cSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 469d9a9e74cSStefano Zampini } else { 470d9a9e74cSStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 471d9a9e74cSStefano Zampini local_mat = matis->A; 472d9a9e74cSStefano Zampini } 473686e3a49SStefano Zampini 474b7ce53b6SStefano Zampini /* Set values */ 475ecf5a873SStefano Zampini ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 476b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 477ecf5a873SStefano Zampini PetscInt i,*dummy_rows,*dummy_cols; 478ecf5a873SStefano Zampini 479ecf5a873SStefano Zampini if (local_rows != local_cols) { 480ecf5a873SStefano Zampini ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 481ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 482ecf5a873SStefano Zampini for (i=0;i<local_cols;i++) dummy_cols[i] = i; 483ecf5a873SStefano Zampini } else { 484ecf5a873SStefano Zampini ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 485ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 486ecf5a873SStefano Zampini dummy_cols = dummy_rows; 487ecf5a873SStefano Zampini } 488b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 489d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 490ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 491d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 492ecf5a873SStefano Zampini if (dummy_rows != dummy_cols) { 493ecf5a873SStefano Zampini ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 494ecf5a873SStefano Zampini } else { 495ecf5a873SStefano Zampini ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 496ecf5a873SStefano Zampini } 497686e3a49SStefano Zampini } else if (isseqaij) { 498ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 499686e3a49SStefano Zampini PetscBool done; 500686e3a49SStefano Zampini 501d9a9e74cSStefano Zampini ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 5026c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 503d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 504686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 505ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 506686e3a49SStefano Zampini } 507d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 5086c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 509d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 510686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 511ecf5a873SStefano Zampini PetscInt i; 512c0962df8SStefano Zampini 513686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 514686e3a49SStefano Zampini PetscInt j; 515ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 516686e3a49SStefano Zampini 517ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 518ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 519ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 520686e3a49SStefano Zampini } 521b7ce53b6SStefano Zampini } 522b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 523d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 524b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 525b7ce53b6SStefano Zampini if (isdense) { 526b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 527b7ce53b6SStefano Zampini } 528b7ce53b6SStefano Zampini PetscFunctionReturn(0); 529b7ce53b6SStefano Zampini } 530b7ce53b6SStefano Zampini 531b7ce53b6SStefano Zampini #undef __FUNCT__ 532b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 533b7ce53b6SStefano Zampini /*@ 534b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 535b7ce53b6SStefano Zampini 536b7ce53b6SStefano Zampini Input Parameter: 537b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 538b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 539b7ce53b6SStefano Zampini 540b7ce53b6SStefano Zampini Output Parameter: 541b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 542b7ce53b6SStefano Zampini 543b7ce53b6SStefano Zampini Level: developer 544b7ce53b6SStefano Zampini 545eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 546b7ce53b6SStefano Zampini 547b7ce53b6SStefano Zampini .seealso: MATIS 548b7ce53b6SStefano Zampini @*/ 549b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 550b7ce53b6SStefano Zampini { 551b7ce53b6SStefano Zampini PetscErrorCode ierr; 552b7ce53b6SStefano Zampini 553b7ce53b6SStefano Zampini PetscFunctionBegin; 554b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 555b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 556b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 557b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 558b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 559b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 5606c4ed002SBarry Smith if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 561b7ce53b6SStefano Zampini } 562b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 563b7ce53b6SStefano Zampini PetscFunctionReturn(0); 564b7ce53b6SStefano Zampini } 565b7ce53b6SStefano Zampini 566b7ce53b6SStefano Zampini #undef __FUNCT__ 567ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 568ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 569ad6194a2SStefano Zampini { 570ad6194a2SStefano Zampini PetscErrorCode ierr; 571ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 572ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 573ad6194a2SStefano Zampini Mat B,localmat; 574ad6194a2SStefano Zampini 575ad6194a2SStefano Zampini PetscFunctionBegin; 576ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 577ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 578ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 579e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 580ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 581ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 582b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 583ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 584ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 585ad6194a2SStefano Zampini *newmat = B; 586ad6194a2SStefano Zampini PetscFunctionReturn(0); 587ad6194a2SStefano Zampini } 588ad6194a2SStefano Zampini 589ad6194a2SStefano Zampini #undef __FUNCT__ 59069796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 59169796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 59269796d55SStefano Zampini { 59369796d55SStefano Zampini PetscErrorCode ierr; 59469796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 59569796d55SStefano Zampini PetscBool local_sym; 59669796d55SStefano Zampini 59769796d55SStefano Zampini PetscFunctionBegin; 59869796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 599b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 60069796d55SStefano Zampini PetscFunctionReturn(0); 60169796d55SStefano Zampini } 60269796d55SStefano Zampini 60369796d55SStefano Zampini #undef __FUNCT__ 60469796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 60569796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 60669796d55SStefano Zampini { 60769796d55SStefano Zampini PetscErrorCode ierr; 60869796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 60969796d55SStefano Zampini PetscBool local_sym; 61069796d55SStefano Zampini 61169796d55SStefano Zampini PetscFunctionBegin; 61269796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 613b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 61469796d55SStefano Zampini PetscFunctionReturn(0); 61569796d55SStefano Zampini } 61669796d55SStefano Zampini 61769796d55SStefano Zampini #undef __FUNCT__ 618b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 619dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 620b4319ba4SBarry Smith { 621dfbe8321SBarry Smith PetscErrorCode ierr; 622b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 623b4319ba4SBarry Smith 624b4319ba4SBarry Smith PetscFunctionBegin; 6256bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 626e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 627e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 6286bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 6296bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 63028f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 63128f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 632bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 633dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 634bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 635b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 636b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 6372e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 638b4319ba4SBarry Smith PetscFunctionReturn(0); 639b4319ba4SBarry Smith } 640b4319ba4SBarry Smith 641b4319ba4SBarry Smith #undef __FUNCT__ 642b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 643dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 644b4319ba4SBarry Smith { 645dfbe8321SBarry Smith PetscErrorCode ierr; 646b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 647b4319ba4SBarry Smith PetscScalar zero = 0.0; 648b4319ba4SBarry Smith 649b4319ba4SBarry Smith PetscFunctionBegin; 650b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 651e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 652e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 653b4319ba4SBarry Smith 654b4319ba4SBarry Smith /* multiply the local matrix */ 655b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 656b4319ba4SBarry Smith 657b4319ba4SBarry Smith /* scatter product back into global memory */ 6582dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 659e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 660e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 661b4319ba4SBarry Smith PetscFunctionReturn(0); 662b4319ba4SBarry Smith } 663b4319ba4SBarry Smith 664b4319ba4SBarry Smith #undef __FUNCT__ 6652e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 666b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 6672e74eeadSLisandro Dalcin { 668650997f4SStefano Zampini Vec temp_vec; 6692e74eeadSLisandro Dalcin PetscErrorCode ierr; 6702e74eeadSLisandro Dalcin 6712e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 672650997f4SStefano Zampini if (v3 != v2) { 673650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 674650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 675650997f4SStefano Zampini } else { 676650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 677650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 678650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 679650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 680650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 681650997f4SStefano Zampini } 6822e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6832e74eeadSLisandro Dalcin } 6842e74eeadSLisandro Dalcin 6852e74eeadSLisandro Dalcin #undef __FUNCT__ 6862e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 687e176bc59SStefano Zampini PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 6882e74eeadSLisandro Dalcin { 6892e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 6902e74eeadSLisandro Dalcin PetscErrorCode ierr; 6912e74eeadSLisandro Dalcin 692e176bc59SStefano Zampini PetscFunctionBegin; 6932e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 694e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 695e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6962e74eeadSLisandro Dalcin 6972e74eeadSLisandro Dalcin /* multiply the local matrix */ 698e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 6992e74eeadSLisandro Dalcin 7002e74eeadSLisandro Dalcin /* scatter product back into global vector */ 701e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 702e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 703e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 7042e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7052e74eeadSLisandro Dalcin } 7062e74eeadSLisandro Dalcin 7072e74eeadSLisandro Dalcin #undef __FUNCT__ 7082e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 7092e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 7102e74eeadSLisandro Dalcin { 711650997f4SStefano Zampini Vec temp_vec; 7122e74eeadSLisandro Dalcin PetscErrorCode ierr; 7132e74eeadSLisandro Dalcin 7142e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 715650997f4SStefano Zampini if (v3 != v2) { 716650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 717650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 718650997f4SStefano Zampini } else { 719650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 720650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 721650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 722650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 723650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 724650997f4SStefano Zampini } 7252e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7262e74eeadSLisandro Dalcin } 7272e74eeadSLisandro Dalcin 7282e74eeadSLisandro Dalcin #undef __FUNCT__ 729b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 730dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 731b4319ba4SBarry Smith { 732b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 733dfbe8321SBarry Smith PetscErrorCode ierr; 734b4319ba4SBarry Smith PetscViewer sviewer; 735b4319ba4SBarry Smith 736b4319ba4SBarry Smith PetscFunctionBegin; 7373f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 738b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 7393f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 7406e520ac8SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 741b4319ba4SBarry Smith PetscFunctionReturn(0); 742b4319ba4SBarry Smith } 743b4319ba4SBarry Smith 744b4319ba4SBarry Smith #undef __FUNCT__ 745b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 746784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 747b4319ba4SBarry Smith { 748dfbe8321SBarry Smith PetscErrorCode ierr; 749e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 750b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 751b4319ba4SBarry Smith IS from,to; 752e176bc59SStefano Zampini Vec cglobal,rglobal; 753b4319ba4SBarry Smith 754b4319ba4SBarry Smith PetscFunctionBegin; 755784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 756e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 7573bbff08aSStefano Zampini /* Destroy any previous data */ 75870cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 75970cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 760e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 761e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 7621c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 76328f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 76428f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 7653bbff08aSStefano Zampini 7663bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 767fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 768fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 769fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 770fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 771b4319ba4SBarry Smith 772b4319ba4SBarry Smith /* Create the local matrix A */ 773e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 774e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 775e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 776e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 777f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 778e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 779e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 780e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 781ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 782ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 783b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 784b4319ba4SBarry Smith 785b4319ba4SBarry Smith /* Create the local work vectors */ 7863bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 787b4319ba4SBarry Smith 788e176bc59SStefano Zampini /* setup the global to local scatters */ 789e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 790e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 791784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 792e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 793e176bc59SStefano Zampini if (rmapping != cmapping) { 794e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 795e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 796e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 797e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 798e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 799e176bc59SStefano Zampini } else { 800e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 801e176bc59SStefano Zampini is->cctx = is->rctx; 802e176bc59SStefano Zampini } 803e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 804e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 8056bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 8066bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 80748ff6bf3SStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 808b4319ba4SBarry Smith PetscFunctionReturn(0); 809b4319ba4SBarry Smith } 810b4319ba4SBarry Smith 8112e74eeadSLisandro Dalcin #undef __FUNCT__ 8122e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 8132e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 8142e74eeadSLisandro Dalcin { 8152e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 8162e74eeadSLisandro Dalcin PetscErrorCode ierr; 81797563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 81897563a80SStefano Zampini PetscInt i,zm,zn; 81997563a80SStefano Zampini #endif 82097563a80SStefano Zampini PetscInt rows_l[2048],cols_l[2048]; 8212e74eeadSLisandro Dalcin 8222e74eeadSLisandro Dalcin PetscFunctionBegin; 8232e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 824f23aa3ddSBarry 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); 82597563a80SStefano Zampini /* count negative indices */ 82697563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 82797563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 8282e74eeadSLisandro Dalcin #endif 82997563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 83097563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 83197563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 83297563a80SStefano Zampini /* count negative indices (should be the same as before) */ 83397563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 83497563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 83597563a80SStefano 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"); 83697563a80SStefano 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"); 83797563a80SStefano Zampini #endif 8382e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 8392e74eeadSLisandro Dalcin PetscFunctionReturn(0); 8402e74eeadSLisandro Dalcin } 8412e74eeadSLisandro Dalcin 842b4319ba4SBarry Smith #undef __FUNCT__ 84397563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS" 84497563a80SStefano Zampini PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 84597563a80SStefano Zampini { 84697563a80SStefano Zampini Mat_IS *is = (Mat_IS*)mat->data; 84797563a80SStefano Zampini PetscErrorCode ierr; 84897563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 84997563a80SStefano Zampini PetscInt i,zm,zn; 85097563a80SStefano Zampini #endif 85197563a80SStefano Zampini PetscInt rows_l[2048],cols_l[2048]; 85297563a80SStefano Zampini 85397563a80SStefano Zampini PetscFunctionBegin; 85497563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 85597563a80SStefano 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); 85697563a80SStefano Zampini /* count negative indices */ 85797563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 85897563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 85997563a80SStefano Zampini #endif 86097563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 86197563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 86297563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 86397563a80SStefano Zampini /* count negative indices (should be the same as before) */ 86497563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 86597563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 86697563a80SStefano 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"); 86797563a80SStefano 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"); 86897563a80SStefano Zampini #endif 86997563a80SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 87097563a80SStefano Zampini PetscFunctionReturn(0); 87197563a80SStefano Zampini } 87297563a80SStefano Zampini 87397563a80SStefano Zampini #undef __FUNCT__ 874b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 87513f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 876b4319ba4SBarry Smith { 877dfbe8321SBarry Smith PetscErrorCode ierr; 878b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 879b4319ba4SBarry Smith 880b4319ba4SBarry Smith PetscFunctionBegin; 881b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 882b4319ba4SBarry Smith PetscFunctionReturn(0); 883b4319ba4SBarry Smith } 884b4319ba4SBarry Smith 885b4319ba4SBarry Smith #undef __FUNCT__ 886f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 887f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 888f0006bf2SLisandro Dalcin { 889f0006bf2SLisandro Dalcin PetscErrorCode ierr; 890f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 891f0006bf2SLisandro Dalcin 892f0006bf2SLisandro Dalcin PetscFunctionBegin; 893f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 894f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 895f0006bf2SLisandro Dalcin } 896f0006bf2SLisandro Dalcin 897f0006bf2SLisandro Dalcin #undef __FUNCT__ 8982e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 8992b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 9002e74eeadSLisandro Dalcin { 9016e520ac8SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 9026e520ac8SStefano Zampini PetscInt nr,nl,len,i; 9036e520ac8SStefano Zampini PetscInt *lrows; 9042e74eeadSLisandro Dalcin PetscErrorCode ierr; 9052e74eeadSLisandro Dalcin 9062e74eeadSLisandro Dalcin PetscFunctionBegin; 9076e520ac8SStefano Zampini /* get locally owned rows */ 9086e520ac8SStefano Zampini ierr = MatZeroRowsMapLocal_Private(A,n,rows,&len,&lrows);CHKERRQ(ierr); 9096e520ac8SStefano Zampini /* fix right hand side if needed */ 9106e520ac8SStefano Zampini if (x && b) { 9116e520ac8SStefano Zampini const PetscScalar *xx; 9126e520ac8SStefano Zampini PetscScalar *bb; 9136e520ac8SStefano Zampini 9146e520ac8SStefano Zampini ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr); 9156e520ac8SStefano Zampini ierr = VecGetArray(b, &bb);CHKERRQ(ierr); 9166e520ac8SStefano Zampini for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]]; 9176e520ac8SStefano Zampini ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr); 9186e520ac8SStefano Zampini ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr); 9192e74eeadSLisandro Dalcin } 9206e520ac8SStefano Zampini /* get rows associated to the local matrices */ 9216e520ac8SStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 9226e520ac8SStefano Zampini ierr = MatISComputeSF_Private(A);CHKERRQ(ierr); 9236e520ac8SStefano Zampini } 9246e520ac8SStefano Zampini ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr); 9256e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr); 9266e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 9276e520ac8SStefano Zampini for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1; 9286e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 9296e520ac8SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 9306e520ac8SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 9316e520ac8SStefano Zampini ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr); 9326e520ac8SStefano Zampini for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i; 9337230de76SStefano Zampini ierr = MatISZeroRowsLocal_Private(A,nr,lrows,diag);CHKERRQ(ierr); 9346e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 9352e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9362e74eeadSLisandro Dalcin } 9372e74eeadSLisandro Dalcin 9382e74eeadSLisandro Dalcin #undef __FUNCT__ 9397230de76SStefano Zampini #define __FUNCT__ "MatISZeroRowsLocal_Private" 9407230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag) 941b4319ba4SBarry Smith { 942b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 943dfbe8321SBarry Smith PetscErrorCode ierr; 944b4319ba4SBarry Smith 945b4319ba4SBarry Smith PetscFunctionBegin; 9466e520ac8SStefano Zampini if (diag != 0.) { 947b4319ba4SBarry Smith /* 948b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 949b4319ba4SBarry Smith work properly in the interface nodes. 950b4319ba4SBarry Smith */ 951b4319ba4SBarry Smith Vec counter; 952e176bc59SStefano Zampini ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr); 953e176bc59SStefano Zampini ierr = VecSet(counter,0.);CHKERRQ(ierr); 954e176bc59SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 955e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 956e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 957e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 958e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9596bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 960b4319ba4SBarry Smith } 9612b404112SStefano Zampini 962958c9bccSBarry Smith if (!n) { 963b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 964b4319ba4SBarry Smith } else { 9657230de76SStefano Zampini PetscInt i; 966b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 9672205254eSKarl Rupp 9682b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 9696e520ac8SStefano Zampini if (diag != 0.) { 9706e520ac8SStefano Zampini const PetscScalar *array; 9716e520ac8SStefano Zampini ierr = VecGetArrayRead(is->y,&array);CHKERRQ(ierr); 972b4319ba4SBarry Smith for (i=0; i<n; i++) { 973f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 974b4319ba4SBarry Smith } 9756e520ac8SStefano Zampini ierr = VecRestoreArrayRead(is->y,&array);CHKERRQ(ierr); 9766e520ac8SStefano Zampini } 977b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 978b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 979b4319ba4SBarry Smith } 980b4319ba4SBarry Smith PetscFunctionReturn(0); 981b4319ba4SBarry Smith } 982b4319ba4SBarry Smith 983b4319ba4SBarry Smith #undef __FUNCT__ 984b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 985dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 986b4319ba4SBarry Smith { 987b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 988dfbe8321SBarry Smith PetscErrorCode ierr; 989b4319ba4SBarry Smith 990b4319ba4SBarry Smith PetscFunctionBegin; 991b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 992b4319ba4SBarry Smith PetscFunctionReturn(0); 993b4319ba4SBarry Smith } 994b4319ba4SBarry Smith 995b4319ba4SBarry Smith #undef __FUNCT__ 996b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 997dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 998b4319ba4SBarry Smith { 999b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1000dfbe8321SBarry Smith PetscErrorCode ierr; 1001b4319ba4SBarry Smith 1002b4319ba4SBarry Smith PetscFunctionBegin; 1003b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 1004b4319ba4SBarry Smith PetscFunctionReturn(0); 1005b4319ba4SBarry Smith } 1006b4319ba4SBarry Smith 1007b4319ba4SBarry Smith #undef __FUNCT__ 1008b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 10097087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 1010b4319ba4SBarry Smith { 1011b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 1012b4319ba4SBarry Smith 1013b4319ba4SBarry Smith PetscFunctionBegin; 1014b4319ba4SBarry Smith *local = is->A; 1015b4319ba4SBarry Smith PetscFunctionReturn(0); 1016b4319ba4SBarry Smith } 1017b4319ba4SBarry Smith 1018b4319ba4SBarry Smith #undef __FUNCT__ 1019b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 1020b4319ba4SBarry Smith /*@ 1021b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 1022b4319ba4SBarry Smith 1023b4319ba4SBarry Smith Input Parameter: 1024b4319ba4SBarry Smith . mat - the matrix 1025b4319ba4SBarry Smith 1026b4319ba4SBarry Smith Output Parameter: 1027eb82efa4SStefano Zampini . local - the local matrix 1028b4319ba4SBarry Smith 1029b4319ba4SBarry Smith Level: advanced 1030b4319ba4SBarry Smith 1031b4319ba4SBarry Smith Notes: 1032b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 1033b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 1034b4319ba4SBarry Smith of the MatSetValues() operation. 1035b4319ba4SBarry Smith 1036b4319ba4SBarry Smith .seealso: MATIS 1037b4319ba4SBarry Smith @*/ 10387087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 1039b4319ba4SBarry Smith { 10404ac538c5SBarry Smith PetscErrorCode ierr; 1041b4319ba4SBarry Smith 1042b4319ba4SBarry Smith PetscFunctionBegin; 10430700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1044b4319ba4SBarry Smith PetscValidPointer(local,2); 10454ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 1046b4319ba4SBarry Smith PetscFunctionReturn(0); 1047b4319ba4SBarry Smith } 1048b4319ba4SBarry Smith 10493b03a366Sstefano_zampini #undef __FUNCT__ 10503b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 10513b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 10523b03a366Sstefano_zampini { 10533b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 10543b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 10553b03a366Sstefano_zampini PetscErrorCode ierr; 10563b03a366Sstefano_zampini 10573b03a366Sstefano_zampini PetscFunctionBegin; 10584e4c7dbeSStefano Zampini if (is->A) { 10593b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 10603b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 1061f23aa3ddSBarry 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); 10624e4c7dbeSStefano Zampini } 10633b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 10643b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 10653b03a366Sstefano_zampini is->A = local; 10663b03a366Sstefano_zampini PetscFunctionReturn(0); 10673b03a366Sstefano_zampini } 10683b03a366Sstefano_zampini 10693b03a366Sstefano_zampini #undef __FUNCT__ 10703b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 10713b03a366Sstefano_zampini /*@ 1072eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 10733b03a366Sstefano_zampini 10743b03a366Sstefano_zampini Input Parameter: 10753b03a366Sstefano_zampini . mat - the matrix 1076eb82efa4SStefano Zampini . local - the local matrix 10773b03a366Sstefano_zampini 10783b03a366Sstefano_zampini Output Parameter: 10793b03a366Sstefano_zampini 10803b03a366Sstefano_zampini Level: advanced 10813b03a366Sstefano_zampini 10823b03a366Sstefano_zampini Notes: 10833b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 10843b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 10853b03a366Sstefano_zampini 10863b03a366Sstefano_zampini .seealso: MATIS 10873b03a366Sstefano_zampini @*/ 10883b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 10893b03a366Sstefano_zampini { 10903b03a366Sstefano_zampini PetscErrorCode ierr; 10913b03a366Sstefano_zampini 10923b03a366Sstefano_zampini PetscFunctionBegin; 10933b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1094b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 10953b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 10963b03a366Sstefano_zampini PetscFunctionReturn(0); 10973b03a366Sstefano_zampini } 10983b03a366Sstefano_zampini 10996726f965SBarry Smith #undef __FUNCT__ 11006726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 11016726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 11026726f965SBarry Smith { 11036726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 11046726f965SBarry Smith PetscErrorCode ierr; 11056726f965SBarry Smith 11066726f965SBarry Smith PetscFunctionBegin; 11076726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 11086726f965SBarry Smith PetscFunctionReturn(0); 11096726f965SBarry Smith } 11106726f965SBarry Smith 11116726f965SBarry Smith #undef __FUNCT__ 11122e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 11132e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 11142e74eeadSLisandro Dalcin { 11152e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 11162e74eeadSLisandro Dalcin PetscErrorCode ierr; 11172e74eeadSLisandro Dalcin 11182e74eeadSLisandro Dalcin PetscFunctionBegin; 11192e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 11202e74eeadSLisandro Dalcin PetscFunctionReturn(0); 11212e74eeadSLisandro Dalcin } 11222e74eeadSLisandro Dalcin 11232e74eeadSLisandro Dalcin #undef __FUNCT__ 11242e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 11252e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 11262e74eeadSLisandro Dalcin { 11272e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 11282e74eeadSLisandro Dalcin PetscErrorCode ierr; 11292e74eeadSLisandro Dalcin 11302e74eeadSLisandro Dalcin PetscFunctionBegin; 11312e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 1132e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 11332e74eeadSLisandro Dalcin 11342e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 11352e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 1136e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1137e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11382e74eeadSLisandro Dalcin PetscFunctionReturn(0); 11392e74eeadSLisandro Dalcin } 11402e74eeadSLisandro Dalcin 11412e74eeadSLisandro Dalcin #undef __FUNCT__ 11426726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 1143ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 11446726f965SBarry Smith { 11456726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 11466726f965SBarry Smith PetscErrorCode ierr; 11476726f965SBarry Smith 11486726f965SBarry Smith PetscFunctionBegin; 11494e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 11506726f965SBarry Smith PetscFunctionReturn(0); 11516726f965SBarry Smith } 11526726f965SBarry Smith 1153284134d9SBarry Smith #undef __FUNCT__ 1154284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 1155284134d9SBarry Smith /*@ 11563c212e90SHong Zhang MatCreateIS - Creates a "process" unassembled matrix, assembled on each 1157284134d9SBarry Smith process but not across processes. 1158284134d9SBarry Smith 1159284134d9SBarry Smith Input Parameters: 1160284134d9SBarry Smith + comm - MPI communicator that will share the matrix 1161e176bc59SStefano Zampini . bs - block size of the matrix 1162df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 1163e176bc59SStefano Zampini . rmap - local to global map for rows 1164e176bc59SStefano Zampini - cmap - local to global map for cols 1165284134d9SBarry Smith 1166284134d9SBarry Smith Output Parameter: 1167284134d9SBarry Smith . A - the resulting matrix 1168284134d9SBarry Smith 11698e6c10adSSatish Balay Level: advanced 11708e6c10adSSatish Balay 11713c212e90SHong Zhang Notes: See MATIS for more details. 11723c212e90SHong Zhang m and n are NOT related to the size of the map; they are the size of the part of the vector owned 1173e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 11743c212e90SHong Zhang If either rmap or cmap are NULL, then the matrix is assumed to be square. 1175284134d9SBarry Smith 1176284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 1177284134d9SBarry Smith @*/ 1178e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 1179284134d9SBarry Smith { 1180284134d9SBarry Smith PetscErrorCode ierr; 1181284134d9SBarry Smith 1182284134d9SBarry Smith PetscFunctionBegin; 1183e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 1184284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1185d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 1186284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1187284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 1188e176bc59SStefano Zampini if (rmap && cmap) { 1189e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 1190e176bc59SStefano Zampini } else if (!rmap) { 1191e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 1192e176bc59SStefano Zampini } else { 1193e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 1194e176bc59SStefano Zampini } 1195284134d9SBarry Smith PetscFunctionReturn(0); 1196284134d9SBarry Smith } 1197284134d9SBarry Smith 1198b4319ba4SBarry Smith /*MC 1199eb82efa4SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC). 1200b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 1201b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 1202b4319ba4SBarry Smith product is handled "implicitly". 1203b4319ba4SBarry Smith 1204b4319ba4SBarry Smith Operations Provided: 12056726f965SBarry Smith + MatMult() 12062e74eeadSLisandro Dalcin . MatMultAdd() 12072e74eeadSLisandro Dalcin . MatMultTranspose() 12082e74eeadSLisandro Dalcin . MatMultTransposeAdd() 12096726f965SBarry Smith . MatZeroEntries() 12106726f965SBarry Smith . MatSetOption() 12112e74eeadSLisandro Dalcin . MatZeroRows() 12122e74eeadSLisandro Dalcin . MatSetValues() 121397563a80SStefano Zampini . MatSetValuesBlocked() 12146726f965SBarry Smith . MatSetValuesLocal() 121597563a80SStefano Zampini . MatSetValuesBlockedLocal() 12162e74eeadSLisandro Dalcin . MatScale() 12172e74eeadSLisandro Dalcin . MatGetDiagonal() 12182b404112SStefano Zampini . MatMissingDiagonal() 12192b404112SStefano Zampini . MatDuplicate() 12202b404112SStefano Zampini . MatCopy() 12216726f965SBarry Smith - MatSetLocalToGlobalMapping() 1222b4319ba4SBarry Smith 1223b4319ba4SBarry Smith Options Database Keys: 1224b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1225b4319ba4SBarry Smith 1226b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1227b4319ba4SBarry Smith 1228b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1229b4319ba4SBarry Smith 1230b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1231eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1232b4319ba4SBarry Smith 1233b4319ba4SBarry Smith Level: advanced 1234b4319ba4SBarry Smith 12353c212e90SHong Zhang .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC, MatCreateIS() 1236b4319ba4SBarry Smith 1237b4319ba4SBarry Smith M*/ 1238b4319ba4SBarry Smith 1239b4319ba4SBarry Smith #undef __FUNCT__ 1240b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 12418cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1242b4319ba4SBarry Smith { 1243dfbe8321SBarry Smith PetscErrorCode ierr; 1244b4319ba4SBarry Smith Mat_IS *b; 1245b4319ba4SBarry Smith 1246b4319ba4SBarry Smith PetscFunctionBegin; 1247b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1248b4319ba4SBarry Smith A->data = (void*)b; 1249b4319ba4SBarry Smith 1250e176bc59SStefano Zampini /* matrix ops */ 1251e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1252b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 12532e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 12542e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 12552e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1256b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1257b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 12582e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 125998921651SStefano Zampini A->ops->setvaluesblocked = MatSetValuesBlocked_IS; 1260b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1261f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 12622e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1263b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1264b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1265b4319ba4SBarry Smith A->ops->view = MatView_IS; 12666726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 12672e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 12682e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 12696726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 127069796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 127169796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1272ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 12736bd84002SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_IS; 12742b404112SStefano Zampini A->ops->copy = MatCopy_IS; 1275659959c5SStefano Zampini A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_IS; 1276b4319ba4SBarry Smith 1277b7ce53b6SStefano Zampini /* special MATIS functions */ 1278bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1279bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1280b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 12812e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 128217667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1283b4319ba4SBarry Smith PetscFunctionReturn(0); 1284b4319ba4SBarry Smith } 1285