1b4319ba4SBarry Smith /* 2b4319ba4SBarry Smith Creates a matrix class for using the Neumann-Neumann type preconditioners. 3b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 4b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 5b4319ba4SBarry Smith product is handled "implicitly". 6b4319ba4SBarry Smith 7b4319ba4SBarry Smith Currently this allows for only one subdomain per processor. 8b4319ba4SBarry Smith */ 9b4319ba4SBarry Smith 10c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/ 1128f4e0baSStefano Zampini 12f26d0771SStefano Zampini #define MATIS_MAX_ENTRIES_INSERTION 2048 13f26d0771SStefano Zampini 14a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat); 15a72627d2SStefano Zampini 1628f4e0baSStefano Zampini #undef __FUNCT__ 17*d7f69cd0SStefano Zampini #define __FUNCT__ "MatTranspose_IS" 18*d7f69cd0SStefano Zampini PetscErrorCode MatTranspose_IS(Mat A,MatReuse reuse,Mat *B) 19*d7f69cd0SStefano Zampini { 20*d7f69cd0SStefano Zampini Mat C,lC,lA; 21*d7f69cd0SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 22*d7f69cd0SStefano Zampini PetscBool notset = PETSC_FALSE,cong = PETSC_TRUE; 23*d7f69cd0SStefano Zampini PetscErrorCode ierr; 24*d7f69cd0SStefano Zampini 25*d7f69cd0SStefano Zampini PetscFunctionBegin; 26*d7f69cd0SStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 27*d7f69cd0SStefano Zampini PetscBool rcong,ccong; 28*d7f69cd0SStefano Zampini 29*d7f69cd0SStefano Zampini ierr = PetscLayoutCompare((*B)->cmap,A->rmap,&rcong);CHKERRQ(ierr); 30*d7f69cd0SStefano Zampini ierr = PetscLayoutCompare((*B)->rmap,A->cmap,&ccong);CHKERRQ(ierr); 31*d7f69cd0SStefano Zampini cong = rcong || ccong; 32*d7f69cd0SStefano Zampini } 33*d7f69cd0SStefano Zampini if (reuse == MAT_INITIAL_MATRIX || *B == A || !cong) { 34*d7f69cd0SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr); 35*d7f69cd0SStefano Zampini } else { 36*d7f69cd0SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*B,MATIS,¬set);CHKERRQ(ierr); 37*d7f69cd0SStefano Zampini C = *B; 38*d7f69cd0SStefano Zampini } 39*d7f69cd0SStefano Zampini 40*d7f69cd0SStefano Zampini if (!notset) { 41*d7f69cd0SStefano Zampini ierr = MatSetSizes(C,A->cmap->n,A->rmap->n,A->cmap->N,A->rmap->N);CHKERRQ(ierr); 42*d7f69cd0SStefano Zampini ierr = MatSetBlockSizes(C,PetscAbs(A->cmap->bs),PetscAbs(A->rmap->bs));CHKERRQ(ierr); 43*d7f69cd0SStefano Zampini ierr = MatSetType(C,MATIS);CHKERRQ(ierr); 44*d7f69cd0SStefano Zampini ierr = MatGetLocalToGlobalMapping(A,&rl2g,&cl2g);CHKERRQ(ierr); 45*d7f69cd0SStefano Zampini ierr = MatSetLocalToGlobalMapping(C,cl2g,rl2g);CHKERRQ(ierr); 46*d7f69cd0SStefano Zampini } 47*d7f69cd0SStefano Zampini 48*d7f69cd0SStefano Zampini /* perform local transposition */ 49*d7f69cd0SStefano Zampini ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr); 50*d7f69cd0SStefano Zampini ierr = MatTranspose(lA,MAT_INITIAL_MATRIX,&lC);CHKERRQ(ierr); 51*d7f69cd0SStefano Zampini ierr = MatISSetLocalMat(C,lC);CHKERRQ(ierr); 52*d7f69cd0SStefano Zampini ierr = MatDestroy(&lC);CHKERRQ(ierr); 53*d7f69cd0SStefano Zampini 54*d7f69cd0SStefano Zampini ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 55*d7f69cd0SStefano Zampini ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 56*d7f69cd0SStefano Zampini 57*d7f69cd0SStefano Zampini if (reuse == MAT_INITIAL_MATRIX || *B != A) { 58*d7f69cd0SStefano Zampini if (!cong) { 59*d7f69cd0SStefano Zampini ierr = MatDestroy(B);CHKERRQ(ierr); 60*d7f69cd0SStefano Zampini } 61*d7f69cd0SStefano Zampini *B = C; 62*d7f69cd0SStefano Zampini } else { 63*d7f69cd0SStefano Zampini ierr = MatHeaderMerge(A,&C);CHKERRQ(ierr); 64*d7f69cd0SStefano Zampini } 65*d7f69cd0SStefano Zampini PetscFunctionReturn(0); 66*d7f69cd0SStefano Zampini } 67*d7f69cd0SStefano Zampini 68*d7f69cd0SStefano Zampini #undef __FUNCT__ 693fd1c9e7SStefano Zampini #define __FUNCT__ "MatDiagonalSet_IS" 703fd1c9e7SStefano Zampini PetscErrorCode MatDiagonalSet_IS(Mat A,Vec D,InsertMode insmode) 713fd1c9e7SStefano Zampini { 723fd1c9e7SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 733fd1c9e7SStefano Zampini PetscErrorCode ierr; 743fd1c9e7SStefano Zampini 753fd1c9e7SStefano Zampini PetscFunctionBegin; 763fd1c9e7SStefano Zampini if (!D) { /* this code branch is used by MatShift_IS */ 773fd1c9e7SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 783fd1c9e7SStefano Zampini } else { 793fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 803fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 813fd1c9e7SStefano Zampini } 823fd1c9e7SStefano Zampini ierr = VecPointwiseDivide(is->y,is->y,is->counter);CHKERRQ(ierr); 833fd1c9e7SStefano Zampini ierr = MatDiagonalSet(is->A,is->y,insmode);CHKERRQ(ierr); 843fd1c9e7SStefano Zampini PetscFunctionReturn(0); 853fd1c9e7SStefano Zampini } 863fd1c9e7SStefano Zampini 873fd1c9e7SStefano Zampini #undef __FUNCT__ 883fd1c9e7SStefano Zampini #define __FUNCT__ "MatShift_IS" 893fd1c9e7SStefano Zampini PetscErrorCode MatShift_IS(Mat A,PetscScalar a) 903fd1c9e7SStefano Zampini { 913fd1c9e7SStefano Zampini PetscErrorCode ierr; 923fd1c9e7SStefano Zampini 933fd1c9e7SStefano Zampini PetscFunctionBegin; 943fd1c9e7SStefano Zampini ierr = MatDiagonalSet_IS(A,NULL,ADD_VALUES);CHKERRQ(ierr); 953fd1c9e7SStefano Zampini PetscFunctionReturn(0); 963fd1c9e7SStefano Zampini } 973fd1c9e7SStefano Zampini 983fd1c9e7SStefano Zampini #undef __FUNCT__ 99f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesLocal_SubMat_IS" 100f26d0771SStefano Zampini static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 101f26d0771SStefano Zampini { 102f26d0771SStefano Zampini PetscErrorCode ierr; 103f26d0771SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 104f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 105f26d0771SStefano Zampini 106f26d0771SStefano Zampini PetscFunctionBegin; 107f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 108f26d0771SStefano Zampini if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n); 109f26d0771SStefano Zampini #endif 110f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 111f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 112f26d0771SStefano Zampini ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 113f26d0771SStefano Zampini PetscFunctionReturn(0); 114f26d0771SStefano Zampini } 115f26d0771SStefano Zampini 116f26d0771SStefano Zampini #undef __FUNCT__ 117f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesBlockedLocal_SubMat_IS" 118f26d0771SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 119f26d0771SStefano Zampini { 120f26d0771SStefano Zampini PetscErrorCode ierr; 121f26d0771SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 122f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 123f26d0771SStefano Zampini 124f26d0771SStefano Zampini PetscFunctionBegin; 125f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 126f26d0771SStefano Zampini if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column block indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n); 127f26d0771SStefano Zampini #endif 128f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 129f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 130f26d0771SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 131f26d0771SStefano Zampini PetscFunctionReturn(0); 132f26d0771SStefano Zampini } 133f26d0771SStefano Zampini 134f26d0771SStefano Zampini #undef __FUNCT__ 135a8116848SStefano Zampini #define __FUNCT__ "PetscLayoutMapLocal_Private" 136f0ae7da4SStefano Zampini static PetscErrorCode PetscLayoutMapLocal_Private(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs) 137a8116848SStefano Zampini { 138a8116848SStefano Zampini PetscInt *owners = map->range; 139a8116848SStefano Zampini PetscInt n = map->n; 140a8116848SStefano Zampini PetscSF sf; 141a8116848SStefano Zampini PetscInt *lidxs,*work = NULL; 142a8116848SStefano Zampini PetscSFNode *ridxs; 143a8116848SStefano Zampini PetscMPIInt rank; 144a8116848SStefano Zampini PetscInt r, p = 0, len = 0; 145a8116848SStefano Zampini PetscErrorCode ierr; 146a8116848SStefano Zampini 147a8116848SStefano Zampini PetscFunctionBegin; 148a8116848SStefano Zampini /* Create SF where leaves are input idxs and roots are owned idxs (code adapted from MatZeroRowsMapLocal_Private) */ 149a8116848SStefano Zampini ierr = MPI_Comm_rank(map->comm,&rank);CHKERRQ(ierr); 150a8116848SStefano Zampini ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr); 151a8116848SStefano Zampini for (r = 0; r < n; ++r) lidxs[r] = -1; 152a8116848SStefano Zampini ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr); 153a8116848SStefano Zampini for (r = 0; r < N; ++r) { 154a8116848SStefano Zampini const PetscInt idx = idxs[r]; 155a8116848SStefano Zampini if (idx < 0 || map->N <= idx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D out of range [0,%D)",idx,map->N); 156a8116848SStefano Zampini if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 157a8116848SStefano Zampini ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr); 158a8116848SStefano Zampini } 159a8116848SStefano Zampini ridxs[r].rank = p; 160a8116848SStefano Zampini ridxs[r].index = idxs[r] - owners[p]; 161a8116848SStefano Zampini } 162a8116848SStefano Zampini ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr); 163a8116848SStefano Zampini ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr); 164a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 165a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 166f0ae7da4SStefano Zampini if (ogidxs) { /* communicate global idxs */ 167a8116848SStefano Zampini PetscInt cum = 0,start,*work2; 168f0ae7da4SStefano Zampini 169f0ae7da4SStefano Zampini ierr = PetscMalloc1(n,&work);CHKERRQ(ierr); 170a8116848SStefano Zampini ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr); 171a8116848SStefano Zampini for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++; 172a8116848SStefano Zampini ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRQ(ierr); 173a8116848SStefano Zampini start -= cum; 174a8116848SStefano Zampini cum = 0; 175a8116848SStefano Zampini for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++; 176a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 177a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 178a8116848SStefano Zampini ierr = PetscFree(work2);CHKERRQ(ierr); 179a8116848SStefano Zampini } 180a8116848SStefano Zampini ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 181a8116848SStefano Zampini /* Compress and put in indices */ 182a8116848SStefano Zampini for (r = 0; r < n; ++r) 183a8116848SStefano Zampini if (lidxs[r] >= 0) { 184a8116848SStefano Zampini if (work) work[len] = work[r]; 185a8116848SStefano Zampini lidxs[len++] = r; 186a8116848SStefano Zampini } 187a8116848SStefano Zampini if (on) *on = len; 188a8116848SStefano Zampini if (oidxs) *oidxs = lidxs; 189a8116848SStefano Zampini if (ogidxs) *ogidxs = work; 190a8116848SStefano Zampini PetscFunctionReturn(0); 191a8116848SStefano Zampini } 192a8116848SStefano Zampini 193a8116848SStefano Zampini #undef __FUNCT__ 194a8116848SStefano Zampini #define __FUNCT__ "MatGetSubMatrix_IS" 195a8116848SStefano Zampini static PetscErrorCode MatGetSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat) 196a8116848SStefano Zampini { 197a8116848SStefano Zampini Mat locmat,newlocmat; 198a8116848SStefano Zampini Mat_IS *newmatis; 199a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG) 200a8116848SStefano Zampini Vec rtest,ltest; 201a8116848SStefano Zampini const PetscScalar *array; 202a8116848SStefano Zampini #endif 203a8116848SStefano Zampini const PetscInt *idxs; 204a8116848SStefano Zampini PetscInt i,m,n; 205a8116848SStefano Zampini PetscErrorCode ierr; 206a8116848SStefano Zampini 207a8116848SStefano Zampini PetscFunctionBegin; 208a8116848SStefano Zampini if (scall == MAT_REUSE_MATRIX) { 209a8116848SStefano Zampini PetscBool ismatis; 210a8116848SStefano Zampini 211a8116848SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr); 212a8116848SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type"); 213a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 214a8116848SStefano Zampini if (!newmatis->getsub_ris) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS"); 215a8116848SStefano Zampini if (!newmatis->getsub_cis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS"); 216a8116848SStefano Zampini } 217a8116848SStefano Zampini /* irow and icol may not have duplicate entries */ 218a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG) 219a8116848SStefano Zampini ierr = MatCreateVecs(mat,<est,&rtest);CHKERRQ(ierr); 220a8116848SStefano Zampini ierr = ISGetLocalSize(irow,&n);CHKERRQ(ierr); 221a8116848SStefano Zampini ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 222a8116848SStefano Zampini for (i=0;i<n;i++) { 223a8116848SStefano Zampini ierr = VecSetValue(rtest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 224a8116848SStefano Zampini } 225a8116848SStefano Zampini ierr = VecAssemblyBegin(rtest);CHKERRQ(ierr); 226a8116848SStefano Zampini ierr = VecAssemblyEnd(rtest);CHKERRQ(ierr); 227a8116848SStefano Zampini ierr = VecGetLocalSize(rtest,&n);CHKERRQ(ierr); 228a8116848SStefano Zampini ierr = VecGetOwnershipRange(rtest,&m,NULL);CHKERRQ(ierr); 229a8116848SStefano Zampini ierr = VecGetArrayRead(rtest,&array);CHKERRQ(ierr); 230fd479f66SMatthew G. Knepley for (i=0;i<n;i++) if (array[i] != 0. && array[i] != 1.) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Index %D counted %D times! Irow may not have duplicate entries",i+m,(PetscInt)PetscRealPart(array[i])); 231a8116848SStefano Zampini ierr = VecRestoreArrayRead(rtest,&array);CHKERRQ(ierr); 232a8116848SStefano Zampini ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 233a8116848SStefano Zampini ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 234a8116848SStefano Zampini ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 235a8116848SStefano Zampini for (i=0;i<n;i++) { 236a8116848SStefano Zampini ierr = VecSetValue(ltest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 237a8116848SStefano Zampini } 238a8116848SStefano Zampini ierr = VecAssemblyBegin(ltest);CHKERRQ(ierr); 239a8116848SStefano Zampini ierr = VecAssemblyEnd(ltest);CHKERRQ(ierr); 240a8116848SStefano Zampini ierr = VecGetLocalSize(ltest,&n);CHKERRQ(ierr); 241a8116848SStefano Zampini ierr = VecGetOwnershipRange(ltest,&m,NULL);CHKERRQ(ierr); 242a8116848SStefano Zampini ierr = VecGetArrayRead(ltest,&array);CHKERRQ(ierr); 243fd479f66SMatthew G. Knepley for (i=0;i<n;i++) if (array[i] != 0. && array[i] != 1.) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Index %D counted %D times! Icol may not have duplicate entries",i+m,(PetscInt)PetscRealPart(array[i])); 244a8116848SStefano Zampini ierr = VecRestoreArrayRead(ltest,&array);CHKERRQ(ierr); 245a8116848SStefano Zampini ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 246a8116848SStefano Zampini ierr = VecDestroy(&rtest);CHKERRQ(ierr); 247a8116848SStefano Zampini ierr = VecDestroy(<est);CHKERRQ(ierr); 248a8116848SStefano Zampini #endif 249a8116848SStefano Zampini if (scall == MAT_INITIAL_MATRIX) { 250a8116848SStefano Zampini Mat_IS *matis = (Mat_IS*)mat->data; 251a8116848SStefano Zampini ISLocalToGlobalMapping rl2g; 252a8116848SStefano Zampini IS is; 253a8116848SStefano Zampini PetscInt *lidxs,*lgidxs,*newgidxs; 2546dd40735SStefano Zampini PetscInt ll,newloc; 255a8116848SStefano Zampini MPI_Comm comm; 256a8116848SStefano Zampini 257a8116848SStefano Zampini ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 258a8116848SStefano Zampini ierr = ISGetLocalSize(irow,&m);CHKERRQ(ierr); 259a8116848SStefano Zampini ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 260a8116848SStefano Zampini ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 261a8116848SStefano Zampini ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr); 262a8116848SStefano Zampini ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 263a8116848SStefano Zampini /* communicate irow to their owners in the layout */ 264a8116848SStefano Zampini ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 265f0ae7da4SStefano Zampini ierr = PetscLayoutMapLocal_Private(mat->rmap,m,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 266a8116848SStefano Zampini ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 267a8116848SStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 268a8116848SStefano Zampini ierr = MatISComputeSF_Private(mat);CHKERRQ(ierr); 269a8116848SStefano Zampini } 270a8116848SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,matis->sf_nroots*sizeof(PetscInt));CHKERRQ(ierr); 271a8116848SStefano Zampini for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1; 272a8116848SStefano Zampini ierr = PetscFree(lidxs);CHKERRQ(ierr); 273a8116848SStefano Zampini ierr = PetscFree(lgidxs);CHKERRQ(ierr); 274a8116848SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 275a8116848SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 276a8116848SStefano Zampini for (i=0,newloc=0;i<matis->sf_nleaves;i++) if (matis->sf_leafdata[i]) newloc++; 277a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 278a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 279a8116848SStefano Zampini for (i=0,newloc=0;i<matis->sf_nleaves;i++) 280a8116848SStefano Zampini if (matis->sf_leafdata[i]) { 281a8116848SStefano Zampini lidxs[newloc] = i; 282a8116848SStefano Zampini newgidxs[newloc++] = matis->sf_leafdata[i]-1; 283a8116848SStefano Zampini } 284a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 285a8116848SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 286a8116848SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 287a8116848SStefano Zampini /* local is to extract local submatrix */ 288a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 289a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris);CHKERRQ(ierr); 290a8116848SStefano Zampini if (mat->congruentlayouts == -1) { /* first time we compare rows and cols layouts */ 291a8116848SStefano Zampini PetscBool cong; 292a8116848SStefano Zampini ierr = PetscLayoutCompare(mat->rmap,mat->cmap,&cong);CHKERRQ(ierr); 293a8116848SStefano Zampini if (cong) mat->congruentlayouts = 1; 294a8116848SStefano Zampini else mat->congruentlayouts = 0; 295a8116848SStefano Zampini } 296a8116848SStefano Zampini if (mat->congruentlayouts && irow == icol) { 297a8116848SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g);CHKERRQ(ierr); 298a8116848SStefano Zampini ierr = PetscObjectReference((PetscObject)newmatis->getsub_ris);CHKERRQ(ierr); 299a8116848SStefano Zampini newmatis->getsub_cis = newmatis->getsub_ris; 300a8116848SStefano Zampini } else { 301a8116848SStefano Zampini ISLocalToGlobalMapping cl2g; 302a8116848SStefano Zampini 303a8116848SStefano Zampini /* communicate icol to their owners in the layout */ 304a8116848SStefano Zampini ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 305f0ae7da4SStefano Zampini ierr = PetscLayoutMapLocal_Private(mat->cmap,n,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 306a8116848SStefano Zampini ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 307a8116848SStefano Zampini ierr = PetscMemzero(matis->csf_rootdata,matis->csf_nroots*sizeof(PetscInt));CHKERRQ(ierr); 308a8116848SStefano Zampini for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1; 309a8116848SStefano Zampini ierr = PetscFree(lidxs);CHKERRQ(ierr); 310a8116848SStefano Zampini ierr = PetscFree(lgidxs);CHKERRQ(ierr); 311a8116848SStefano Zampini ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 312a8116848SStefano Zampini ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 313a8116848SStefano Zampini for (i=0,newloc=0;i<matis->csf_nleaves;i++) if (matis->csf_leafdata[i]) newloc++; 314a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 315a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 316a8116848SStefano Zampini for (i=0,newloc=0;i<matis->csf_nleaves;i++) 317a8116848SStefano Zampini if (matis->csf_leafdata[i]) { 318a8116848SStefano Zampini lidxs[newloc] = i; 319a8116848SStefano Zampini newgidxs[newloc++] = matis->csf_leafdata[i]-1; 320a8116848SStefano Zampini } 321a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 322a8116848SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 323a8116848SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 324a8116848SStefano Zampini /* local is to extract local submatrix */ 325a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis);CHKERRQ(ierr); 326a8116848SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr); 327a8116848SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 328a8116848SStefano Zampini } 329a8116848SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 330a8116848SStefano Zampini } else { 331a8116848SStefano Zampini ierr = MatISGetLocalMat(*newmat,&newlocmat);CHKERRQ(ierr); 332a8116848SStefano Zampini } 333a8116848SStefano Zampini ierr = MatISGetLocalMat(mat,&locmat);CHKERRQ(ierr); 334a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 335a8116848SStefano Zampini ierr = MatGetSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat);CHKERRQ(ierr); 336a8116848SStefano Zampini if (scall == MAT_INITIAL_MATRIX) { 337a8116848SStefano Zampini ierr = MatISSetLocalMat(*newmat,newlocmat);CHKERRQ(ierr); 338a8116848SStefano Zampini ierr = MatDestroy(&newlocmat);CHKERRQ(ierr); 339a8116848SStefano Zampini } 340a8116848SStefano Zampini ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 341a8116848SStefano Zampini ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 342a8116848SStefano Zampini PetscFunctionReturn(0); 343a8116848SStefano Zampini } 344a8116848SStefano Zampini 345a8116848SStefano Zampini #undef __FUNCT__ 3462b404112SStefano Zampini #define __FUNCT__ "MatCopy_IS" 347a8116848SStefano Zampini static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str) 3482b404112SStefano Zampini { 3492b404112SStefano Zampini Mat_IS *a = (Mat_IS*)A->data,*b; 3502b404112SStefano Zampini PetscBool ismatis; 3512b404112SStefano Zampini PetscErrorCode ierr; 3522b404112SStefano Zampini 3532b404112SStefano Zampini PetscFunctionBegin; 3542b404112SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr); 3552b404112SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented"); 3562b404112SStefano Zampini b = (Mat_IS*)B->data; 3572b404112SStefano Zampini ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr); 3582b404112SStefano Zampini PetscFunctionReturn(0); 3592b404112SStefano Zampini } 3602b404112SStefano Zampini 3612b404112SStefano Zampini #undef __FUNCT__ 3626bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS" 363a8116848SStefano Zampini static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool *missing,PetscInt *d) 3646bd84002SStefano Zampini { 365527b2640SStefano Zampini Vec v; 366527b2640SStefano Zampini const PetscScalar *array; 367527b2640SStefano Zampini PetscInt i,n; 3686bd84002SStefano Zampini PetscErrorCode ierr; 3696bd84002SStefano Zampini 3706bd84002SStefano Zampini PetscFunctionBegin; 371527b2640SStefano Zampini *missing = PETSC_FALSE; 372527b2640SStefano Zampini ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr); 373527b2640SStefano Zampini ierr = MatGetDiagonal(A,v);CHKERRQ(ierr); 374527b2640SStefano Zampini ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 375527b2640SStefano Zampini ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr); 376527b2640SStefano Zampini for (i=0;i<n;i++) if (array[i] == 0.) break; 377527b2640SStefano Zampini ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr); 378527b2640SStefano Zampini ierr = VecDestroy(&v);CHKERRQ(ierr); 379527b2640SStefano Zampini if (i != n) *missing = PETSC_TRUE; 380527b2640SStefano Zampini if (d) { 381527b2640SStefano Zampini *d = -1; 382527b2640SStefano Zampini if (*missing) { 383527b2640SStefano Zampini PetscInt rstart; 384527b2640SStefano Zampini ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 385527b2640SStefano Zampini *d = i+rstart; 386527b2640SStefano Zampini } 387527b2640SStefano Zampini } 3886bd84002SStefano Zampini PetscFunctionReturn(0); 3896bd84002SStefano Zampini } 3906bd84002SStefano Zampini 3916bd84002SStefano Zampini #undef __FUNCT__ 39228f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private" 393a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B) 39428f4e0baSStefano Zampini { 39528f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 39628f4e0baSStefano Zampini const PetscInt *gidxs; 39728f4e0baSStefano Zampini PetscErrorCode ierr; 39828f4e0baSStefano Zampini 39928f4e0baSStefano Zampini PetscFunctionBegin; 40028f4e0baSStefano Zampini ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr); 40128f4e0baSStefano Zampini ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr); 40228f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 4033bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 40428f4e0baSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 4053bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 40628f4e0baSStefano Zampini ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 407a8116848SStefano Zampini if (B->rmap->mapping != B->cmap->mapping) { /* setup SF for columns */ 408a8116848SStefano Zampini ierr = MatGetSize(matis->A,NULL,&matis->csf_nleaves);CHKERRQ(ierr); 409a8116848SStefano Zampini ierr = MatGetLocalSize(B,NULL,&matis->csf_nroots);CHKERRQ(ierr); 410a8116848SStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf);CHKERRQ(ierr); 411a8116848SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 412a8116848SStefano Zampini ierr = PetscSFSetGraphLayout(matis->csf,B->cmap,matis->csf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 413a8116848SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 414a8116848SStefano Zampini ierr = PetscMalloc2(matis->csf_nroots,&matis->csf_rootdata,matis->csf_nleaves,&matis->csf_leafdata);CHKERRQ(ierr); 415a8116848SStefano Zampini } else { 416a8116848SStefano Zampini matis->csf = matis->sf; 417a8116848SStefano Zampini matis->csf_nleaves = matis->sf_nleaves; 418a8116848SStefano Zampini matis->csf_nroots = matis->sf_nroots; 419a8116848SStefano Zampini matis->csf_leafdata = matis->sf_leafdata; 420a8116848SStefano Zampini matis->csf_rootdata = matis->sf_rootdata; 421a8116848SStefano Zampini } 42228f4e0baSStefano Zampini PetscFunctionReturn(0); 42328f4e0baSStefano Zampini } 4242e1947a5SStefano Zampini 4252e1947a5SStefano Zampini #undef __FUNCT__ 4262e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 427eb82efa4SStefano Zampini /*@ 428a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 429a88811baSStefano Zampini 430a88811baSStefano Zampini Collective on MPI_Comm 431a88811baSStefano Zampini 432a88811baSStefano Zampini Input Parameters: 433a88811baSStefano Zampini + B - the matrix 434a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 435a88811baSStefano Zampini (same value is used for all local rows) 436a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 437a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 438a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 439a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 440a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 441a88811baSStefano Zampini the diagonal entry even if it is zero. 442a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 443a88811baSStefano Zampini submatrix (same value is used for all local rows). 444a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 445a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 446a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 447a88811baSStefano Zampini structure. The size of this array is equal to the number 448a88811baSStefano Zampini of local rows, i.e 'm'. 449a88811baSStefano Zampini 450a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 451a88811baSStefano Zampini 452a88811baSStefano Zampini Level: intermediate 453a88811baSStefano Zampini 454a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 455a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 456a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 457a88811baSStefano Zampini 458a88811baSStefano Zampini .keywords: matrix 459a88811baSStefano Zampini 4603c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS 461a88811baSStefano Zampini @*/ 4622e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 4632e1947a5SStefano Zampini { 4642e1947a5SStefano Zampini PetscErrorCode ierr; 4652e1947a5SStefano Zampini 4662e1947a5SStefano Zampini PetscFunctionBegin; 4672e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 4682e1947a5SStefano Zampini PetscValidType(B,1); 4692e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 4702e1947a5SStefano Zampini PetscFunctionReturn(0); 4712e1947a5SStefano Zampini } 4722e1947a5SStefano Zampini 4732e1947a5SStefano Zampini #undef __FUNCT__ 4742e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 4757230de76SStefano Zampini static PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 4762e1947a5SStefano Zampini { 4772e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 47828f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 4792e1947a5SStefano Zampini PetscErrorCode ierr; 4802e1947a5SStefano Zampini 4812e1947a5SStefano Zampini PetscFunctionBegin; 4826c4ed002SBarry Smith if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 48328f4e0baSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 48428f4e0baSStefano Zampini ierr = MatISComputeSF_Private(B);CHKERRQ(ierr); 48528f4e0baSStefano Zampini } 4862e1947a5SStefano Zampini if (!d_nnz) { 48728f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz; 4882e1947a5SStefano Zampini } else { 48928f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 4902e1947a5SStefano Zampini } 4912e1947a5SStefano Zampini if (!o_nnz) { 49228f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz; 4932e1947a5SStefano Zampini } else { 49428f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 4952e1947a5SStefano Zampini } 49628f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 49728f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 49828f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 49928f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 50028f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves;i++) { 50128f4e0baSStefano Zampini matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 5022e1947a5SStefano Zampini } 50328f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 50428f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 50528f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 5062e1947a5SStefano Zampini } 50728f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 50828f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 50928f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 5102e1947a5SStefano Zampini } 51128f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 5122e1947a5SStefano Zampini PetscFunctionReturn(0); 5132e1947a5SStefano Zampini } 514b4319ba4SBarry Smith 515b4319ba4SBarry Smith #undef __FUNCT__ 5163927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private" 5173927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 5183927de2eSStefano Zampini { 5193927de2eSStefano Zampini Mat_IS *matis = (Mat_IS*)(A->data); 5203927de2eSStefano Zampini PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 521ecf5a873SStefano Zampini const PetscInt *global_indices_r,*global_indices_c; 5223927de2eSStefano Zampini PetscInt i,j,bs,rows,cols; 5233927de2eSStefano Zampini PetscInt lrows,lcols; 5243927de2eSStefano Zampini PetscInt local_rows,local_cols; 5253927de2eSStefano Zampini PetscMPIInt nsubdomains; 5263927de2eSStefano Zampini PetscBool isdense,issbaij; 5273927de2eSStefano Zampini PetscErrorCode ierr; 5283927de2eSStefano Zampini 5293927de2eSStefano Zampini PetscFunctionBegin; 5303927de2eSStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr); 5313927de2eSStefano Zampini ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr); 5323927de2eSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 5333927de2eSStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 5343927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 5353927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 536ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 537ecf5a873SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 5387230de76SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 539ecf5a873SStefano Zampini } else { 540ecf5a873SStefano Zampini global_indices_c = global_indices_r; 541ecf5a873SStefano Zampini } 542ecf5a873SStefano Zampini 5433927de2eSStefano Zampini if (issbaij) { 5443927de2eSStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 5453927de2eSStefano Zampini } 5463927de2eSStefano Zampini /* 547ecf5a873SStefano Zampini An SF reduce is needed to sum up properly on shared rows. 5483927de2eSStefano Zampini Note that generally preallocation is not exact, since it overestimates nonzeros 5493927de2eSStefano Zampini */ 5503927de2eSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 5513927de2eSStefano Zampini ierr = MatISComputeSF_Private(A);CHKERRQ(ierr); 5523927de2eSStefano Zampini } 5533927de2eSStefano Zampini ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr); 5543927de2eSStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr); 5553927de2eSStefano Zampini /* All processes need to compute entire row ownership */ 5563927de2eSStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 5573927de2eSStefano Zampini ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 5583927de2eSStefano Zampini for (i=0;i<nsubdomains;i++) { 5593927de2eSStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 5603927de2eSStefano Zampini row_ownership[j] = i; 5613927de2eSStefano Zampini } 5623927de2eSStefano Zampini } 5637230de76SStefano Zampini ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 5643927de2eSStefano Zampini 5653927de2eSStefano Zampini /* 5663927de2eSStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 5673927de2eSStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 5683927de2eSStefano Zampini */ 5693927de2eSStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 5703927de2eSStefano Zampini /* preallocation as a MATAIJ */ 5713927de2eSStefano Zampini if (isdense) { /* special case for dense local matrices */ 5723927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 573ecf5a873SStefano Zampini PetscInt index_row = global_indices_r[i]; 5743927de2eSStefano Zampini for (j=i;j<local_rows;j++) { 5753927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 576ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[j]; 5773927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 5783927de2eSStefano Zampini my_dnz[i] += 1; 5793927de2eSStefano Zampini } else { /* offdiag block */ 5803927de2eSStefano Zampini my_onz[i] += 1; 5813927de2eSStefano Zampini } 5823927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 5833927de2eSStefano Zampini if (i != j) { 5843927de2eSStefano Zampini owner = row_ownership[index_col]; 5853927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 5863927de2eSStefano Zampini my_dnz[j] += 1; 5873927de2eSStefano Zampini } else { 5883927de2eSStefano Zampini my_onz[j] += 1; 5893927de2eSStefano Zampini } 5903927de2eSStefano Zampini } 5913927de2eSStefano Zampini } 5923927de2eSStefano Zampini } 593bb1015c3SStefano Zampini } else if (matis->A->ops->getrowij) { 594bb1015c3SStefano Zampini const PetscInt *ii,*jj,*jptr; 595bb1015c3SStefano Zampini PetscBool done; 596bb1015c3SStefano Zampini ierr = MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr); 597bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__); 598bb1015c3SStefano Zampini jptr = jj; 599bb1015c3SStefano Zampini for (i=0;i<local_rows;i++) { 600bb1015c3SStefano Zampini PetscInt index_row = global_indices_r[i]; 601bb1015c3SStefano Zampini for (j=0;j<ii[i+1]-ii[i];j++,jptr++) { 602bb1015c3SStefano Zampini PetscInt owner = row_ownership[index_row]; 603bb1015c3SStefano Zampini PetscInt index_col = global_indices_c[*jptr]; 604bb1015c3SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 605bb1015c3SStefano Zampini my_dnz[i] += 1; 606bb1015c3SStefano Zampini } else { /* offdiag block */ 607bb1015c3SStefano Zampini my_onz[i] += 1; 608bb1015c3SStefano Zampini } 609bb1015c3SStefano Zampini /* same as before, interchanging rows and cols */ 610bb1015c3SStefano Zampini if (issbaij && index_col != index_row) { 611bb1015c3SStefano Zampini owner = row_ownership[index_col]; 612bb1015c3SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 613bb1015c3SStefano Zampini my_dnz[*jptr] += 1; 614bb1015c3SStefano Zampini } else { 615bb1015c3SStefano Zampini my_onz[*jptr] += 1; 616bb1015c3SStefano Zampini } 617bb1015c3SStefano Zampini } 618bb1015c3SStefano Zampini } 619bb1015c3SStefano Zampini } 620bb1015c3SStefano Zampini ierr = MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr); 621bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__); 622bb1015c3SStefano Zampini } else { /* loop over rows and use MatGetRow */ 6233927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 6243927de2eSStefano Zampini const PetscInt *cols; 625ecf5a873SStefano Zampini PetscInt ncols,index_row = global_indices_r[i]; 6263927de2eSStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 6273927de2eSStefano Zampini for (j=0;j<ncols;j++) { 6283927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 629ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[cols[j]]; 6303927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 6313927de2eSStefano Zampini my_dnz[i] += 1; 6323927de2eSStefano Zampini } else { /* offdiag block */ 6333927de2eSStefano Zampini my_onz[i] += 1; 6343927de2eSStefano Zampini } 6353927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 636d9a9e74cSStefano Zampini if (issbaij && index_col != index_row) { 6373927de2eSStefano Zampini owner = row_ownership[index_col]; 6383927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 639d9a9e74cSStefano Zampini my_dnz[cols[j]] += 1; 6403927de2eSStefano Zampini } else { 641d9a9e74cSStefano Zampini my_onz[cols[j]] += 1; 6423927de2eSStefano Zampini } 6433927de2eSStefano Zampini } 6443927de2eSStefano Zampini } 6453927de2eSStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 6463927de2eSStefano Zampini } 6473927de2eSStefano Zampini } 648ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 649ecf5a873SStefano Zampini if (global_indices_c != global_indices_r) { 6507230de76SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 651ecf5a873SStefano Zampini } 6523927de2eSStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 653ecf5a873SStefano Zampini 654ecf5a873SStefano Zampini /* Reduce my_dnz and my_onz */ 6553927de2eSStefano Zampini if (maxreduce) { 6563927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 6573927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 658bb1015c3SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 6593927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 6603927de2eSStefano Zampini } else { 6613927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 6623927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 663bb1015c3SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 6643927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 6653927de2eSStefano Zampini } 6663927de2eSStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 6673927de2eSStefano Zampini 6683927de2eSStefano Zampini /* Resize preallocation if overestimated */ 6693927de2eSStefano Zampini for (i=0;i<lrows;i++) { 6703927de2eSStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 6713927de2eSStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 6723927de2eSStefano Zampini } 6733927de2eSStefano Zampini /* set preallocation */ 6743927de2eSStefano Zampini ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 6753927de2eSStefano Zampini for (i=0;i<lrows/bs;i++) { 6763927de2eSStefano Zampini dnz[i] = dnz[i*bs]/bs; 6773927de2eSStefano Zampini onz[i] = onz[i*bs]/bs; 6783927de2eSStefano Zampini } 6793927de2eSStefano Zampini ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 6803927de2eSStefano Zampini ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 6813927de2eSStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 6823927de2eSStefano Zampini if (issbaij) { 6833927de2eSStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 6843927de2eSStefano Zampini } 6853927de2eSStefano Zampini PetscFunctionReturn(0); 6863927de2eSStefano Zampini } 6873927de2eSStefano Zampini 6883927de2eSStefano Zampini #undef __FUNCT__ 689b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 6907230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 691b7ce53b6SStefano Zampini { 692b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 693d9a9e74cSStefano Zampini Mat local_mat; 694b7ce53b6SStefano Zampini /* info on mat */ 6953cfa4ea4SStefano Zampini PetscInt bs,rows,cols,lrows,lcols; 696b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 697686e3a49SStefano Zampini PetscBool isdense,issbaij,isseqaij; 6987c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 699b7ce53b6SStefano Zampini /* values insertion */ 700b7ce53b6SStefano Zampini PetscScalar *array; 701b7ce53b6SStefano Zampini /* work */ 702b7ce53b6SStefano Zampini PetscErrorCode ierr; 703b7ce53b6SStefano Zampini 704b7ce53b6SStefano Zampini PetscFunctionBegin; 705b7ce53b6SStefano Zampini /* get info from mat */ 7067c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 7077c03b4e8SStefano Zampini if (nsubdomains == 1) { 7087c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 7097c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 7107c03b4e8SStefano Zampini } else { 7117c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 7127c03b4e8SStefano Zampini } 7137c03b4e8SStefano Zampini PetscFunctionReturn(0); 7147c03b4e8SStefano Zampini } 715b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 716b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 7173cfa4ea4SStefano Zampini ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 718b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 719b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 720686e3a49SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 721b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 722b7ce53b6SStefano Zampini 723b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 724b7ce53b6SStefano Zampini MatType new_mat_type; 7253927de2eSStefano Zampini PetscBool issbaij_red; 726b7ce53b6SStefano Zampini 727b7ce53b6SStefano Zampini /* determining new matrix type */ 728b2566f29SBarry Smith ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 729b7ce53b6SStefano Zampini if (issbaij_red) { 730b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 731b7ce53b6SStefano Zampini } else { 732b7ce53b6SStefano Zampini if (bs>1) { 733b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 734b7ce53b6SStefano Zampini } else { 735b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 736b7ce53b6SStefano Zampini } 737b7ce53b6SStefano Zampini } 738b7ce53b6SStefano Zampini 7393927de2eSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 7403cfa4ea4SStefano Zampini ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 7413927de2eSStefano Zampini ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 7423927de2eSStefano Zampini ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr); 7433927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 744b7ce53b6SStefano Zampini } else { 7453cfa4ea4SStefano Zampini PetscInt mbs,mrows,mcols,mlrows,mlcols; 746b7ce53b6SStefano Zampini /* some checks */ 747b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 748b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 7493cfa4ea4SStefano Zampini ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 7506c4ed002SBarry Smith if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 7516c4ed002SBarry Smith if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 7526c4ed002SBarry Smith if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 7536c4ed002SBarry Smith if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 7546c4ed002SBarry Smith if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 755b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 756b7ce53b6SStefano Zampini } 757d9a9e74cSStefano Zampini 758d9a9e74cSStefano Zampini if (issbaij) { 759d9a9e74cSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 760d9a9e74cSStefano Zampini } else { 761d9a9e74cSStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 762d9a9e74cSStefano Zampini local_mat = matis->A; 763d9a9e74cSStefano Zampini } 764686e3a49SStefano Zampini 765b7ce53b6SStefano Zampini /* Set values */ 766ecf5a873SStefano Zampini ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 767b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 768ecf5a873SStefano Zampini PetscInt i,*dummy_rows,*dummy_cols; 769ecf5a873SStefano Zampini 770ecf5a873SStefano Zampini if (local_rows != local_cols) { 771ecf5a873SStefano Zampini ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 772ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 773ecf5a873SStefano Zampini for (i=0;i<local_cols;i++) dummy_cols[i] = i; 774ecf5a873SStefano Zampini } else { 775ecf5a873SStefano Zampini ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 776ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 777ecf5a873SStefano Zampini dummy_cols = dummy_rows; 778ecf5a873SStefano Zampini } 779b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 780d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 781ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 782d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 783ecf5a873SStefano Zampini if (dummy_rows != dummy_cols) { 784ecf5a873SStefano Zampini ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 785ecf5a873SStefano Zampini } else { 786ecf5a873SStefano Zampini ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 787ecf5a873SStefano Zampini } 788686e3a49SStefano Zampini } else if (isseqaij) { 789ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 790686e3a49SStefano Zampini PetscBool done; 791686e3a49SStefano Zampini 792d9a9e74cSStefano Zampini ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 793bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__); 794d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 795686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 796ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 797686e3a49SStefano Zampini } 798d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 799bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__); 800d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 801686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 802ecf5a873SStefano Zampini PetscInt i; 803c0962df8SStefano Zampini 804686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 805686e3a49SStefano Zampini PetscInt j; 806ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 807686e3a49SStefano Zampini 808ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 809ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 810ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 811686e3a49SStefano Zampini } 812b7ce53b6SStefano Zampini } 813b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 814d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 815b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 816b7ce53b6SStefano Zampini if (isdense) { 817b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 818b7ce53b6SStefano Zampini } 819b7ce53b6SStefano Zampini PetscFunctionReturn(0); 820b7ce53b6SStefano Zampini } 821b7ce53b6SStefano Zampini 822b7ce53b6SStefano Zampini #undef __FUNCT__ 823b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 824b7ce53b6SStefano Zampini /*@ 825b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 826b7ce53b6SStefano Zampini 827b7ce53b6SStefano Zampini Input Parameter: 828b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 829b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 830b7ce53b6SStefano Zampini 831b7ce53b6SStefano Zampini Output Parameter: 832b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 833b7ce53b6SStefano Zampini 834b7ce53b6SStefano Zampini Level: developer 835b7ce53b6SStefano Zampini 836eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 837b7ce53b6SStefano Zampini 838b7ce53b6SStefano Zampini .seealso: MATIS 839b7ce53b6SStefano Zampini @*/ 840b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 841b7ce53b6SStefano Zampini { 842b7ce53b6SStefano Zampini PetscErrorCode ierr; 843b7ce53b6SStefano Zampini 844b7ce53b6SStefano Zampini PetscFunctionBegin; 845b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 846b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 847b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 848b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 849b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 850b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 8516c4ed002SBarry Smith if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 852b7ce53b6SStefano Zampini } 853b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 854b7ce53b6SStefano Zampini PetscFunctionReturn(0); 855b7ce53b6SStefano Zampini } 856b7ce53b6SStefano Zampini 857b7ce53b6SStefano Zampini #undef __FUNCT__ 858ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 859ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 860ad6194a2SStefano Zampini { 861ad6194a2SStefano Zampini PetscErrorCode ierr; 862ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 863ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 864ad6194a2SStefano Zampini Mat B,localmat; 865ad6194a2SStefano Zampini 866ad6194a2SStefano Zampini PetscFunctionBegin; 867ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 868ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 869ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 870e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 871ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 872ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 873b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 874ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 875ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 876ad6194a2SStefano Zampini *newmat = B; 877ad6194a2SStefano Zampini PetscFunctionReturn(0); 878ad6194a2SStefano Zampini } 879ad6194a2SStefano Zampini 880ad6194a2SStefano Zampini #undef __FUNCT__ 88169796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 882a8116848SStefano Zampini static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 88369796d55SStefano Zampini { 88469796d55SStefano Zampini PetscErrorCode ierr; 88569796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 88669796d55SStefano Zampini PetscBool local_sym; 88769796d55SStefano Zampini 88869796d55SStefano Zampini PetscFunctionBegin; 88969796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 890b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 89169796d55SStefano Zampini PetscFunctionReturn(0); 89269796d55SStefano Zampini } 89369796d55SStefano Zampini 89469796d55SStefano Zampini #undef __FUNCT__ 89569796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 896a8116848SStefano Zampini static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 89769796d55SStefano Zampini { 89869796d55SStefano Zampini PetscErrorCode ierr; 89969796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 90069796d55SStefano Zampini PetscBool local_sym; 90169796d55SStefano Zampini 90269796d55SStefano Zampini PetscFunctionBegin; 90369796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 904b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 90569796d55SStefano Zampini PetscFunctionReturn(0); 90669796d55SStefano Zampini } 90769796d55SStefano Zampini 90869796d55SStefano Zampini #undef __FUNCT__ 909b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 910a8116848SStefano Zampini static PetscErrorCode MatDestroy_IS(Mat A) 911b4319ba4SBarry Smith { 912dfbe8321SBarry Smith PetscErrorCode ierr; 913b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 914b4319ba4SBarry Smith 915b4319ba4SBarry Smith PetscFunctionBegin; 9166bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 917e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 918e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 9196bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 9206bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 9213fd1c9e7SStefano Zampini ierr = VecDestroy(&b->counter);CHKERRQ(ierr); 922a8116848SStefano Zampini ierr = ISDestroy(&b->getsub_ris);CHKERRQ(ierr); 923a8116848SStefano Zampini ierr = ISDestroy(&b->getsub_cis);CHKERRQ(ierr); 924a8116848SStefano Zampini if (b->sf != b->csf) { 925a8116848SStefano Zampini ierr = PetscSFDestroy(&b->csf);CHKERRQ(ierr); 926a8116848SStefano Zampini ierr = PetscFree2(b->csf_rootdata,b->csf_leafdata);CHKERRQ(ierr); 927a8116848SStefano Zampini } 92828f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 92928f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 930bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 931dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 932bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 933b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 934b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 9352e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 936b4319ba4SBarry Smith PetscFunctionReturn(0); 937b4319ba4SBarry Smith } 938b4319ba4SBarry Smith 939b4319ba4SBarry Smith #undef __FUNCT__ 940b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 941a8116848SStefano Zampini static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 942b4319ba4SBarry Smith { 943dfbe8321SBarry Smith PetscErrorCode ierr; 944b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 945b4319ba4SBarry Smith PetscScalar zero = 0.0; 946b4319ba4SBarry Smith 947b4319ba4SBarry Smith PetscFunctionBegin; 948b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 949e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 950e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 951b4319ba4SBarry Smith 952b4319ba4SBarry Smith /* multiply the local matrix */ 953b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 954b4319ba4SBarry Smith 955b4319ba4SBarry Smith /* scatter product back into global memory */ 9562dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 957e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 958e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 959b4319ba4SBarry Smith PetscFunctionReturn(0); 960b4319ba4SBarry Smith } 961b4319ba4SBarry Smith 962b4319ba4SBarry Smith #undef __FUNCT__ 9632e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 964a8116848SStefano Zampini static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 9652e74eeadSLisandro Dalcin { 966650997f4SStefano Zampini Vec temp_vec; 9672e74eeadSLisandro Dalcin PetscErrorCode ierr; 9682e74eeadSLisandro Dalcin 9692e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 970650997f4SStefano Zampini if (v3 != v2) { 971650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 972650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 973650997f4SStefano Zampini } else { 974650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 975650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 976650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 977650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 978650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 979650997f4SStefano Zampini } 9802e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9812e74eeadSLisandro Dalcin } 9822e74eeadSLisandro Dalcin 9832e74eeadSLisandro Dalcin #undef __FUNCT__ 9842e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 985a8116848SStefano Zampini static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 9862e74eeadSLisandro Dalcin { 9872e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9882e74eeadSLisandro Dalcin PetscErrorCode ierr; 9892e74eeadSLisandro Dalcin 990e176bc59SStefano Zampini PetscFunctionBegin; 9912e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 992e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 993e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9942e74eeadSLisandro Dalcin 9952e74eeadSLisandro Dalcin /* multiply the local matrix */ 996e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 9972e74eeadSLisandro Dalcin 9982e74eeadSLisandro Dalcin /* scatter product back into global vector */ 999e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 1000e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1001e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10022e74eeadSLisandro Dalcin PetscFunctionReturn(0); 10032e74eeadSLisandro Dalcin } 10042e74eeadSLisandro Dalcin 10052e74eeadSLisandro Dalcin #undef __FUNCT__ 10062e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 1007a8116848SStefano Zampini static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 10082e74eeadSLisandro Dalcin { 1009650997f4SStefano Zampini Vec temp_vec; 10102e74eeadSLisandro Dalcin PetscErrorCode ierr; 10112e74eeadSLisandro Dalcin 10122e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 1013650997f4SStefano Zampini if (v3 != v2) { 1014650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 1015650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 1016650997f4SStefano Zampini } else { 1017650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 1018650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 1019650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 1020650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 1021650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 1022650997f4SStefano Zampini } 10232e74eeadSLisandro Dalcin PetscFunctionReturn(0); 10242e74eeadSLisandro Dalcin } 10252e74eeadSLisandro Dalcin 10262e74eeadSLisandro Dalcin #undef __FUNCT__ 1027b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 1028a8116848SStefano Zampini static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 1029b4319ba4SBarry Smith { 1030b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 1031dfbe8321SBarry Smith PetscErrorCode ierr; 1032b4319ba4SBarry Smith PetscViewer sviewer; 1033b4319ba4SBarry Smith 1034b4319ba4SBarry Smith PetscFunctionBegin; 10353f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 1036b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 10373f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 10386e520ac8SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1039b4319ba4SBarry Smith PetscFunctionReturn(0); 1040b4319ba4SBarry Smith } 1041b4319ba4SBarry Smith 1042b4319ba4SBarry Smith #undef __FUNCT__ 1043b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 1044a8116848SStefano Zampini static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 1045b4319ba4SBarry Smith { 1046dfbe8321SBarry Smith PetscErrorCode ierr; 1047e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 1048b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1049b4319ba4SBarry Smith IS from,to; 1050e176bc59SStefano Zampini Vec cglobal,rglobal; 1051b4319ba4SBarry Smith 1052b4319ba4SBarry Smith PetscFunctionBegin; 1053784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 1054e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 10553bbff08aSStefano Zampini /* Destroy any previous data */ 105670cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 105770cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 10583fd1c9e7SStefano Zampini ierr = VecDestroy(&is->counter);CHKERRQ(ierr); 1059e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 1060e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 10611c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 106228f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 106328f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 10643bbff08aSStefano Zampini 10653bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 1066fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1067fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1068fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 1069fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 1070b4319ba4SBarry Smith 1071b4319ba4SBarry Smith /* Create the local matrix A */ 1072e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 1073e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 1074e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 1075e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 1076f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 1077e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 1078e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 1079e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 1080ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 1081ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 1082b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 1083b4319ba4SBarry Smith 1084f26d0771SStefano Zampini if (!is->islocalref) { /* setup scatters and local vectors for MatMult */ 1085b4319ba4SBarry Smith /* Create the local work vectors */ 10863bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 1087b4319ba4SBarry Smith 1088e176bc59SStefano Zampini /* setup the global to local scatters */ 1089e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 1090e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 1091784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 1092e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 1093e176bc59SStefano Zampini if (rmapping != cmapping) { 1094e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 1095e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 1096e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 1097e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 1098e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 1099e176bc59SStefano Zampini } else { 1100e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 1101e176bc59SStefano Zampini is->cctx = is->rctx; 1102e176bc59SStefano Zampini } 11033fd1c9e7SStefano Zampini 11043fd1c9e7SStefano Zampini /* interface counter vector (local) */ 11053fd1c9e7SStefano Zampini ierr = VecDuplicate(is->y,&is->counter);CHKERRQ(ierr); 11063fd1c9e7SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 11073fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11083fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 11093fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 11103fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 11113fd1c9e7SStefano Zampini 11123fd1c9e7SStefano Zampini /* free workspace */ 1113e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 1114e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 11156bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 11166bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1117f26d0771SStefano Zampini } 111848ff6bf3SStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 1119b4319ba4SBarry Smith PetscFunctionReturn(0); 1120b4319ba4SBarry Smith } 1121b4319ba4SBarry Smith 11222e74eeadSLisandro Dalcin #undef __FUNCT__ 11232e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 1124a8116848SStefano Zampini static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 11252e74eeadSLisandro Dalcin { 11262e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 11272e74eeadSLisandro Dalcin PetscErrorCode ierr; 112897563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 112997563a80SStefano Zampini PetscInt i,zm,zn; 113097563a80SStefano Zampini #endif 1131f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 11322e74eeadSLisandro Dalcin 11332e74eeadSLisandro Dalcin PetscFunctionBegin; 11342e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 1135f26d0771SStefano Zampini if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n); 113697563a80SStefano Zampini /* count negative indices */ 113797563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 113897563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 11392e74eeadSLisandro Dalcin #endif 114097563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 114197563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 114297563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 114397563a80SStefano Zampini /* count negative indices (should be the same as before) */ 114497563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 114597563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 114697563a80SStefano 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"); 114797563a80SStefano 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"); 114897563a80SStefano Zampini #endif 11492e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 11502e74eeadSLisandro Dalcin PetscFunctionReturn(0); 11512e74eeadSLisandro Dalcin } 11522e74eeadSLisandro Dalcin 1153b4319ba4SBarry Smith #undef __FUNCT__ 115497563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS" 1155a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 115697563a80SStefano Zampini { 115797563a80SStefano Zampini Mat_IS *is = (Mat_IS*)mat->data; 115897563a80SStefano Zampini PetscErrorCode ierr; 115997563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 116097563a80SStefano Zampini PetscInt i,zm,zn; 116197563a80SStefano Zampini #endif 1162f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 116397563a80SStefano Zampini 116497563a80SStefano Zampini PetscFunctionBegin; 116597563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 1166f26d0771SStefano Zampini if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column block indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n); 116797563a80SStefano Zampini /* count negative indices */ 116897563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 116997563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 117097563a80SStefano Zampini #endif 117197563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 117297563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 117397563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 117497563a80SStefano Zampini /* count negative indices (should be the same as before) */ 117597563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 117697563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 117797563a80SStefano 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"); 117897563a80SStefano 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"); 117997563a80SStefano Zampini #endif 118097563a80SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 118197563a80SStefano Zampini PetscFunctionReturn(0); 118297563a80SStefano Zampini } 118397563a80SStefano Zampini 118497563a80SStefano Zampini #undef __FUNCT__ 1185b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 1186a8116848SStefano Zampini static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1187b4319ba4SBarry Smith { 1188dfbe8321SBarry Smith PetscErrorCode ierr; 1189b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1190b4319ba4SBarry Smith 1191b4319ba4SBarry Smith PetscFunctionBegin; 1192b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1193b4319ba4SBarry Smith PetscFunctionReturn(0); 1194b4319ba4SBarry Smith } 1195b4319ba4SBarry Smith 1196b4319ba4SBarry Smith #undef __FUNCT__ 1197f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 1198a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1199f0006bf2SLisandro Dalcin { 1200f0006bf2SLisandro Dalcin PetscErrorCode ierr; 1201f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 1202f0006bf2SLisandro Dalcin 1203f0006bf2SLisandro Dalcin PetscFunctionBegin; 1204f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1205f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 1206f0006bf2SLisandro Dalcin } 1207f0006bf2SLisandro Dalcin 1208f0006bf2SLisandro Dalcin #undef __FUNCT__ 1209f0ae7da4SStefano Zampini #define __FUNCT__ "MatISZeroRowsColumnsLocal_Private" 1210f0ae7da4SStefano Zampini static PetscErrorCode MatISZeroRowsColumnsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,PetscBool columns) 1211f0ae7da4SStefano Zampini { 1212f0ae7da4SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 1213f0ae7da4SStefano Zampini PetscErrorCode ierr; 1214f0ae7da4SStefano Zampini 1215f0ae7da4SStefano Zampini PetscFunctionBegin; 1216f0ae7da4SStefano Zampini if (!n) { 1217f0ae7da4SStefano Zampini is->pure_neumann = PETSC_TRUE; 1218f0ae7da4SStefano Zampini } else { 1219f0ae7da4SStefano Zampini PetscInt i; 1220f0ae7da4SStefano Zampini is->pure_neumann = PETSC_FALSE; 1221f0ae7da4SStefano Zampini 1222f0ae7da4SStefano Zampini if (columns) { 1223f0ae7da4SStefano Zampini ierr = MatZeroRowsColumns(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 1224f0ae7da4SStefano Zampini } else { 1225f0ae7da4SStefano Zampini ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 1226f0ae7da4SStefano Zampini } 1227f0ae7da4SStefano Zampini if (diag != 0.) { 1228f0ae7da4SStefano Zampini const PetscScalar *array; 1229f0ae7da4SStefano Zampini ierr = VecGetArrayRead(is->counter,&array);CHKERRQ(ierr); 1230f0ae7da4SStefano Zampini for (i=0; i<n; i++) { 1231f0ae7da4SStefano Zampini ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 1232f0ae7da4SStefano Zampini } 1233f0ae7da4SStefano Zampini ierr = VecRestoreArrayRead(is->counter,&array);CHKERRQ(ierr); 1234f0ae7da4SStefano Zampini } 1235f0ae7da4SStefano Zampini ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1236f0ae7da4SStefano Zampini ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1237f0ae7da4SStefano Zampini } 1238f0ae7da4SStefano Zampini PetscFunctionReturn(0); 1239f0ae7da4SStefano Zampini } 1240f0ae7da4SStefano Zampini 1241f0ae7da4SStefano Zampini #undef __FUNCT__ 1242f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_Private_IS" 1243f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Private_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b,PetscBool columns) 12442e74eeadSLisandro Dalcin { 12456e520ac8SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 12466e520ac8SStefano Zampini PetscInt nr,nl,len,i; 12476e520ac8SStefano Zampini PetscInt *lrows; 12482e74eeadSLisandro Dalcin PetscErrorCode ierr; 12492e74eeadSLisandro Dalcin 12502e74eeadSLisandro Dalcin PetscFunctionBegin; 1251f0ae7da4SStefano Zampini #if defined(PETSC_USE_DEBUG) 1252f0ae7da4SStefano Zampini if (columns || diag != 0. || (x && b)) { 1253f0ae7da4SStefano Zampini PetscBool cong; 1254f0ae7da4SStefano Zampini ierr = PetscLayoutCompare(A->rmap,A->cmap,&cong);CHKERRQ(ierr); 1255f0ae7da4SStefano Zampini if (!cong && columns) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Columns can be zeroed if and only if A->rmap and A->cmap are congruent for MATIS"); 1256f0ae7da4SStefano Zampini if (!cong && diag != 0.) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Nonzero diagonal value supported if and only if A->rmap and A->cmap are congruent for MATIS"); 1257f0ae7da4SStefano Zampini if (!cong && x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"A->rmap and A->cmap need to be congruent"); 1258f0ae7da4SStefano Zampini } 1259f0ae7da4SStefano Zampini #endif 12606e520ac8SStefano Zampini /* get locally owned rows */ 1261f0ae7da4SStefano Zampini ierr = PetscLayoutMapLocal_Private(A->rmap,n,rows,&len,&lrows,NULL);CHKERRQ(ierr); 12626e520ac8SStefano Zampini /* fix right hand side if needed */ 12636e520ac8SStefano Zampini if (x && b) { 12646e520ac8SStefano Zampini const PetscScalar *xx; 12656e520ac8SStefano Zampini PetscScalar *bb; 12666e520ac8SStefano Zampini 12676e520ac8SStefano Zampini ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr); 12686e520ac8SStefano Zampini ierr = VecGetArray(b, &bb);CHKERRQ(ierr); 12696e520ac8SStefano Zampini for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]]; 12706e520ac8SStefano Zampini ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr); 12716e520ac8SStefano Zampini ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr); 12722e74eeadSLisandro Dalcin } 12736e520ac8SStefano Zampini /* get rows associated to the local matrices */ 12746e520ac8SStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 12756e520ac8SStefano Zampini ierr = MatISComputeSF_Private(A);CHKERRQ(ierr); 12766e520ac8SStefano Zampini } 12776e520ac8SStefano Zampini ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr); 12786e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr); 12796e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 12806e520ac8SStefano Zampini for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1; 12816e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 12826e520ac8SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 12836e520ac8SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 12846e520ac8SStefano Zampini ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr); 12856e520ac8SStefano Zampini for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i; 1286f0ae7da4SStefano Zampini ierr = MatISZeroRowsColumnsLocal_Private(A,nr,lrows,diag,columns);CHKERRQ(ierr); 12876e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 12882e74eeadSLisandro Dalcin PetscFunctionReturn(0); 12892e74eeadSLisandro Dalcin } 12902e74eeadSLisandro Dalcin 12912e74eeadSLisandro Dalcin #undef __FUNCT__ 1292f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRows_IS" 1293f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 1294b4319ba4SBarry Smith { 1295dfbe8321SBarry Smith PetscErrorCode ierr; 1296b4319ba4SBarry Smith 1297b4319ba4SBarry Smith PetscFunctionBegin; 1298f0ae7da4SStefano Zampini ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_FALSE);CHKERRQ(ierr); 1299f0ae7da4SStefano Zampini PetscFunctionReturn(0); 1300f0ae7da4SStefano Zampini } 13012205254eSKarl Rupp 1302f0ae7da4SStefano Zampini #undef __FUNCT__ 1303f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_IS" 1304f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 1305f0ae7da4SStefano Zampini { 1306f0ae7da4SStefano Zampini PetscErrorCode ierr; 1307f0ae7da4SStefano Zampini 1308f0ae7da4SStefano Zampini PetscFunctionBegin; 1309f0ae7da4SStefano Zampini ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_TRUE);CHKERRQ(ierr); 1310b4319ba4SBarry Smith PetscFunctionReturn(0); 1311b4319ba4SBarry Smith } 1312b4319ba4SBarry Smith 1313b4319ba4SBarry Smith #undef __FUNCT__ 1314b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 1315a8116848SStefano Zampini static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 1316b4319ba4SBarry Smith { 1317b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1318dfbe8321SBarry Smith PetscErrorCode ierr; 1319b4319ba4SBarry Smith 1320b4319ba4SBarry Smith PetscFunctionBegin; 1321b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 1322b4319ba4SBarry Smith PetscFunctionReturn(0); 1323b4319ba4SBarry Smith } 1324b4319ba4SBarry Smith 1325b4319ba4SBarry Smith #undef __FUNCT__ 1326b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 1327a8116848SStefano Zampini static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 1328b4319ba4SBarry Smith { 1329b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1330dfbe8321SBarry Smith PetscErrorCode ierr; 1331b4319ba4SBarry Smith 1332b4319ba4SBarry Smith PetscFunctionBegin; 1333b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 1334b4319ba4SBarry Smith PetscFunctionReturn(0); 1335b4319ba4SBarry Smith } 1336b4319ba4SBarry Smith 1337b4319ba4SBarry Smith #undef __FUNCT__ 1338b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 1339a8116848SStefano Zampini static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 1340b4319ba4SBarry Smith { 1341b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 1342b4319ba4SBarry Smith 1343b4319ba4SBarry Smith PetscFunctionBegin; 1344b4319ba4SBarry Smith *local = is->A; 1345b4319ba4SBarry Smith PetscFunctionReturn(0); 1346b4319ba4SBarry Smith } 1347b4319ba4SBarry Smith 1348b4319ba4SBarry Smith #undef __FUNCT__ 1349b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 1350b4319ba4SBarry Smith /*@ 1351b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 1352b4319ba4SBarry Smith 1353b4319ba4SBarry Smith Input Parameter: 1354b4319ba4SBarry Smith . mat - the matrix 1355b4319ba4SBarry Smith 1356b4319ba4SBarry Smith Output Parameter: 1357eb82efa4SStefano Zampini . local - the local matrix 1358b4319ba4SBarry Smith 1359b4319ba4SBarry Smith Level: advanced 1360b4319ba4SBarry Smith 1361b4319ba4SBarry Smith Notes: 1362b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 1363b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 1364b4319ba4SBarry Smith of the MatSetValues() operation. 1365b4319ba4SBarry Smith 1366b4319ba4SBarry Smith .seealso: MATIS 1367b4319ba4SBarry Smith @*/ 13687087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 1369b4319ba4SBarry Smith { 13704ac538c5SBarry Smith PetscErrorCode ierr; 1371b4319ba4SBarry Smith 1372b4319ba4SBarry Smith PetscFunctionBegin; 13730700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1374b4319ba4SBarry Smith PetscValidPointer(local,2); 13754ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 1376b4319ba4SBarry Smith PetscFunctionReturn(0); 1377b4319ba4SBarry Smith } 1378b4319ba4SBarry Smith 13793b03a366Sstefano_zampini #undef __FUNCT__ 13803b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 1381a8116848SStefano Zampini static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 13823b03a366Sstefano_zampini { 13833b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 13843b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 13853b03a366Sstefano_zampini PetscErrorCode ierr; 13863b03a366Sstefano_zampini 13873b03a366Sstefano_zampini PetscFunctionBegin; 13884e4c7dbeSStefano Zampini if (is->A) { 13893b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 13903b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 1391f0ae7da4SStefano Zampini 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)",orows,ocols,nrows,ncols); 13924e4c7dbeSStefano Zampini } 13933b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 13943b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 13953b03a366Sstefano_zampini is->A = local; 13963b03a366Sstefano_zampini PetscFunctionReturn(0); 13973b03a366Sstefano_zampini } 13983b03a366Sstefano_zampini 13993b03a366Sstefano_zampini #undef __FUNCT__ 14003b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 14013b03a366Sstefano_zampini /*@ 1402eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 14033b03a366Sstefano_zampini 14043b03a366Sstefano_zampini Input Parameter: 14053b03a366Sstefano_zampini . mat - the matrix 1406eb82efa4SStefano Zampini . local - the local matrix 14073b03a366Sstefano_zampini 14083b03a366Sstefano_zampini Output Parameter: 14093b03a366Sstefano_zampini 14103b03a366Sstefano_zampini Level: advanced 14113b03a366Sstefano_zampini 14123b03a366Sstefano_zampini Notes: 14133b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 14143b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 14153b03a366Sstefano_zampini 14163b03a366Sstefano_zampini .seealso: MATIS 14173b03a366Sstefano_zampini @*/ 14183b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 14193b03a366Sstefano_zampini { 14203b03a366Sstefano_zampini PetscErrorCode ierr; 14213b03a366Sstefano_zampini 14223b03a366Sstefano_zampini PetscFunctionBegin; 14233b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1424b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 14253b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 14263b03a366Sstefano_zampini PetscFunctionReturn(0); 14273b03a366Sstefano_zampini } 14283b03a366Sstefano_zampini 14296726f965SBarry Smith #undef __FUNCT__ 14306726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 1431a8116848SStefano Zampini static PetscErrorCode MatZeroEntries_IS(Mat A) 14326726f965SBarry Smith { 14336726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 14346726f965SBarry Smith PetscErrorCode ierr; 14356726f965SBarry Smith 14366726f965SBarry Smith PetscFunctionBegin; 14376726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 14386726f965SBarry Smith PetscFunctionReturn(0); 14396726f965SBarry Smith } 14406726f965SBarry Smith 14416726f965SBarry Smith #undef __FUNCT__ 14422e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 1443a8116848SStefano Zampini static PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 14442e74eeadSLisandro Dalcin { 14452e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 14462e74eeadSLisandro Dalcin PetscErrorCode ierr; 14472e74eeadSLisandro Dalcin 14482e74eeadSLisandro Dalcin PetscFunctionBegin; 14492e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 14502e74eeadSLisandro Dalcin PetscFunctionReturn(0); 14512e74eeadSLisandro Dalcin } 14522e74eeadSLisandro Dalcin 14532e74eeadSLisandro Dalcin #undef __FUNCT__ 14542e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 1455a8116848SStefano Zampini static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 14562e74eeadSLisandro Dalcin { 14572e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 14582e74eeadSLisandro Dalcin PetscErrorCode ierr; 14592e74eeadSLisandro Dalcin 14602e74eeadSLisandro Dalcin PetscFunctionBegin; 14612e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 1462e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 14632e74eeadSLisandro Dalcin 14642e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 14652e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 1466e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1467e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 14682e74eeadSLisandro Dalcin PetscFunctionReturn(0); 14692e74eeadSLisandro Dalcin } 14702e74eeadSLisandro Dalcin 14712e74eeadSLisandro Dalcin #undef __FUNCT__ 14726726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 1473a8116848SStefano Zampini static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 14746726f965SBarry Smith { 14756726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 14766726f965SBarry Smith PetscErrorCode ierr; 14776726f965SBarry Smith 14786726f965SBarry Smith PetscFunctionBegin; 14794e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 14806726f965SBarry Smith PetscFunctionReturn(0); 14816726f965SBarry Smith } 14826726f965SBarry Smith 1483284134d9SBarry Smith #undef __FUNCT__ 1484f26d0771SStefano Zampini #define __FUNCT__ "MatAXPY_IS" 1485f26d0771SStefano Zampini static PetscErrorCode MatAXPY_IS(Mat Y,PetscScalar a,Mat X,MatStructure str) 1486f26d0771SStefano Zampini { 1487f26d0771SStefano Zampini Mat_IS *y = (Mat_IS*)Y->data; 1488f26d0771SStefano Zampini Mat_IS *x; 1489f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1490f26d0771SStefano Zampini PetscBool ismatis; 1491f26d0771SStefano Zampini #endif 1492f26d0771SStefano Zampini PetscErrorCode ierr; 1493f26d0771SStefano Zampini 1494f26d0771SStefano Zampini PetscFunctionBegin; 1495f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1496f26d0771SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)X,MATIS,&ismatis);CHKERRQ(ierr); 1497f26d0771SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS"); 1498f26d0771SStefano Zampini #endif 1499f26d0771SStefano Zampini x = (Mat_IS*)X->data; 1500f26d0771SStefano Zampini ierr = MatAXPY(y->A,a,x->A,str);CHKERRQ(ierr); 1501f26d0771SStefano Zampini PetscFunctionReturn(0); 1502f26d0771SStefano Zampini } 1503f26d0771SStefano Zampini 1504f26d0771SStefano Zampini #undef __FUNCT__ 1505f26d0771SStefano Zampini #define __FUNCT__ "MatGetLocalSubMatrix_IS" 1506f26d0771SStefano Zampini static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat) 1507f26d0771SStefano Zampini { 1508f26d0771SStefano Zampini Mat lA; 1509f26d0771SStefano Zampini Mat_IS *matis; 1510f26d0771SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 1511f26d0771SStefano Zampini IS is; 1512f26d0771SStefano Zampini const PetscInt *rg,*rl; 1513f26d0771SStefano Zampini PetscInt nrg; 1514f26d0771SStefano Zampini PetscInt N,M,nrl,i,*idxs; 1515f26d0771SStefano Zampini PetscErrorCode ierr; 1516f26d0771SStefano Zampini 1517f26d0771SStefano Zampini PetscFunctionBegin; 1518f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 1519f26d0771SStefano Zampini ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr); 1520f26d0771SStefano Zampini ierr = ISGetIndices(row,&rl);CHKERRQ(ierr); 1521f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr); 1522f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1523f0ae7da4SStefano 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",i,rl[i],nrg); 1524f26d0771SStefano Zampini #endif 1525f26d0771SStefano Zampini ierr = PetscMalloc1(nrg,&idxs);CHKERRQ(ierr); 1526f26d0771SStefano Zampini /* map from [0,nrl) to row */ 1527f26d0771SStefano Zampini for (i=0;i<nrl;i++) idxs[i] = rl[i]; 1528f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1529f26d0771SStefano Zampini for (i=nrl;i<nrg;i++) idxs[i] = nrg; 1530f26d0771SStefano Zampini #else 1531f26d0771SStefano Zampini for (i=nrl;i<nrg;i++) idxs[i] = -1; 1532f26d0771SStefano Zampini #endif 1533f26d0771SStefano Zampini ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr); 1534f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 1535f26d0771SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 1536f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 1537f26d0771SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 1538f26d0771SStefano Zampini /* compute new l2g map for columns */ 1539f26d0771SStefano Zampini if (col != row || A->rmap->mapping != A->cmap->mapping) { 1540f26d0771SStefano Zampini const PetscInt *cg,*cl; 1541f26d0771SStefano Zampini PetscInt ncg; 1542f26d0771SStefano Zampini PetscInt ncl; 1543f26d0771SStefano Zampini 1544f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 1545f26d0771SStefano Zampini ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr); 1546f26d0771SStefano Zampini ierr = ISGetIndices(col,&cl);CHKERRQ(ierr); 1547f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr); 1548f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1549f0ae7da4SStefano 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",i,cl[i],ncg); 1550f26d0771SStefano Zampini #endif 1551f26d0771SStefano Zampini ierr = PetscMalloc1(ncg,&idxs);CHKERRQ(ierr); 1552f26d0771SStefano Zampini /* map from [0,ncl) to col */ 1553f26d0771SStefano Zampini for (i=0;i<ncl;i++) idxs[i] = cl[i]; 1554f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1555f26d0771SStefano Zampini for (i=ncl;i<ncg;i++) idxs[i] = ncg; 1556f26d0771SStefano Zampini #else 1557f26d0771SStefano Zampini for (i=ncl;i<ncg;i++) idxs[i] = -1; 1558f26d0771SStefano Zampini #endif 1559f26d0771SStefano Zampini ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr); 1560f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 1561f26d0771SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 1562f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 1563f26d0771SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 1564f26d0771SStefano Zampini } else { 1565f26d0771SStefano Zampini ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr); 1566f26d0771SStefano Zampini cl2g = rl2g; 1567f26d0771SStefano Zampini } 1568f26d0771SStefano Zampini /* create the MATIS submatrix */ 1569f26d0771SStefano Zampini ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 1570f26d0771SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr); 1571f26d0771SStefano Zampini ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1572f26d0771SStefano Zampini ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr); 1573b0aa3428SStefano Zampini matis = (Mat_IS*)((*submat)->data); 1574f26d0771SStefano Zampini matis->islocalref = PETSC_TRUE; 1575f26d0771SStefano Zampini ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr); 1576f26d0771SStefano Zampini ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr); 1577f26d0771SStefano Zampini ierr = MatISSetLocalMat(*submat,lA);CHKERRQ(ierr); 1578f26d0771SStefano Zampini ierr = MatSetUp(*submat);CHKERRQ(ierr); 1579f26d0771SStefano Zampini ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1580f26d0771SStefano Zampini ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1581f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 1582f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 1583f26d0771SStefano Zampini /* remove unsupported ops */ 1584f26d0771SStefano Zampini ierr = PetscMemzero((*submat)->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1585f26d0771SStefano Zampini (*submat)->ops->destroy = MatDestroy_IS; 1586f26d0771SStefano Zampini (*submat)->ops->setvalueslocal = MatSetValuesLocal_SubMat_IS; 1587f26d0771SStefano Zampini (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS; 1588f26d0771SStefano Zampini (*submat)->ops->assemblybegin = MatAssemblyBegin_IS; 1589f26d0771SStefano Zampini (*submat)->ops->assemblyend = MatAssemblyEnd_IS; 1590f26d0771SStefano Zampini PetscFunctionReturn(0); 1591f26d0771SStefano Zampini } 1592f26d0771SStefano Zampini 1593f26d0771SStefano Zampini #undef __FUNCT__ 1594284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 1595284134d9SBarry Smith /*@ 15963c212e90SHong Zhang MatCreateIS - Creates a "process" unassembled matrix, assembled on each 1597284134d9SBarry Smith process but not across processes. 1598284134d9SBarry Smith 1599284134d9SBarry Smith Input Parameters: 1600284134d9SBarry Smith + comm - MPI communicator that will share the matrix 1601e176bc59SStefano Zampini . bs - block size of the matrix 1602df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 1603e176bc59SStefano Zampini . rmap - local to global map for rows 1604e176bc59SStefano Zampini - cmap - local to global map for cols 1605284134d9SBarry Smith 1606284134d9SBarry Smith Output Parameter: 1607284134d9SBarry Smith . A - the resulting matrix 1608284134d9SBarry Smith 16098e6c10adSSatish Balay Level: advanced 16108e6c10adSSatish Balay 16113c212e90SHong Zhang Notes: See MATIS for more details. 16123c212e90SHong Zhang m and n are NOT related to the size of the map; they are the size of the part of the vector owned 1613e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 16143c212e90SHong Zhang If either rmap or cmap are NULL, then the matrix is assumed to be square. 1615284134d9SBarry Smith 1616284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 1617284134d9SBarry Smith @*/ 1618e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 1619284134d9SBarry Smith { 1620284134d9SBarry Smith PetscErrorCode ierr; 1621284134d9SBarry Smith 1622284134d9SBarry Smith PetscFunctionBegin; 1623e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 1624284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1625d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 1626284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1627284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 1628e176bc59SStefano Zampini if (rmap && cmap) { 1629e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 1630e176bc59SStefano Zampini } else if (!rmap) { 1631e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 1632e176bc59SStefano Zampini } else { 1633e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 1634e176bc59SStefano Zampini } 1635284134d9SBarry Smith PetscFunctionReturn(0); 1636284134d9SBarry Smith } 1637284134d9SBarry Smith 1638b4319ba4SBarry Smith /*MC 1639f26d0771SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition methods (e.g. PCBDDC or KSPFETIDP). 1640b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 1641b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 1642b4319ba4SBarry Smith product is handled "implicitly". 1643b4319ba4SBarry Smith 1644b4319ba4SBarry Smith Operations Provided: 16456726f965SBarry Smith + MatMult() 16462e74eeadSLisandro Dalcin . MatMultAdd() 16472e74eeadSLisandro Dalcin . MatMultTranspose() 16482e74eeadSLisandro Dalcin . MatMultTransposeAdd() 16496726f965SBarry Smith . MatZeroEntries() 16506726f965SBarry Smith . MatSetOption() 16512e74eeadSLisandro Dalcin . MatZeroRows() 16522e74eeadSLisandro Dalcin . MatSetValues() 165397563a80SStefano Zampini . MatSetValuesBlocked() 16546726f965SBarry Smith . MatSetValuesLocal() 165597563a80SStefano Zampini . MatSetValuesBlockedLocal() 16562e74eeadSLisandro Dalcin . MatScale() 16572e74eeadSLisandro Dalcin . MatGetDiagonal() 16582b404112SStefano Zampini . MatMissingDiagonal() 16592b404112SStefano Zampini . MatDuplicate() 16602b404112SStefano Zampini . MatCopy() 1661f26d0771SStefano Zampini . MatAXPY() 1662f26d0771SStefano Zampini . MatGetSubMatrix() 1663f26d0771SStefano Zampini . MatGetLocalSubMatrix() 1664*d7f69cd0SStefano Zampini . MatTranspose() 16656726f965SBarry Smith - MatSetLocalToGlobalMapping() 1666b4319ba4SBarry Smith 1667b4319ba4SBarry Smith Options Database Keys: 1668b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1669b4319ba4SBarry Smith 1670b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1671b4319ba4SBarry Smith 1672b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1673b4319ba4SBarry Smith 1674b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1675eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1676b4319ba4SBarry Smith 1677b4319ba4SBarry Smith Level: advanced 1678b4319ba4SBarry Smith 1679f26d0771SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), MatCreateIS(), PCBDDC, KSPFETIDP 1680b4319ba4SBarry Smith 1681b4319ba4SBarry Smith M*/ 1682b4319ba4SBarry Smith 1683b4319ba4SBarry Smith #undef __FUNCT__ 1684b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 16858cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1686b4319ba4SBarry Smith { 1687dfbe8321SBarry Smith PetscErrorCode ierr; 1688b4319ba4SBarry Smith Mat_IS *b; 1689b4319ba4SBarry Smith 1690b4319ba4SBarry Smith PetscFunctionBegin; 1691b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1692b4319ba4SBarry Smith A->data = (void*)b; 1693b4319ba4SBarry Smith 1694e176bc59SStefano Zampini /* matrix ops */ 1695e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1696b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 16972e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 16982e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 16992e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1700b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1701b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 17022e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 170398921651SStefano Zampini A->ops->setvaluesblocked = MatSetValuesBlocked_IS; 1704b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1705f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 17062e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1707f0ae7da4SStefano Zampini A->ops->zerorowscolumns = MatZeroRowsColumns_IS; 1708b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1709b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1710b4319ba4SBarry Smith A->ops->view = MatView_IS; 17116726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 17122e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 17132e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 17146726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 171569796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 171669796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1717ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 17186bd84002SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_IS; 17192b404112SStefano Zampini A->ops->copy = MatCopy_IS; 1720659959c5SStefano Zampini A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_IS; 1721a8116848SStefano Zampini A->ops->getsubmatrix = MatGetSubMatrix_IS; 1722f26d0771SStefano Zampini A->ops->axpy = MatAXPY_IS; 17233fd1c9e7SStefano Zampini A->ops->diagonalset = MatDiagonalSet_IS; 17243fd1c9e7SStefano Zampini A->ops->shift = MatShift_IS; 1725*d7f69cd0SStefano Zampini A->ops->transpose = MatTranspose_IS; 1726b4319ba4SBarry Smith 1727b7ce53b6SStefano Zampini /* special MATIS functions */ 1728bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1729bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1730b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 17312e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 173217667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1733b4319ba4SBarry Smith PetscFunctionReturn(0); 1734b4319ba4SBarry Smith } 1735