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*/ 114f2d7cafSStefano Zampini #include <petsc/private/sfimpl.h> 1228f4e0baSStefano Zampini 13f26d0771SStefano Zampini #define MATIS_MAX_ENTRIES_INSERTION 2048 14f26d0771SStefano Zampini 15cf0a3239SStefano Zampini #undef __FUNCT__ 16cf0a3239SStefano Zampini #define __FUNCT__ "MatISSetUpSF" 17cf0a3239SStefano Zampini /*@ 18*3d996552SStefano Zampini MatISSetUpSF - Setup star forest objects used by MatIS. 19cf0a3239SStefano Zampini 20cf0a3239SStefano Zampini Collective on MPI_Comm 21cf0a3239SStefano Zampini 22cf0a3239SStefano Zampini Input Parameters: 23cf0a3239SStefano Zampini + A - the matrix 24cf0a3239SStefano Zampini 25cf0a3239SStefano Zampini Level: advanced 26cf0a3239SStefano Zampini 27*3d996552SStefano Zampini Notes: This function does not need to be called by the user. 28cf0a3239SStefano Zampini 29cf0a3239SStefano Zampini .keywords: matrix 30cf0a3239SStefano Zampini 31cf0a3239SStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatISSetPreallocation(), MatISGetLocalMat() 32cf0a3239SStefano Zampini @*/ 33cf0a3239SStefano Zampini PetscErrorCode MatISSetUpSF(Mat A) 34cf0a3239SStefano Zampini { 35cf0a3239SStefano Zampini PetscErrorCode ierr; 36cf0a3239SStefano Zampini 37cf0a3239SStefano Zampini PetscFunctionBegin; 38cf0a3239SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 39cf0a3239SStefano Zampini PetscValidType(A,1); 40cf0a3239SStefano Zampini ierr = PetscTryMethod(A,"MatISSetUpSF_C",(Mat),(A));CHKERRQ(ierr); 41cf0a3239SStefano Zampini PetscFunctionReturn(0); 42cf0a3239SStefano Zampini } 43a72627d2SStefano Zampini 4428f4e0baSStefano Zampini #undef __FUNCT__ 453fd1c9e7SStefano Zampini #define __FUNCT__ "MatDiagonalSet_IS" 463fd1c9e7SStefano Zampini PetscErrorCode MatDiagonalSet_IS(Mat A,Vec D,InsertMode insmode) 473fd1c9e7SStefano Zampini { 483fd1c9e7SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 493fd1c9e7SStefano Zampini PetscErrorCode ierr; 503fd1c9e7SStefano Zampini 513fd1c9e7SStefano Zampini PetscFunctionBegin; 523fd1c9e7SStefano Zampini if (!D) { /* this code branch is used by MatShift_IS */ 533fd1c9e7SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 543fd1c9e7SStefano Zampini } else { 553fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 563fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 573fd1c9e7SStefano Zampini } 583fd1c9e7SStefano Zampini ierr = VecPointwiseDivide(is->y,is->y,is->counter);CHKERRQ(ierr); 593fd1c9e7SStefano Zampini ierr = MatDiagonalSet(is->A,is->y,insmode);CHKERRQ(ierr); 603fd1c9e7SStefano Zampini PetscFunctionReturn(0); 613fd1c9e7SStefano Zampini } 623fd1c9e7SStefano Zampini 633fd1c9e7SStefano Zampini #undef __FUNCT__ 643fd1c9e7SStefano Zampini #define __FUNCT__ "MatShift_IS" 653fd1c9e7SStefano Zampini PetscErrorCode MatShift_IS(Mat A,PetscScalar a) 663fd1c9e7SStefano Zampini { 673fd1c9e7SStefano Zampini PetscErrorCode ierr; 683fd1c9e7SStefano Zampini 693fd1c9e7SStefano Zampini PetscFunctionBegin; 703fd1c9e7SStefano Zampini ierr = MatDiagonalSet_IS(A,NULL,ADD_VALUES);CHKERRQ(ierr); 713fd1c9e7SStefano Zampini PetscFunctionReturn(0); 723fd1c9e7SStefano Zampini } 733fd1c9e7SStefano Zampini 743fd1c9e7SStefano Zampini #undef __FUNCT__ 75f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesLocal_SubMat_IS" 76f26d0771SStefano Zampini static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 77f26d0771SStefano Zampini { 78f26d0771SStefano Zampini PetscErrorCode ierr; 79f26d0771SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 80f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 81f26d0771SStefano Zampini 82f26d0771SStefano Zampini PetscFunctionBegin; 83f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 84f26d0771SStefano 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); 85f26d0771SStefano Zampini #endif 86f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 87f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 88f26d0771SStefano Zampini ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 89f26d0771SStefano Zampini PetscFunctionReturn(0); 90f26d0771SStefano Zampini } 91f26d0771SStefano Zampini 92f26d0771SStefano Zampini #undef __FUNCT__ 93f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesBlockedLocal_SubMat_IS" 94f26d0771SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 95f26d0771SStefano Zampini { 96f26d0771SStefano Zampini PetscErrorCode ierr; 97f26d0771SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 98f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 99f26d0771SStefano Zampini 100f26d0771SStefano Zampini PetscFunctionBegin; 101f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 102f26d0771SStefano 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); 103f26d0771SStefano Zampini #endif 104f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 105f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 106f26d0771SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 107f26d0771SStefano Zampini PetscFunctionReturn(0); 108f26d0771SStefano Zampini } 109f26d0771SStefano Zampini 110f26d0771SStefano Zampini #undef __FUNCT__ 111a8116848SStefano Zampini #define __FUNCT__ "PetscLayoutMapLocal_Private" 112f0ae7da4SStefano Zampini static PetscErrorCode PetscLayoutMapLocal_Private(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs) 113a8116848SStefano Zampini { 114a8116848SStefano Zampini PetscInt *owners = map->range; 115a8116848SStefano Zampini PetscInt n = map->n; 116a8116848SStefano Zampini PetscSF sf; 117a8116848SStefano Zampini PetscInt *lidxs,*work = NULL; 118a8116848SStefano Zampini PetscSFNode *ridxs; 119a8116848SStefano Zampini PetscMPIInt rank; 120a8116848SStefano Zampini PetscInt r, p = 0, len = 0; 121a8116848SStefano Zampini PetscErrorCode ierr; 122a8116848SStefano Zampini 123a8116848SStefano Zampini PetscFunctionBegin; 124a8116848SStefano Zampini /* Create SF where leaves are input idxs and roots are owned idxs (code adapted from MatZeroRowsMapLocal_Private) */ 125a8116848SStefano Zampini ierr = MPI_Comm_rank(map->comm,&rank);CHKERRQ(ierr); 126a8116848SStefano Zampini ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr); 127a8116848SStefano Zampini for (r = 0; r < n; ++r) lidxs[r] = -1; 128a8116848SStefano Zampini ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr); 129a8116848SStefano Zampini for (r = 0; r < N; ++r) { 130a8116848SStefano Zampini const PetscInt idx = idxs[r]; 131a8116848SStefano 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); 132a8116848SStefano Zampini if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 133a8116848SStefano Zampini ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr); 134a8116848SStefano Zampini } 135a8116848SStefano Zampini ridxs[r].rank = p; 136a8116848SStefano Zampini ridxs[r].index = idxs[r] - owners[p]; 137a8116848SStefano Zampini } 138a8116848SStefano Zampini ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr); 139a8116848SStefano Zampini ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr); 140a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 141a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 142f0ae7da4SStefano Zampini if (ogidxs) { /* communicate global idxs */ 143a8116848SStefano Zampini PetscInt cum = 0,start,*work2; 144f0ae7da4SStefano Zampini 145f0ae7da4SStefano Zampini ierr = PetscMalloc1(n,&work);CHKERRQ(ierr); 146a8116848SStefano Zampini ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr); 147a8116848SStefano Zampini for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++; 148a8116848SStefano Zampini ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRQ(ierr); 149a8116848SStefano Zampini start -= cum; 150a8116848SStefano Zampini cum = 0; 151a8116848SStefano Zampini for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++; 152a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 153a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 154a8116848SStefano Zampini ierr = PetscFree(work2);CHKERRQ(ierr); 155a8116848SStefano Zampini } 156a8116848SStefano Zampini ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 157a8116848SStefano Zampini /* Compress and put in indices */ 158a8116848SStefano Zampini for (r = 0; r < n; ++r) 159a8116848SStefano Zampini if (lidxs[r] >= 0) { 160a8116848SStefano Zampini if (work) work[len] = work[r]; 161a8116848SStefano Zampini lidxs[len++] = r; 162a8116848SStefano Zampini } 163a8116848SStefano Zampini if (on) *on = len; 164a8116848SStefano Zampini if (oidxs) *oidxs = lidxs; 165a8116848SStefano Zampini if (ogidxs) *ogidxs = work; 166a8116848SStefano Zampini PetscFunctionReturn(0); 167a8116848SStefano Zampini } 168a8116848SStefano Zampini 169a8116848SStefano Zampini #undef __FUNCT__ 170a8116848SStefano Zampini #define __FUNCT__ "MatGetSubMatrix_IS" 171a8116848SStefano Zampini static PetscErrorCode MatGetSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat) 172a8116848SStefano Zampini { 173a8116848SStefano Zampini Mat locmat,newlocmat; 174a8116848SStefano Zampini Mat_IS *newmatis; 175a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG) 176a8116848SStefano Zampini Vec rtest,ltest; 177a8116848SStefano Zampini const PetscScalar *array; 178a8116848SStefano Zampini #endif 179a8116848SStefano Zampini const PetscInt *idxs; 180a8116848SStefano Zampini PetscInt i,m,n; 181a8116848SStefano Zampini PetscErrorCode ierr; 182a8116848SStefano Zampini 183a8116848SStefano Zampini PetscFunctionBegin; 184a8116848SStefano Zampini if (scall == MAT_REUSE_MATRIX) { 185a8116848SStefano Zampini PetscBool ismatis; 186a8116848SStefano Zampini 187a8116848SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr); 188a8116848SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type"); 189a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 190a8116848SStefano Zampini if (!newmatis->getsub_ris) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS"); 191a8116848SStefano Zampini if (!newmatis->getsub_cis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS"); 192a8116848SStefano Zampini } 193a8116848SStefano Zampini /* irow and icol may not have duplicate entries */ 194a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG) 195a8116848SStefano Zampini ierr = MatCreateVecs(mat,<est,&rtest);CHKERRQ(ierr); 196a8116848SStefano Zampini ierr = ISGetLocalSize(irow,&n);CHKERRQ(ierr); 197a8116848SStefano Zampini ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 198a8116848SStefano Zampini for (i=0;i<n;i++) { 199a8116848SStefano Zampini ierr = VecSetValue(rtest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 200a8116848SStefano Zampini } 201a8116848SStefano Zampini ierr = VecAssemblyBegin(rtest);CHKERRQ(ierr); 202a8116848SStefano Zampini ierr = VecAssemblyEnd(rtest);CHKERRQ(ierr); 203a8116848SStefano Zampini ierr = VecGetLocalSize(rtest,&n);CHKERRQ(ierr); 204a8116848SStefano Zampini ierr = VecGetOwnershipRange(rtest,&m,NULL);CHKERRQ(ierr); 205a8116848SStefano Zampini ierr = VecGetArrayRead(rtest,&array);CHKERRQ(ierr); 206fd479f66SMatthew 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])); 207a8116848SStefano Zampini ierr = VecRestoreArrayRead(rtest,&array);CHKERRQ(ierr); 208a8116848SStefano Zampini ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 209a8116848SStefano Zampini ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 210a8116848SStefano Zampini ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 211a8116848SStefano Zampini for (i=0;i<n;i++) { 212a8116848SStefano Zampini ierr = VecSetValue(ltest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 213a8116848SStefano Zampini } 214a8116848SStefano Zampini ierr = VecAssemblyBegin(ltest);CHKERRQ(ierr); 215a8116848SStefano Zampini ierr = VecAssemblyEnd(ltest);CHKERRQ(ierr); 216a8116848SStefano Zampini ierr = VecGetLocalSize(ltest,&n);CHKERRQ(ierr); 217a8116848SStefano Zampini ierr = VecGetOwnershipRange(ltest,&m,NULL);CHKERRQ(ierr); 218a8116848SStefano Zampini ierr = VecGetArrayRead(ltest,&array);CHKERRQ(ierr); 219fd479f66SMatthew 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])); 220a8116848SStefano Zampini ierr = VecRestoreArrayRead(ltest,&array);CHKERRQ(ierr); 221a8116848SStefano Zampini ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 222a8116848SStefano Zampini ierr = VecDestroy(&rtest);CHKERRQ(ierr); 223a8116848SStefano Zampini ierr = VecDestroy(<est);CHKERRQ(ierr); 224a8116848SStefano Zampini #endif 225a8116848SStefano Zampini if (scall == MAT_INITIAL_MATRIX) { 226a8116848SStefano Zampini Mat_IS *matis = (Mat_IS*)mat->data; 227a8116848SStefano Zampini ISLocalToGlobalMapping rl2g; 228a8116848SStefano Zampini IS is; 229a8116848SStefano Zampini PetscInt *lidxs,*lgidxs,*newgidxs; 2306dd40735SStefano Zampini PetscInt ll,newloc; 231a8116848SStefano Zampini MPI_Comm comm; 232a8116848SStefano Zampini 233a8116848SStefano Zampini ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 234a8116848SStefano Zampini ierr = ISGetLocalSize(irow,&m);CHKERRQ(ierr); 235a8116848SStefano Zampini ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 236a8116848SStefano Zampini ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 237a8116848SStefano Zampini ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr); 238a8116848SStefano Zampini ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 239a8116848SStefano Zampini /* communicate irow to their owners in the layout */ 240a8116848SStefano Zampini ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 241f0ae7da4SStefano Zampini ierr = PetscLayoutMapLocal_Private(mat->rmap,m,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 242a8116848SStefano Zampini ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 243*3d996552SStefano Zampini ierr = MatISSetUpSF(mat);CHKERRQ(ierr); 244*3d996552SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,matis->sf->nroots*sizeof(PetscInt));CHKERRQ(ierr); 245a8116848SStefano Zampini for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1; 246a8116848SStefano Zampini ierr = PetscFree(lidxs);CHKERRQ(ierr); 247a8116848SStefano Zampini ierr = PetscFree(lgidxs);CHKERRQ(ierr); 248a8116848SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 249a8116848SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 250*3d996552SStefano Zampini for (i=0,newloc=0;i<matis->sf->nleaves;i++) if (matis->sf_leafdata[i]) newloc++; 251a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 252a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 253*3d996552SStefano Zampini for (i=0,newloc=0;i<matis->sf->nleaves;i++) 254a8116848SStefano Zampini if (matis->sf_leafdata[i]) { 255a8116848SStefano Zampini lidxs[newloc] = i; 256a8116848SStefano Zampini newgidxs[newloc++] = matis->sf_leafdata[i]-1; 257a8116848SStefano Zampini } 258a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 259a8116848SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 260a8116848SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 261a8116848SStefano Zampini /* local is to extract local submatrix */ 262a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 263a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris);CHKERRQ(ierr); 264a8116848SStefano Zampini if (mat->congruentlayouts == -1) { /* first time we compare rows and cols layouts */ 265a8116848SStefano Zampini PetscBool cong; 266a8116848SStefano Zampini ierr = PetscLayoutCompare(mat->rmap,mat->cmap,&cong);CHKERRQ(ierr); 267a8116848SStefano Zampini if (cong) mat->congruentlayouts = 1; 268a8116848SStefano Zampini else mat->congruentlayouts = 0; 269a8116848SStefano Zampini } 270a8116848SStefano Zampini if (mat->congruentlayouts && irow == icol) { 271a8116848SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g);CHKERRQ(ierr); 272a8116848SStefano Zampini ierr = PetscObjectReference((PetscObject)newmatis->getsub_ris);CHKERRQ(ierr); 273a8116848SStefano Zampini newmatis->getsub_cis = newmatis->getsub_ris; 274a8116848SStefano Zampini } else { 275a8116848SStefano Zampini ISLocalToGlobalMapping cl2g; 276a8116848SStefano Zampini 277a8116848SStefano Zampini /* communicate icol to their owners in the layout */ 278a8116848SStefano Zampini ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 279f0ae7da4SStefano Zampini ierr = PetscLayoutMapLocal_Private(mat->cmap,n,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 280a8116848SStefano Zampini ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 281*3d996552SStefano Zampini ierr = PetscMemzero(matis->csf_rootdata,matis->csf->nroots*sizeof(PetscInt));CHKERRQ(ierr); 282a8116848SStefano Zampini for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1; 283a8116848SStefano Zampini ierr = PetscFree(lidxs);CHKERRQ(ierr); 284a8116848SStefano Zampini ierr = PetscFree(lgidxs);CHKERRQ(ierr); 285a8116848SStefano Zampini ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 286a8116848SStefano Zampini ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 287*3d996552SStefano Zampini for (i=0,newloc=0;i<matis->csf->nleaves;i++) if (matis->csf_leafdata[i]) newloc++; 288a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 289a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 290*3d996552SStefano Zampini for (i=0,newloc=0;i<matis->csf->nleaves;i++) 291a8116848SStefano Zampini if (matis->csf_leafdata[i]) { 292a8116848SStefano Zampini lidxs[newloc] = i; 293a8116848SStefano Zampini newgidxs[newloc++] = matis->csf_leafdata[i]-1; 294a8116848SStefano Zampini } 295a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 296a8116848SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 297a8116848SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 298a8116848SStefano Zampini /* local is to extract local submatrix */ 299a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis);CHKERRQ(ierr); 300a8116848SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr); 301a8116848SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 302a8116848SStefano Zampini } 303a8116848SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 304a8116848SStefano Zampini } else { 305a8116848SStefano Zampini ierr = MatISGetLocalMat(*newmat,&newlocmat);CHKERRQ(ierr); 306a8116848SStefano Zampini } 307a8116848SStefano Zampini ierr = MatISGetLocalMat(mat,&locmat);CHKERRQ(ierr); 308a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 309a8116848SStefano Zampini ierr = MatGetSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat);CHKERRQ(ierr); 310a8116848SStefano Zampini if (scall == MAT_INITIAL_MATRIX) { 311a8116848SStefano Zampini ierr = MatISSetLocalMat(*newmat,newlocmat);CHKERRQ(ierr); 312a8116848SStefano Zampini ierr = MatDestroy(&newlocmat);CHKERRQ(ierr); 313a8116848SStefano Zampini } 314a8116848SStefano Zampini ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 315a8116848SStefano Zampini ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 316a8116848SStefano Zampini PetscFunctionReturn(0); 317a8116848SStefano Zampini } 318a8116848SStefano Zampini 319a8116848SStefano Zampini #undef __FUNCT__ 3202b404112SStefano Zampini #define __FUNCT__ "MatCopy_IS" 321a8116848SStefano Zampini static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str) 3222b404112SStefano Zampini { 3232b404112SStefano Zampini Mat_IS *a = (Mat_IS*)A->data,*b; 3242b404112SStefano Zampini PetscBool ismatis; 3252b404112SStefano Zampini PetscErrorCode ierr; 3262b404112SStefano Zampini 3272b404112SStefano Zampini PetscFunctionBegin; 3282b404112SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr); 3292b404112SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented"); 3302b404112SStefano Zampini b = (Mat_IS*)B->data; 3312b404112SStefano Zampini ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr); 3322b404112SStefano Zampini PetscFunctionReturn(0); 3332b404112SStefano Zampini } 3342b404112SStefano Zampini 3352b404112SStefano Zampini #undef __FUNCT__ 3366bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS" 337a8116848SStefano Zampini static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool *missing,PetscInt *d) 3386bd84002SStefano Zampini { 339527b2640SStefano Zampini Vec v; 340527b2640SStefano Zampini const PetscScalar *array; 341527b2640SStefano Zampini PetscInt i,n; 3426bd84002SStefano Zampini PetscErrorCode ierr; 3436bd84002SStefano Zampini 3446bd84002SStefano Zampini PetscFunctionBegin; 345527b2640SStefano Zampini *missing = PETSC_FALSE; 346527b2640SStefano Zampini ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr); 347527b2640SStefano Zampini ierr = MatGetDiagonal(A,v);CHKERRQ(ierr); 348527b2640SStefano Zampini ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 349527b2640SStefano Zampini ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr); 350527b2640SStefano Zampini for (i=0;i<n;i++) if (array[i] == 0.) break; 351527b2640SStefano Zampini ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr); 352527b2640SStefano Zampini ierr = VecDestroy(&v);CHKERRQ(ierr); 353527b2640SStefano Zampini if (i != n) *missing = PETSC_TRUE; 354527b2640SStefano Zampini if (d) { 355527b2640SStefano Zampini *d = -1; 356527b2640SStefano Zampini if (*missing) { 357527b2640SStefano Zampini PetscInt rstart; 358527b2640SStefano Zampini ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 359527b2640SStefano Zampini *d = i+rstart; 360527b2640SStefano Zampini } 361527b2640SStefano Zampini } 3626bd84002SStefano Zampini PetscFunctionReturn(0); 3636bd84002SStefano Zampini } 3646bd84002SStefano Zampini 3656bd84002SStefano Zampini #undef __FUNCT__ 366cf0a3239SStefano Zampini #define __FUNCT__ "MatISSetUpSF_IS" 367cf0a3239SStefano Zampini static PetscErrorCode MatISSetUpSF_IS(Mat B) 36828f4e0baSStefano Zampini { 36928f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 37028f4e0baSStefano Zampini const PetscInt *gidxs; 3714f2d7cafSStefano Zampini PetscInt nleaves; 37228f4e0baSStefano Zampini PetscErrorCode ierr; 37328f4e0baSStefano Zampini 37428f4e0baSStefano Zampini PetscFunctionBegin; 3754f2d7cafSStefano Zampini if (matis->sf) PetscFunctionReturn(0); 37628f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 3773bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 3784f2d7cafSStefano Zampini ierr = ISLocalToGlobalMappingGetSize(B->rmap->mapping,&nleaves);CHKERRQ(ierr); 3794f2d7cafSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 3803bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 3814f2d7cafSStefano Zampini ierr = PetscMalloc2(matis->sf->nroots,&matis->sf_rootdata,matis->sf->nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 382a8116848SStefano Zampini if (B->rmap->mapping != B->cmap->mapping) { /* setup SF for columns */ 383*3d996552SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(B->cmap->mapping,&nleaves);CHKERRQ(ierr); 384a8116848SStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf);CHKERRQ(ierr); 385a8116848SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 386*3d996552SStefano Zampini ierr = PetscSFSetGraphLayout(matis->csf,B->cmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 387a8116848SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 388*3d996552SStefano Zampini ierr = PetscMalloc2(matis->csf->nroots,&matis->csf_rootdata,matis->csf->nleaves,&matis->csf_leafdata);CHKERRQ(ierr); 389a8116848SStefano Zampini } else { 390a8116848SStefano Zampini matis->csf = matis->sf; 391a8116848SStefano Zampini matis->csf_leafdata = matis->sf_leafdata; 392a8116848SStefano Zampini matis->csf_rootdata = matis->sf_rootdata; 393a8116848SStefano Zampini } 39428f4e0baSStefano Zampini PetscFunctionReturn(0); 39528f4e0baSStefano Zampini } 3962e1947a5SStefano Zampini 3972e1947a5SStefano Zampini #undef __FUNCT__ 3982e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 399eb82efa4SStefano Zampini /*@ 400a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 401a88811baSStefano Zampini 402a88811baSStefano Zampini Collective on MPI_Comm 403a88811baSStefano Zampini 404a88811baSStefano Zampini Input Parameters: 405a88811baSStefano Zampini + B - the matrix 406a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 407a88811baSStefano Zampini (same value is used for all local rows) 408a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 409a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 410a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 411a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 412a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 413a88811baSStefano Zampini the diagonal entry even if it is zero. 414a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 415a88811baSStefano Zampini submatrix (same value is used for all local rows). 416a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 417a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 418a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 419a88811baSStefano Zampini structure. The size of this array is equal to the number 420a88811baSStefano Zampini of local rows, i.e 'm'. 421a88811baSStefano Zampini 422a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 423a88811baSStefano Zampini 424a88811baSStefano Zampini Level: intermediate 425a88811baSStefano Zampini 426a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 427a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 428a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 429a88811baSStefano Zampini 430a88811baSStefano Zampini .keywords: matrix 431a88811baSStefano Zampini 4323c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS 433a88811baSStefano Zampini @*/ 4342e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 4352e1947a5SStefano Zampini { 4362e1947a5SStefano Zampini PetscErrorCode ierr; 4372e1947a5SStefano Zampini 4382e1947a5SStefano Zampini PetscFunctionBegin; 4392e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 4402e1947a5SStefano Zampini PetscValidType(B,1); 4412e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 4422e1947a5SStefano Zampini PetscFunctionReturn(0); 4432e1947a5SStefano Zampini } 4442e1947a5SStefano Zampini 4452e1947a5SStefano Zampini #undef __FUNCT__ 4462e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 4477230de76SStefano Zampini static PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 4482e1947a5SStefano Zampini { 4492e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 45028f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 4512e1947a5SStefano Zampini PetscErrorCode ierr; 4522e1947a5SStefano Zampini 4532e1947a5SStefano Zampini PetscFunctionBegin; 4546c4ed002SBarry Smith if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 455cf0a3239SStefano Zampini ierr = MatISSetUpSF(B);CHKERRQ(ierr); 4564f2d7cafSStefano Zampini 4574f2d7cafSStefano Zampini if (!d_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nz; 4584f2d7cafSStefano Zampini else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 4594f2d7cafSStefano Zampini 4604f2d7cafSStefano Zampini if (!o_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nz; 4614f2d7cafSStefano Zampini else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 4624f2d7cafSStefano Zampini 46328f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 46428f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 46528f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 46628f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 4674f2d7cafSStefano Zampini 4684f2d7cafSStefano Zampini for (i=0;i<matis->sf->nleaves;i++) matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 46928f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 4704f2d7cafSStefano Zampini 4714f2d7cafSStefano Zampini for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 47228f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 4734f2d7cafSStefano Zampini 4744f2d7cafSStefano Zampini for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 47528f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 4762e1947a5SStefano Zampini PetscFunctionReturn(0); 4772e1947a5SStefano Zampini } 478b4319ba4SBarry Smith 479b4319ba4SBarry Smith #undef __FUNCT__ 4803927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private" 4813927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 4823927de2eSStefano Zampini { 4833927de2eSStefano Zampini Mat_IS *matis = (Mat_IS*)(A->data); 4843927de2eSStefano Zampini PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 485ecf5a873SStefano Zampini const PetscInt *global_indices_r,*global_indices_c; 4863927de2eSStefano Zampini PetscInt i,j,bs,rows,cols; 4873927de2eSStefano Zampini PetscInt lrows,lcols; 4883927de2eSStefano Zampini PetscInt local_rows,local_cols; 4893927de2eSStefano Zampini PetscMPIInt nsubdomains; 4903927de2eSStefano Zampini PetscBool isdense,issbaij; 4913927de2eSStefano Zampini PetscErrorCode ierr; 4923927de2eSStefano Zampini 4933927de2eSStefano Zampini PetscFunctionBegin; 4943927de2eSStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr); 4953927de2eSStefano Zampini ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr); 4963927de2eSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 4973927de2eSStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 4983927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 4993927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 500ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 501ecf5a873SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 5027230de76SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 503ecf5a873SStefano Zampini } else { 504ecf5a873SStefano Zampini global_indices_c = global_indices_r; 505ecf5a873SStefano Zampini } 506ecf5a873SStefano Zampini 5073927de2eSStefano Zampini if (issbaij) { 5083927de2eSStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 5093927de2eSStefano Zampini } 5103927de2eSStefano Zampini /* 511ecf5a873SStefano Zampini An SF reduce is needed to sum up properly on shared rows. 5123927de2eSStefano Zampini Note that generally preallocation is not exact, since it overestimates nonzeros 5133927de2eSStefano Zampini */ 514cf0a3239SStefano Zampini ierr = MatISSetUpSF(A);CHKERRQ(ierr); 5153927de2eSStefano Zampini ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr); 5163927de2eSStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr); 5173927de2eSStefano Zampini /* All processes need to compute entire row ownership */ 5183927de2eSStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 5193927de2eSStefano Zampini ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 5203927de2eSStefano Zampini for (i=0;i<nsubdomains;i++) { 5213927de2eSStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 5223927de2eSStefano Zampini row_ownership[j] = i; 5233927de2eSStefano Zampini } 5243927de2eSStefano Zampini } 5257230de76SStefano Zampini ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 5263927de2eSStefano Zampini 5273927de2eSStefano Zampini /* 5283927de2eSStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 5293927de2eSStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 5303927de2eSStefano Zampini */ 5313927de2eSStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 5323927de2eSStefano Zampini /* preallocation as a MATAIJ */ 5333927de2eSStefano Zampini if (isdense) { /* special case for dense local matrices */ 5343927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 535ecf5a873SStefano Zampini PetscInt index_row = global_indices_r[i]; 5363927de2eSStefano Zampini for (j=i;j<local_rows;j++) { 5373927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 538ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[j]; 5393927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 5403927de2eSStefano Zampini my_dnz[i] += 1; 5413927de2eSStefano Zampini } else { /* offdiag block */ 5423927de2eSStefano Zampini my_onz[i] += 1; 5433927de2eSStefano Zampini } 5443927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 5453927de2eSStefano Zampini if (i != j) { 5463927de2eSStefano Zampini owner = row_ownership[index_col]; 5473927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 5483927de2eSStefano Zampini my_dnz[j] += 1; 5493927de2eSStefano Zampini } else { 5503927de2eSStefano Zampini my_onz[j] += 1; 5513927de2eSStefano Zampini } 5523927de2eSStefano Zampini } 5533927de2eSStefano Zampini } 5543927de2eSStefano Zampini } 555bb1015c3SStefano Zampini } else if (matis->A->ops->getrowij) { 556bb1015c3SStefano Zampini const PetscInt *ii,*jj,*jptr; 557bb1015c3SStefano Zampini PetscBool done; 558bb1015c3SStefano Zampini ierr = MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr); 559bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__); 560bb1015c3SStefano Zampini jptr = jj; 561bb1015c3SStefano Zampini for (i=0;i<local_rows;i++) { 562bb1015c3SStefano Zampini PetscInt index_row = global_indices_r[i]; 563bb1015c3SStefano Zampini for (j=0;j<ii[i+1]-ii[i];j++,jptr++) { 564bb1015c3SStefano Zampini PetscInt owner = row_ownership[index_row]; 565bb1015c3SStefano Zampini PetscInt index_col = global_indices_c[*jptr]; 566bb1015c3SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 567bb1015c3SStefano Zampini my_dnz[i] += 1; 568bb1015c3SStefano Zampini } else { /* offdiag block */ 569bb1015c3SStefano Zampini my_onz[i] += 1; 570bb1015c3SStefano Zampini } 571bb1015c3SStefano Zampini /* same as before, interchanging rows and cols */ 572bb1015c3SStefano Zampini if (issbaij && index_col != index_row) { 573bb1015c3SStefano Zampini owner = row_ownership[index_col]; 574bb1015c3SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 575bb1015c3SStefano Zampini my_dnz[*jptr] += 1; 576bb1015c3SStefano Zampini } else { 577bb1015c3SStefano Zampini my_onz[*jptr] += 1; 578bb1015c3SStefano Zampini } 579bb1015c3SStefano Zampini } 580bb1015c3SStefano Zampini } 581bb1015c3SStefano Zampini } 582bb1015c3SStefano Zampini ierr = MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr); 583bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__); 584bb1015c3SStefano Zampini } else { /* loop over rows and use MatGetRow */ 5853927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 5863927de2eSStefano Zampini const PetscInt *cols; 587ecf5a873SStefano Zampini PetscInt ncols,index_row = global_indices_r[i]; 5883927de2eSStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 5893927de2eSStefano Zampini for (j=0;j<ncols;j++) { 5903927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 591ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[cols[j]]; 5923927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 5933927de2eSStefano Zampini my_dnz[i] += 1; 5943927de2eSStefano Zampini } else { /* offdiag block */ 5953927de2eSStefano Zampini my_onz[i] += 1; 5963927de2eSStefano Zampini } 5973927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 598d9a9e74cSStefano Zampini if (issbaij && index_col != index_row) { 5993927de2eSStefano Zampini owner = row_ownership[index_col]; 6003927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 601d9a9e74cSStefano Zampini my_dnz[cols[j]] += 1; 6023927de2eSStefano Zampini } else { 603d9a9e74cSStefano Zampini my_onz[cols[j]] += 1; 6043927de2eSStefano Zampini } 6053927de2eSStefano Zampini } 6063927de2eSStefano Zampini } 6073927de2eSStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 6083927de2eSStefano Zampini } 6093927de2eSStefano Zampini } 610ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 611ecf5a873SStefano Zampini if (global_indices_c != global_indices_r) { 6127230de76SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 613ecf5a873SStefano Zampini } 6143927de2eSStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 615ecf5a873SStefano Zampini 616ecf5a873SStefano Zampini /* Reduce my_dnz and my_onz */ 6173927de2eSStefano Zampini if (maxreduce) { 6183927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 6193927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 620bb1015c3SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 6213927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 6223927de2eSStefano Zampini } else { 6233927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 6243927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 625bb1015c3SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 6263927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 6273927de2eSStefano Zampini } 6283927de2eSStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 6293927de2eSStefano Zampini 6303927de2eSStefano Zampini /* Resize preallocation if overestimated */ 6313927de2eSStefano Zampini for (i=0;i<lrows;i++) { 6323927de2eSStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 6333927de2eSStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 6343927de2eSStefano Zampini } 6353927de2eSStefano Zampini /* set preallocation */ 6363927de2eSStefano Zampini ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 6373927de2eSStefano Zampini for (i=0;i<lrows/bs;i++) { 6383927de2eSStefano Zampini dnz[i] = dnz[i*bs]/bs; 6393927de2eSStefano Zampini onz[i] = onz[i*bs]/bs; 6403927de2eSStefano Zampini } 6413927de2eSStefano Zampini ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 6423927de2eSStefano Zampini ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 6433927de2eSStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 6443927de2eSStefano Zampini if (issbaij) { 6453927de2eSStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 6463927de2eSStefano Zampini } 6473927de2eSStefano Zampini PetscFunctionReturn(0); 6483927de2eSStefano Zampini } 6493927de2eSStefano Zampini 6503927de2eSStefano Zampini #undef __FUNCT__ 651b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 6527230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 653b7ce53b6SStefano Zampini { 654b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 655d9a9e74cSStefano Zampini Mat local_mat; 656b7ce53b6SStefano Zampini /* info on mat */ 6573cfa4ea4SStefano Zampini PetscInt bs,rows,cols,lrows,lcols; 658b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 659686e3a49SStefano Zampini PetscBool isdense,issbaij,isseqaij; 6607c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 661b7ce53b6SStefano Zampini /* values insertion */ 662b7ce53b6SStefano Zampini PetscScalar *array; 663b7ce53b6SStefano Zampini /* work */ 664b7ce53b6SStefano Zampini PetscErrorCode ierr; 665b7ce53b6SStefano Zampini 666b7ce53b6SStefano Zampini PetscFunctionBegin; 667b7ce53b6SStefano Zampini /* get info from mat */ 6687c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 6697c03b4e8SStefano Zampini if (nsubdomains == 1) { 6707c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 6717c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 6727c03b4e8SStefano Zampini } else { 6737c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 6747c03b4e8SStefano Zampini } 6757c03b4e8SStefano Zampini PetscFunctionReturn(0); 6767c03b4e8SStefano Zampini } 677b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 678b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 6793cfa4ea4SStefano Zampini ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 680b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 681b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 682686e3a49SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 683b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 684b7ce53b6SStefano Zampini 685b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 686b7ce53b6SStefano Zampini MatType new_mat_type; 6873927de2eSStefano Zampini PetscBool issbaij_red; 688b7ce53b6SStefano Zampini 689b7ce53b6SStefano Zampini /* determining new matrix type */ 690b2566f29SBarry Smith ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 691b7ce53b6SStefano Zampini if (issbaij_red) { 692b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 693b7ce53b6SStefano Zampini } else { 694b7ce53b6SStefano Zampini if (bs>1) { 695b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 696b7ce53b6SStefano Zampini } else { 697b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 698b7ce53b6SStefano Zampini } 699b7ce53b6SStefano Zampini } 700b7ce53b6SStefano Zampini 7013927de2eSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 7023cfa4ea4SStefano Zampini ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 7033927de2eSStefano Zampini ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 7043927de2eSStefano Zampini ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr); 7053927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 706b7ce53b6SStefano Zampini } else { 7073cfa4ea4SStefano Zampini PetscInt mbs,mrows,mcols,mlrows,mlcols; 708b7ce53b6SStefano Zampini /* some checks */ 709b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 710b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 7113cfa4ea4SStefano Zampini ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 7126c4ed002SBarry Smith if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 7136c4ed002SBarry Smith if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 7146c4ed002SBarry Smith if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 7156c4ed002SBarry Smith if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 7166c4ed002SBarry Smith if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 717b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 718b7ce53b6SStefano Zampini } 719d9a9e74cSStefano Zampini 720d9a9e74cSStefano Zampini if (issbaij) { 721d9a9e74cSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 722d9a9e74cSStefano Zampini } else { 723d9a9e74cSStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 724d9a9e74cSStefano Zampini local_mat = matis->A; 725d9a9e74cSStefano Zampini } 726686e3a49SStefano Zampini 727b7ce53b6SStefano Zampini /* Set values */ 728ecf5a873SStefano Zampini ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 729b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 730ecf5a873SStefano Zampini PetscInt i,*dummy_rows,*dummy_cols; 731ecf5a873SStefano Zampini 732ecf5a873SStefano Zampini if (local_rows != local_cols) { 733ecf5a873SStefano Zampini ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 734ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 735ecf5a873SStefano Zampini for (i=0;i<local_cols;i++) dummy_cols[i] = i; 736ecf5a873SStefano Zampini } else { 737ecf5a873SStefano Zampini ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 738ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 739ecf5a873SStefano Zampini dummy_cols = dummy_rows; 740ecf5a873SStefano Zampini } 741b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 742d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 743ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 744d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 745ecf5a873SStefano Zampini if (dummy_rows != dummy_cols) { 746ecf5a873SStefano Zampini ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 747ecf5a873SStefano Zampini } else { 748ecf5a873SStefano Zampini ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 749ecf5a873SStefano Zampini } 750686e3a49SStefano Zampini } else if (isseqaij) { 751ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 752686e3a49SStefano Zampini PetscBool done; 753686e3a49SStefano Zampini 754d9a9e74cSStefano Zampini ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 755bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__); 756d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 757686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 758ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 759686e3a49SStefano Zampini } 760d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 761bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__); 762d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 763686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 764ecf5a873SStefano Zampini PetscInt i; 765c0962df8SStefano Zampini 766686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 767686e3a49SStefano Zampini PetscInt j; 768ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 769686e3a49SStefano Zampini 770ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 771ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 772ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 773686e3a49SStefano Zampini } 774b7ce53b6SStefano Zampini } 775b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 776d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 777b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 778b7ce53b6SStefano Zampini if (isdense) { 779b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 780b7ce53b6SStefano Zampini } 781b7ce53b6SStefano Zampini PetscFunctionReturn(0); 782b7ce53b6SStefano Zampini } 783b7ce53b6SStefano Zampini 784b7ce53b6SStefano Zampini #undef __FUNCT__ 785b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 786b7ce53b6SStefano Zampini /*@ 787b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 788b7ce53b6SStefano Zampini 789b7ce53b6SStefano Zampini Input Parameter: 790b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 791b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 792b7ce53b6SStefano Zampini 793b7ce53b6SStefano Zampini Output Parameter: 794b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 795b7ce53b6SStefano Zampini 796b7ce53b6SStefano Zampini Level: developer 797b7ce53b6SStefano Zampini 798eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 799b7ce53b6SStefano Zampini 800b7ce53b6SStefano Zampini .seealso: MATIS 801b7ce53b6SStefano Zampini @*/ 802b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 803b7ce53b6SStefano Zampini { 804b7ce53b6SStefano Zampini PetscErrorCode ierr; 805b7ce53b6SStefano Zampini 806b7ce53b6SStefano Zampini PetscFunctionBegin; 807b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 808b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 809b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 810b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 811b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 812b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 8136c4ed002SBarry Smith if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 814b7ce53b6SStefano Zampini } 815b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 816b7ce53b6SStefano Zampini PetscFunctionReturn(0); 817b7ce53b6SStefano Zampini } 818b7ce53b6SStefano Zampini 819b7ce53b6SStefano Zampini #undef __FUNCT__ 820ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 821ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 822ad6194a2SStefano Zampini { 823ad6194a2SStefano Zampini PetscErrorCode ierr; 824ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 825ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 826ad6194a2SStefano Zampini Mat B,localmat; 827ad6194a2SStefano Zampini 828ad6194a2SStefano Zampini PetscFunctionBegin; 829ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 830ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 831ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 832e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 833ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 834ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 835b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 836ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 837ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 838ad6194a2SStefano Zampini *newmat = B; 839ad6194a2SStefano Zampini PetscFunctionReturn(0); 840ad6194a2SStefano Zampini } 841ad6194a2SStefano Zampini 842ad6194a2SStefano Zampini #undef __FUNCT__ 84369796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 844a8116848SStefano Zampini static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 84569796d55SStefano Zampini { 84669796d55SStefano Zampini PetscErrorCode ierr; 84769796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 84869796d55SStefano Zampini PetscBool local_sym; 84969796d55SStefano Zampini 85069796d55SStefano Zampini PetscFunctionBegin; 85169796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 852b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 85369796d55SStefano Zampini PetscFunctionReturn(0); 85469796d55SStefano Zampini } 85569796d55SStefano Zampini 85669796d55SStefano Zampini #undef __FUNCT__ 85769796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 858a8116848SStefano Zampini static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 85969796d55SStefano Zampini { 86069796d55SStefano Zampini PetscErrorCode ierr; 86169796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 86269796d55SStefano Zampini PetscBool local_sym; 86369796d55SStefano Zampini 86469796d55SStefano Zampini PetscFunctionBegin; 86569796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 866b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 86769796d55SStefano Zampini PetscFunctionReturn(0); 86869796d55SStefano Zampini } 86969796d55SStefano Zampini 87069796d55SStefano Zampini #undef __FUNCT__ 871b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 872a8116848SStefano Zampini static PetscErrorCode MatDestroy_IS(Mat A) 873b4319ba4SBarry Smith { 874dfbe8321SBarry Smith PetscErrorCode ierr; 875b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 876b4319ba4SBarry Smith 877b4319ba4SBarry Smith PetscFunctionBegin; 8786bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 879e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 880e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 8816bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 8826bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 8833fd1c9e7SStefano Zampini ierr = VecDestroy(&b->counter);CHKERRQ(ierr); 884a8116848SStefano Zampini ierr = ISDestroy(&b->getsub_ris);CHKERRQ(ierr); 885a8116848SStefano Zampini ierr = ISDestroy(&b->getsub_cis);CHKERRQ(ierr); 886a8116848SStefano Zampini if (b->sf != b->csf) { 887a8116848SStefano Zampini ierr = PetscSFDestroy(&b->csf);CHKERRQ(ierr); 888a8116848SStefano Zampini ierr = PetscFree2(b->csf_rootdata,b->csf_leafdata);CHKERRQ(ierr); 889a8116848SStefano Zampini } 89028f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 89128f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 892bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 893dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 894bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 895b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 896b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 8972e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 898cf0a3239SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",NULL);CHKERRQ(ierr); 899b4319ba4SBarry Smith PetscFunctionReturn(0); 900b4319ba4SBarry Smith } 901b4319ba4SBarry Smith 902b4319ba4SBarry Smith #undef __FUNCT__ 903b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 904a8116848SStefano Zampini static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 905b4319ba4SBarry Smith { 906dfbe8321SBarry Smith PetscErrorCode ierr; 907b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 908b4319ba4SBarry Smith PetscScalar zero = 0.0; 909b4319ba4SBarry Smith 910b4319ba4SBarry Smith PetscFunctionBegin; 911b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 912e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 913e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 914b4319ba4SBarry Smith 915b4319ba4SBarry Smith /* multiply the local matrix */ 916b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 917b4319ba4SBarry Smith 918b4319ba4SBarry Smith /* scatter product back into global memory */ 9192dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 920e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 921e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 922b4319ba4SBarry Smith PetscFunctionReturn(0); 923b4319ba4SBarry Smith } 924b4319ba4SBarry Smith 925b4319ba4SBarry Smith #undef __FUNCT__ 9262e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 927a8116848SStefano Zampini static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 9282e74eeadSLisandro Dalcin { 929650997f4SStefano Zampini Vec temp_vec; 9302e74eeadSLisandro Dalcin PetscErrorCode ierr; 9312e74eeadSLisandro Dalcin 9322e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 933650997f4SStefano Zampini if (v3 != v2) { 934650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 935650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 936650997f4SStefano Zampini } else { 937650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 938650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 939650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 940650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 941650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 942650997f4SStefano Zampini } 9432e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9442e74eeadSLisandro Dalcin } 9452e74eeadSLisandro Dalcin 9462e74eeadSLisandro Dalcin #undef __FUNCT__ 9472e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 948a8116848SStefano Zampini static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 9492e74eeadSLisandro Dalcin { 9502e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9512e74eeadSLisandro Dalcin PetscErrorCode ierr; 9522e74eeadSLisandro Dalcin 953e176bc59SStefano Zampini PetscFunctionBegin; 9542e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 955e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 956e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9572e74eeadSLisandro Dalcin 9582e74eeadSLisandro Dalcin /* multiply the local matrix */ 959e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 9602e74eeadSLisandro Dalcin 9612e74eeadSLisandro Dalcin /* scatter product back into global vector */ 962e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 963e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 964e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9652e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9662e74eeadSLisandro Dalcin } 9672e74eeadSLisandro Dalcin 9682e74eeadSLisandro Dalcin #undef __FUNCT__ 9692e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 970a8116848SStefano Zampini static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 9712e74eeadSLisandro Dalcin { 972650997f4SStefano Zampini Vec temp_vec; 9732e74eeadSLisandro Dalcin PetscErrorCode ierr; 9742e74eeadSLisandro Dalcin 9752e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 976650997f4SStefano Zampini if (v3 != v2) { 977650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 978650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 979650997f4SStefano Zampini } else { 980650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 981650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 982650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 983650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 984650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 985650997f4SStefano Zampini } 9862e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9872e74eeadSLisandro Dalcin } 9882e74eeadSLisandro Dalcin 9892e74eeadSLisandro Dalcin #undef __FUNCT__ 990b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 991a8116848SStefano Zampini static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 992b4319ba4SBarry Smith { 993b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 994dfbe8321SBarry Smith PetscErrorCode ierr; 995b4319ba4SBarry Smith PetscViewer sviewer; 996ee2491ecSStefano Zampini PetscBool isascii,view = PETSC_TRUE; 997b4319ba4SBarry Smith 998b4319ba4SBarry Smith PetscFunctionBegin; 999ee2491ecSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 1000ee2491ecSStefano Zampini if (isascii) { 1001ee2491ecSStefano Zampini PetscViewerFormat format; 1002ee2491ecSStefano Zampini 1003ee2491ecSStefano Zampini ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1004ee2491ecSStefano Zampini if (format == PETSC_VIEWER_ASCII_INFO) view = PETSC_FALSE; 1005ee2491ecSStefano Zampini } 1006ee2491ecSStefano Zampini if (!view) PetscFunctionReturn(0); 10073f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 1008b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 10093f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 10106e520ac8SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1011b4319ba4SBarry Smith PetscFunctionReturn(0); 1012b4319ba4SBarry Smith } 1013b4319ba4SBarry Smith 1014b4319ba4SBarry Smith #undef __FUNCT__ 1015b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 1016a8116848SStefano Zampini static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 1017b4319ba4SBarry Smith { 1018dfbe8321SBarry Smith PetscErrorCode ierr; 1019e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 1020b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1021b4319ba4SBarry Smith IS from,to; 1022e176bc59SStefano Zampini Vec cglobal,rglobal; 1023b4319ba4SBarry Smith 1024b4319ba4SBarry Smith PetscFunctionBegin; 1025784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 1026e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 10273bbff08aSStefano Zampini /* Destroy any previous data */ 102870cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 102970cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 10303fd1c9e7SStefano Zampini ierr = VecDestroy(&is->counter);CHKERRQ(ierr); 1031e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 1032e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 10331c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 103428f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 103528f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 10363bbff08aSStefano Zampini 10373bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 1038fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1039fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1040fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 1041fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 1042b4319ba4SBarry Smith 1043b4319ba4SBarry Smith /* Create the local matrix A */ 1044e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 1045e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 1046e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 1047e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 1048f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 1049e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 1050e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 1051e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 1052ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 1053ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 1054b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 1055b4319ba4SBarry Smith 1056f26d0771SStefano Zampini if (!is->islocalref) { /* setup scatters and local vectors for MatMult */ 1057b4319ba4SBarry Smith /* Create the local work vectors */ 10583bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 1059b4319ba4SBarry Smith 1060e176bc59SStefano Zampini /* setup the global to local scatters */ 1061e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 1062e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 1063784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 1064e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 1065e176bc59SStefano Zampini if (rmapping != cmapping) { 1066e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 1067e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 1068e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 1069e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 1070e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 1071e176bc59SStefano Zampini } else { 1072e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 1073e176bc59SStefano Zampini is->cctx = is->rctx; 1074e176bc59SStefano Zampini } 10753fd1c9e7SStefano Zampini 10763fd1c9e7SStefano Zampini /* interface counter vector (local) */ 10773fd1c9e7SStefano Zampini ierr = VecDuplicate(is->y,&is->counter);CHKERRQ(ierr); 10783fd1c9e7SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 10793fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10803fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10813fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10823fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10833fd1c9e7SStefano Zampini 10843fd1c9e7SStefano Zampini /* free workspace */ 1085e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 1086e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 10876bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 10886bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1089f26d0771SStefano Zampini } 10909c0b00dbSStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 1091b4319ba4SBarry Smith PetscFunctionReturn(0); 1092b4319ba4SBarry Smith } 1093b4319ba4SBarry Smith 10942e74eeadSLisandro Dalcin #undef __FUNCT__ 10952e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 1096a8116848SStefano Zampini static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 10972e74eeadSLisandro Dalcin { 10982e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 10992e74eeadSLisandro Dalcin PetscErrorCode ierr; 110097563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 110197563a80SStefano Zampini PetscInt i,zm,zn; 110297563a80SStefano Zampini #endif 1103f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 11042e74eeadSLisandro Dalcin 11052e74eeadSLisandro Dalcin PetscFunctionBegin; 11062e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 1107f26d0771SStefano 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); 110897563a80SStefano Zampini /* count negative indices */ 110997563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 111097563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 11112e74eeadSLisandro Dalcin #endif 111297563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 111397563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 111497563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 111597563a80SStefano Zampini /* count negative indices (should be the same as before) */ 111697563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 111797563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 111897563a80SStefano 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"); 111997563a80SStefano 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"); 112097563a80SStefano Zampini #endif 11212e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 11222e74eeadSLisandro Dalcin PetscFunctionReturn(0); 11232e74eeadSLisandro Dalcin } 11242e74eeadSLisandro Dalcin 1125b4319ba4SBarry Smith #undef __FUNCT__ 112697563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS" 1127a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 112897563a80SStefano Zampini { 112997563a80SStefano Zampini Mat_IS *is = (Mat_IS*)mat->data; 113097563a80SStefano Zampini PetscErrorCode ierr; 113197563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 113297563a80SStefano Zampini PetscInt i,zm,zn; 113397563a80SStefano Zampini #endif 1134f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 113597563a80SStefano Zampini 113697563a80SStefano Zampini PetscFunctionBegin; 113797563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 1138f26d0771SStefano 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); 113997563a80SStefano Zampini /* count negative indices */ 114097563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 114197563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 114297563a80SStefano Zampini #endif 114397563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 114497563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 114597563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 114697563a80SStefano Zampini /* count negative indices (should be the same as before) */ 114797563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 114897563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 114997563a80SStefano 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"); 115097563a80SStefano 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"); 115197563a80SStefano Zampini #endif 115297563a80SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 115397563a80SStefano Zampini PetscFunctionReturn(0); 115497563a80SStefano Zampini } 115597563a80SStefano Zampini 115697563a80SStefano Zampini #undef __FUNCT__ 1157b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 1158a8116848SStefano Zampini static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1159b4319ba4SBarry Smith { 1160dfbe8321SBarry Smith PetscErrorCode ierr; 1161b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1162b4319ba4SBarry Smith 1163b4319ba4SBarry Smith PetscFunctionBegin; 1164b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1165b4319ba4SBarry Smith PetscFunctionReturn(0); 1166b4319ba4SBarry Smith } 1167b4319ba4SBarry Smith 1168b4319ba4SBarry Smith #undef __FUNCT__ 1169f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 1170a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1171f0006bf2SLisandro Dalcin { 1172f0006bf2SLisandro Dalcin PetscErrorCode ierr; 1173f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 1174f0006bf2SLisandro Dalcin 1175f0006bf2SLisandro Dalcin PetscFunctionBegin; 1176f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1177f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 1178f0006bf2SLisandro Dalcin } 1179f0006bf2SLisandro Dalcin 1180f0006bf2SLisandro Dalcin #undef __FUNCT__ 1181f0ae7da4SStefano Zampini #define __FUNCT__ "MatISZeroRowsColumnsLocal_Private" 1182f0ae7da4SStefano Zampini static PetscErrorCode MatISZeroRowsColumnsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,PetscBool columns) 11832e74eeadSLisandro Dalcin { 1184f0ae7da4SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 11852e74eeadSLisandro Dalcin PetscErrorCode ierr; 11862e74eeadSLisandro Dalcin 11872e74eeadSLisandro Dalcin PetscFunctionBegin; 1188f0ae7da4SStefano Zampini if (!n) { 1189f0ae7da4SStefano Zampini is->pure_neumann = PETSC_TRUE; 1190f0ae7da4SStefano Zampini } else { 1191f0ae7da4SStefano Zampini PetscInt i; 1192f0ae7da4SStefano Zampini is->pure_neumann = PETSC_FALSE; 1193f0ae7da4SStefano Zampini 1194f0ae7da4SStefano Zampini if (columns) { 1195f0ae7da4SStefano Zampini ierr = MatZeroRowsColumns(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 1196f0ae7da4SStefano Zampini } else { 1197f0ae7da4SStefano Zampini ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 11982e74eeadSLisandro Dalcin } 1199f0ae7da4SStefano Zampini if (diag != 0.) { 1200f0ae7da4SStefano Zampini const PetscScalar *array; 1201f0ae7da4SStefano Zampini ierr = VecGetArrayRead(is->counter,&array);CHKERRQ(ierr); 1202f0ae7da4SStefano Zampini for (i=0; i<n; i++) { 1203f0ae7da4SStefano Zampini ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 1204f0ae7da4SStefano Zampini } 1205f0ae7da4SStefano Zampini ierr = VecRestoreArrayRead(is->counter,&array);CHKERRQ(ierr); 1206f0ae7da4SStefano Zampini } 1207f0ae7da4SStefano Zampini ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1208f0ae7da4SStefano Zampini ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1209f0ae7da4SStefano Zampini } 12102e74eeadSLisandro Dalcin PetscFunctionReturn(0); 12112e74eeadSLisandro Dalcin } 12122e74eeadSLisandro Dalcin 12132e74eeadSLisandro Dalcin #undef __FUNCT__ 1214f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_Private_IS" 1215f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Private_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b,PetscBool columns) 1216b4319ba4SBarry Smith { 12176e520ac8SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 12186e520ac8SStefano Zampini PetscInt nr,nl,len,i; 12196e520ac8SStefano Zampini PetscInt *lrows; 1220dfbe8321SBarry Smith PetscErrorCode ierr; 1221b4319ba4SBarry Smith 1222b4319ba4SBarry Smith PetscFunctionBegin; 1223f0ae7da4SStefano Zampini #if defined(PETSC_USE_DEBUG) 1224f0ae7da4SStefano Zampini if (columns || diag != 0. || (x && b)) { 1225f0ae7da4SStefano Zampini PetscBool cong; 1226f0ae7da4SStefano Zampini ierr = PetscLayoutCompare(A->rmap,A->cmap,&cong);CHKERRQ(ierr); 1227f0ae7da4SStefano 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"); 1228f0ae7da4SStefano 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"); 1229f0ae7da4SStefano Zampini if (!cong && x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"A->rmap and A->cmap need to be congruent"); 1230b4319ba4SBarry Smith } 1231f0ae7da4SStefano Zampini #endif 12326e520ac8SStefano Zampini /* get locally owned rows */ 1233f0ae7da4SStefano Zampini ierr = PetscLayoutMapLocal_Private(A->rmap,n,rows,&len,&lrows,NULL);CHKERRQ(ierr); 12346e520ac8SStefano Zampini /* fix right hand side if needed */ 12356e520ac8SStefano Zampini if (x && b) { 12366e520ac8SStefano Zampini const PetscScalar *xx; 12376e520ac8SStefano Zampini PetscScalar *bb; 12382205254eSKarl Rupp 12396e520ac8SStefano Zampini ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr); 12406e520ac8SStefano Zampini ierr = VecGetArray(b, &bb);CHKERRQ(ierr); 12416e520ac8SStefano Zampini for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]]; 12426e520ac8SStefano Zampini ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr); 12436e520ac8SStefano Zampini ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr); 1244b4319ba4SBarry Smith } 12456e520ac8SStefano Zampini /* get rows associated to the local matrices */ 1246*3d996552SStefano Zampini ierr = MatISSetUpSF(A);CHKERRQ(ierr); 12476e520ac8SStefano Zampini ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr); 12486e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr); 12496e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 12506e520ac8SStefano Zampini for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1; 12516e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 12526e520ac8SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 12536e520ac8SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 12546e520ac8SStefano Zampini ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr); 12556e520ac8SStefano Zampini for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i; 1256f0ae7da4SStefano Zampini ierr = MatISZeroRowsColumnsLocal_Private(A,nr,lrows,diag,columns);CHKERRQ(ierr); 12576e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 1258b4319ba4SBarry Smith PetscFunctionReturn(0); 1259b4319ba4SBarry Smith } 1260b4319ba4SBarry Smith 1261b4319ba4SBarry Smith #undef __FUNCT__ 1262f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRows_IS" 1263f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 1264b4319ba4SBarry Smith { 1265b4319ba4SBarry Smith PetscErrorCode ierr; 1266b4319ba4SBarry Smith 1267b4319ba4SBarry Smith PetscFunctionBegin; 1268f0ae7da4SStefano Zampini ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_FALSE);CHKERRQ(ierr); 1269f0ae7da4SStefano Zampini PetscFunctionReturn(0); 1270f0ae7da4SStefano Zampini } 1271b4319ba4SBarry Smith 1272f0ae7da4SStefano Zampini #undef __FUNCT__ 1273f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_IS" 1274f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 1275f0ae7da4SStefano Zampini { 1276f0ae7da4SStefano Zampini PetscErrorCode ierr; 1277f0ae7da4SStefano Zampini 1278f0ae7da4SStefano Zampini PetscFunctionBegin; 1279f0ae7da4SStefano Zampini ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_TRUE);CHKERRQ(ierr); 1280b4319ba4SBarry Smith PetscFunctionReturn(0); 1281b4319ba4SBarry Smith } 1282b4319ba4SBarry Smith 1283b4319ba4SBarry Smith #undef __FUNCT__ 1284b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 1285a8116848SStefano Zampini static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 1286b4319ba4SBarry Smith { 1287b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1288dfbe8321SBarry Smith PetscErrorCode ierr; 1289b4319ba4SBarry Smith 1290b4319ba4SBarry Smith PetscFunctionBegin; 1291b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 1292b4319ba4SBarry Smith PetscFunctionReturn(0); 1293b4319ba4SBarry Smith } 1294b4319ba4SBarry Smith 1295b4319ba4SBarry Smith #undef __FUNCT__ 1296b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 1297a8116848SStefano Zampini static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 1298b4319ba4SBarry Smith { 1299b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1300dfbe8321SBarry Smith PetscErrorCode ierr; 1301b4319ba4SBarry Smith 1302b4319ba4SBarry Smith PetscFunctionBegin; 1303b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 1304b4319ba4SBarry Smith PetscFunctionReturn(0); 1305b4319ba4SBarry Smith } 1306b4319ba4SBarry Smith 1307b4319ba4SBarry Smith #undef __FUNCT__ 1308b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 1309a8116848SStefano Zampini static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 1310b4319ba4SBarry Smith { 1311b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 1312b4319ba4SBarry Smith 1313b4319ba4SBarry Smith PetscFunctionBegin; 1314b4319ba4SBarry Smith *local = is->A; 1315b4319ba4SBarry Smith PetscFunctionReturn(0); 1316b4319ba4SBarry Smith } 1317b4319ba4SBarry Smith 1318b4319ba4SBarry Smith #undef __FUNCT__ 1319b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 1320b4319ba4SBarry Smith /*@ 1321b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 1322b4319ba4SBarry Smith 1323b4319ba4SBarry Smith Input Parameter: 1324b4319ba4SBarry Smith . mat - the matrix 1325b4319ba4SBarry Smith 1326b4319ba4SBarry Smith Output Parameter: 1327eb82efa4SStefano Zampini . local - the local matrix 1328b4319ba4SBarry Smith 1329b4319ba4SBarry Smith Level: advanced 1330b4319ba4SBarry Smith 1331b4319ba4SBarry Smith Notes: 1332b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 1333b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 1334b4319ba4SBarry Smith of the MatSetValues() operation. 1335b4319ba4SBarry Smith 1336b4319ba4SBarry Smith .seealso: MATIS 1337b4319ba4SBarry Smith @*/ 13387087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 1339b4319ba4SBarry Smith { 13404ac538c5SBarry Smith PetscErrorCode ierr; 1341b4319ba4SBarry Smith 1342b4319ba4SBarry Smith PetscFunctionBegin; 13430700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1344b4319ba4SBarry Smith PetscValidPointer(local,2); 13454ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 1346b4319ba4SBarry Smith PetscFunctionReturn(0); 1347b4319ba4SBarry Smith } 1348b4319ba4SBarry Smith 13493b03a366Sstefano_zampini #undef __FUNCT__ 13503b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 1351a8116848SStefano Zampini static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 13523b03a366Sstefano_zampini { 13533b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 13543b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 13553b03a366Sstefano_zampini PetscErrorCode ierr; 13563b03a366Sstefano_zampini 13573b03a366Sstefano_zampini PetscFunctionBegin; 13584e4c7dbeSStefano Zampini if (is->A) { 13593b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 13603b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 1361f0ae7da4SStefano 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); 13624e4c7dbeSStefano Zampini } 13633b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 13643b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 13653b03a366Sstefano_zampini is->A = local; 13663b03a366Sstefano_zampini PetscFunctionReturn(0); 13673b03a366Sstefano_zampini } 13683b03a366Sstefano_zampini 13693b03a366Sstefano_zampini #undef __FUNCT__ 13703b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 13713b03a366Sstefano_zampini /*@ 1372eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 13733b03a366Sstefano_zampini 13743b03a366Sstefano_zampini Input Parameter: 13753b03a366Sstefano_zampini . mat - the matrix 1376eb82efa4SStefano Zampini . local - the local matrix 13773b03a366Sstefano_zampini 13783b03a366Sstefano_zampini Output Parameter: 13793b03a366Sstefano_zampini 13803b03a366Sstefano_zampini Level: advanced 13813b03a366Sstefano_zampini 13823b03a366Sstefano_zampini Notes: 13833b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 13843b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 13853b03a366Sstefano_zampini 13863b03a366Sstefano_zampini .seealso: MATIS 13873b03a366Sstefano_zampini @*/ 13883b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 13893b03a366Sstefano_zampini { 13903b03a366Sstefano_zampini PetscErrorCode ierr; 13913b03a366Sstefano_zampini 13923b03a366Sstefano_zampini PetscFunctionBegin; 13933b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1394b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 13953b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 13963b03a366Sstefano_zampini PetscFunctionReturn(0); 13973b03a366Sstefano_zampini } 13983b03a366Sstefano_zampini 13996726f965SBarry Smith #undef __FUNCT__ 14006726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 1401a8116848SStefano Zampini static PetscErrorCode MatZeroEntries_IS(Mat A) 14026726f965SBarry Smith { 14036726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 14046726f965SBarry Smith PetscErrorCode ierr; 14056726f965SBarry Smith 14066726f965SBarry Smith PetscFunctionBegin; 14076726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 14086726f965SBarry Smith PetscFunctionReturn(0); 14096726f965SBarry Smith } 14106726f965SBarry Smith 14116726f965SBarry Smith #undef __FUNCT__ 14122e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 1413a8116848SStefano Zampini static PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 14142e74eeadSLisandro Dalcin { 14152e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 14162e74eeadSLisandro Dalcin PetscErrorCode ierr; 14172e74eeadSLisandro Dalcin 14182e74eeadSLisandro Dalcin PetscFunctionBegin; 14192e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 14202e74eeadSLisandro Dalcin PetscFunctionReturn(0); 14212e74eeadSLisandro Dalcin } 14222e74eeadSLisandro Dalcin 14232e74eeadSLisandro Dalcin #undef __FUNCT__ 14242e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 1425a8116848SStefano Zampini static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 14262e74eeadSLisandro Dalcin { 14272e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 14282e74eeadSLisandro Dalcin PetscErrorCode ierr; 14292e74eeadSLisandro Dalcin 14302e74eeadSLisandro Dalcin PetscFunctionBegin; 14312e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 1432e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 14332e74eeadSLisandro Dalcin 14342e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 14352e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 1436e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1437e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 14382e74eeadSLisandro Dalcin PetscFunctionReturn(0); 14392e74eeadSLisandro Dalcin } 14402e74eeadSLisandro Dalcin 14412e74eeadSLisandro Dalcin #undef __FUNCT__ 14426726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 1443a8116848SStefano Zampini static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 14446726f965SBarry Smith { 14456726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 14466726f965SBarry Smith PetscErrorCode ierr; 14476726f965SBarry Smith 14486726f965SBarry Smith PetscFunctionBegin; 14494e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 14506726f965SBarry Smith PetscFunctionReturn(0); 14516726f965SBarry Smith } 14526726f965SBarry Smith 1453284134d9SBarry Smith #undef __FUNCT__ 1454f26d0771SStefano Zampini #define __FUNCT__ "MatAXPY_IS" 1455f26d0771SStefano Zampini static PetscErrorCode MatAXPY_IS(Mat Y,PetscScalar a,Mat X,MatStructure str) 1456f26d0771SStefano Zampini { 1457f26d0771SStefano Zampini Mat_IS *y = (Mat_IS*)Y->data; 1458f26d0771SStefano Zampini Mat_IS *x; 1459f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1460f26d0771SStefano Zampini PetscBool ismatis; 1461f26d0771SStefano Zampini #endif 1462f26d0771SStefano Zampini PetscErrorCode ierr; 1463f26d0771SStefano Zampini 1464f26d0771SStefano Zampini PetscFunctionBegin; 1465f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1466f26d0771SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)X,MATIS,&ismatis);CHKERRQ(ierr); 1467f26d0771SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS"); 1468f26d0771SStefano Zampini #endif 1469f26d0771SStefano Zampini x = (Mat_IS*)X->data; 1470f26d0771SStefano Zampini ierr = MatAXPY(y->A,a,x->A,str);CHKERRQ(ierr); 1471f26d0771SStefano Zampini PetscFunctionReturn(0); 1472f26d0771SStefano Zampini } 1473f26d0771SStefano Zampini 1474f26d0771SStefano Zampini #undef __FUNCT__ 1475f26d0771SStefano Zampini #define __FUNCT__ "MatGetLocalSubMatrix_IS" 1476f26d0771SStefano Zampini static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat) 1477f26d0771SStefano Zampini { 1478f26d0771SStefano Zampini Mat lA; 1479f26d0771SStefano Zampini Mat_IS *matis; 1480f26d0771SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 1481f26d0771SStefano Zampini IS is; 1482f26d0771SStefano Zampini const PetscInt *rg,*rl; 1483f26d0771SStefano Zampini PetscInt nrg; 1484f26d0771SStefano Zampini PetscInt N,M,nrl,i,*idxs; 1485f26d0771SStefano Zampini PetscErrorCode ierr; 1486f26d0771SStefano Zampini 1487f26d0771SStefano Zampini PetscFunctionBegin; 1488f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 1489f26d0771SStefano Zampini ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr); 1490f26d0771SStefano Zampini ierr = ISGetIndices(row,&rl);CHKERRQ(ierr); 1491f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr); 1492f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1493f0ae7da4SStefano 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); 1494f26d0771SStefano Zampini #endif 1495f26d0771SStefano Zampini ierr = PetscMalloc1(nrg,&idxs);CHKERRQ(ierr); 1496f26d0771SStefano Zampini /* map from [0,nrl) to row */ 1497f26d0771SStefano Zampini for (i=0;i<nrl;i++) idxs[i] = rl[i]; 1498f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1499f26d0771SStefano Zampini for (i=nrl;i<nrg;i++) idxs[i] = nrg; 1500f26d0771SStefano Zampini #else 1501f26d0771SStefano Zampini for (i=nrl;i<nrg;i++) idxs[i] = -1; 1502f26d0771SStefano Zampini #endif 1503f26d0771SStefano Zampini ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr); 1504f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 1505f26d0771SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 1506f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 1507f26d0771SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 1508f26d0771SStefano Zampini /* compute new l2g map for columns */ 1509f26d0771SStefano Zampini if (col != row || A->rmap->mapping != A->cmap->mapping) { 1510f26d0771SStefano Zampini const PetscInt *cg,*cl; 1511f26d0771SStefano Zampini PetscInt ncg; 1512f26d0771SStefano Zampini PetscInt ncl; 1513f26d0771SStefano Zampini 1514f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 1515f26d0771SStefano Zampini ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr); 1516f26d0771SStefano Zampini ierr = ISGetIndices(col,&cl);CHKERRQ(ierr); 1517f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr); 1518f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1519f0ae7da4SStefano 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); 1520f26d0771SStefano Zampini #endif 1521f26d0771SStefano Zampini ierr = PetscMalloc1(ncg,&idxs);CHKERRQ(ierr); 1522f26d0771SStefano Zampini /* map from [0,ncl) to col */ 1523f26d0771SStefano Zampini for (i=0;i<ncl;i++) idxs[i] = cl[i]; 1524f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1525f26d0771SStefano Zampini for (i=ncl;i<ncg;i++) idxs[i] = ncg; 1526f26d0771SStefano Zampini #else 1527f26d0771SStefano Zampini for (i=ncl;i<ncg;i++) idxs[i] = -1; 1528f26d0771SStefano Zampini #endif 1529f26d0771SStefano Zampini ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr); 1530f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 1531f26d0771SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 1532f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 1533f26d0771SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 1534f26d0771SStefano Zampini } else { 1535f26d0771SStefano Zampini ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr); 1536f26d0771SStefano Zampini cl2g = rl2g; 1537f26d0771SStefano Zampini } 1538f26d0771SStefano Zampini /* create the MATIS submatrix */ 1539f26d0771SStefano Zampini ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 1540f26d0771SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr); 1541f26d0771SStefano Zampini ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1542f26d0771SStefano Zampini ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr); 1543b0aa3428SStefano Zampini matis = (Mat_IS*)((*submat)->data); 1544f26d0771SStefano Zampini matis->islocalref = PETSC_TRUE; 1545f26d0771SStefano Zampini ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr); 1546f26d0771SStefano Zampini ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr); 1547f26d0771SStefano Zampini ierr = MatISSetLocalMat(*submat,lA);CHKERRQ(ierr); 1548f26d0771SStefano Zampini ierr = MatSetUp(*submat);CHKERRQ(ierr); 1549f26d0771SStefano Zampini ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1550f26d0771SStefano Zampini ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1551f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 1552f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 1553f26d0771SStefano Zampini /* remove unsupported ops */ 1554f26d0771SStefano Zampini ierr = PetscMemzero((*submat)->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1555f26d0771SStefano Zampini (*submat)->ops->destroy = MatDestroy_IS; 1556f26d0771SStefano Zampini (*submat)->ops->setvalueslocal = MatSetValuesLocal_SubMat_IS; 1557f26d0771SStefano Zampini (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS; 1558f26d0771SStefano Zampini (*submat)->ops->assemblybegin = MatAssemblyBegin_IS; 1559f26d0771SStefano Zampini (*submat)->ops->assemblyend = MatAssemblyEnd_IS; 1560f26d0771SStefano Zampini PetscFunctionReturn(0); 1561f26d0771SStefano Zampini } 1562f26d0771SStefano Zampini 1563f26d0771SStefano Zampini #undef __FUNCT__ 1564284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 1565284134d9SBarry Smith /*@ 15663c212e90SHong Zhang MatCreateIS - Creates a "process" unassembled matrix, assembled on each 1567284134d9SBarry Smith process but not across processes. 1568284134d9SBarry Smith 1569284134d9SBarry Smith Input Parameters: 1570284134d9SBarry Smith + comm - MPI communicator that will share the matrix 1571e176bc59SStefano Zampini . bs - block size of the matrix 1572df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 1573e176bc59SStefano Zampini . rmap - local to global map for rows 1574e176bc59SStefano Zampini - cmap - local to global map for cols 1575284134d9SBarry Smith 1576284134d9SBarry Smith Output Parameter: 1577284134d9SBarry Smith . A - the resulting matrix 1578284134d9SBarry Smith 15798e6c10adSSatish Balay Level: advanced 15808e6c10adSSatish Balay 15813c212e90SHong Zhang Notes: See MATIS for more details. 15823c212e90SHong Zhang m and n are NOT related to the size of the map; they are the size of the part of the vector owned 1583e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 15843c212e90SHong Zhang If either rmap or cmap are NULL, then the matrix is assumed to be square. 1585284134d9SBarry Smith 1586284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 1587284134d9SBarry Smith @*/ 1588e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 1589284134d9SBarry Smith { 1590284134d9SBarry Smith PetscErrorCode ierr; 1591284134d9SBarry Smith 1592284134d9SBarry Smith PetscFunctionBegin; 1593e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 1594284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1595d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 1596284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1597284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 1598e176bc59SStefano Zampini if (rmap && cmap) { 1599e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 1600e176bc59SStefano Zampini } else if (!rmap) { 1601e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 1602e176bc59SStefano Zampini } else { 1603e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 1604e176bc59SStefano Zampini } 1605284134d9SBarry Smith PetscFunctionReturn(0); 1606284134d9SBarry Smith } 1607284134d9SBarry Smith 1608b4319ba4SBarry Smith /*MC 1609f26d0771SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition methods (e.g. PCBDDC or KSPFETIDP). 1610b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 1611b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 1612b4319ba4SBarry Smith product is handled "implicitly". 1613b4319ba4SBarry Smith 1614b4319ba4SBarry Smith Operations Provided: 16156726f965SBarry Smith + MatMult() 16162e74eeadSLisandro Dalcin . MatMultAdd() 16172e74eeadSLisandro Dalcin . MatMultTranspose() 16182e74eeadSLisandro Dalcin . MatMultTransposeAdd() 16196726f965SBarry Smith . MatZeroEntries() 16206726f965SBarry Smith . MatSetOption() 16212e74eeadSLisandro Dalcin . MatZeroRows() 16222e74eeadSLisandro Dalcin . MatSetValues() 162397563a80SStefano Zampini . MatSetValuesBlocked() 16246726f965SBarry Smith . MatSetValuesLocal() 162597563a80SStefano Zampini . MatSetValuesBlockedLocal() 16262e74eeadSLisandro Dalcin . MatScale() 16272e74eeadSLisandro Dalcin . MatGetDiagonal() 16282b404112SStefano Zampini . MatMissingDiagonal() 16292b404112SStefano Zampini . MatDuplicate() 16302b404112SStefano Zampini . MatCopy() 1631f26d0771SStefano Zampini . MatAXPY() 1632f26d0771SStefano Zampini . MatGetSubMatrix() 1633f26d0771SStefano Zampini . MatGetLocalSubMatrix() 16346726f965SBarry Smith - MatSetLocalToGlobalMapping() 1635b4319ba4SBarry Smith 1636b4319ba4SBarry Smith Options Database Keys: 1637b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1638b4319ba4SBarry Smith 1639b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1640b4319ba4SBarry Smith 1641b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1642b4319ba4SBarry Smith 1643b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1644eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1645b4319ba4SBarry Smith 1646b4319ba4SBarry Smith Level: advanced 1647b4319ba4SBarry Smith 1648f26d0771SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), MatCreateIS(), PCBDDC, KSPFETIDP 1649b4319ba4SBarry Smith 1650b4319ba4SBarry Smith M*/ 1651b4319ba4SBarry Smith 1652b4319ba4SBarry Smith #undef __FUNCT__ 1653b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 16548cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1655b4319ba4SBarry Smith { 1656dfbe8321SBarry Smith PetscErrorCode ierr; 1657b4319ba4SBarry Smith Mat_IS *b; 1658b4319ba4SBarry Smith 1659b4319ba4SBarry Smith PetscFunctionBegin; 1660b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1661b4319ba4SBarry Smith A->data = (void*)b; 1662b4319ba4SBarry Smith 1663e176bc59SStefano Zampini /* matrix ops */ 1664e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1665b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 16662e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 16672e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 16682e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1669b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1670b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 16712e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 167298921651SStefano Zampini A->ops->setvaluesblocked = MatSetValuesBlocked_IS; 1673b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1674f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 16752e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1676f0ae7da4SStefano Zampini A->ops->zerorowscolumns = MatZeroRowsColumns_IS; 1677b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1678b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1679b4319ba4SBarry Smith A->ops->view = MatView_IS; 16806726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 16812e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 16822e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 16836726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 168469796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 168569796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1686ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 16876bd84002SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_IS; 16882b404112SStefano Zampini A->ops->copy = MatCopy_IS; 1689659959c5SStefano Zampini A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_IS; 1690a8116848SStefano Zampini A->ops->getsubmatrix = MatGetSubMatrix_IS; 1691f26d0771SStefano Zampini A->ops->axpy = MatAXPY_IS; 16923fd1c9e7SStefano Zampini A->ops->diagonalset = MatDiagonalSet_IS; 16933fd1c9e7SStefano Zampini A->ops->shift = MatShift_IS; 1694b4319ba4SBarry Smith 1695b7ce53b6SStefano Zampini /* special MATIS functions */ 1696bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1697bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1698b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 16992e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 1700cf0a3239SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",MatISSetUpSF_IS);CHKERRQ(ierr); 170117667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1702b4319ba4SBarry Smith PetscFunctionReturn(0); 1703b4319ba4SBarry Smith } 1704