1b4319ba4SBarry Smith /* 2b4319ba4SBarry Smith Creates a matrix class for using the Neumann-Neumann type preconditioners. 3b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 4b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 5b4319ba4SBarry Smith product is handled "implicitly". 6b4319ba4SBarry Smith 7b4319ba4SBarry Smith Currently this allows for only one subdomain per processor. 8b4319ba4SBarry Smith */ 9b4319ba4SBarry Smith 10c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/ 1128f4e0baSStefano Zampini 12f26d0771SStefano Zampini #define MATIS_MAX_ENTRIES_INSERTION 2048 13f26d0771SStefano Zampini 14a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat); 157230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat,PetscInt,const PetscInt[],PetscScalar); 16a72627d2SStefano Zampini 1728f4e0baSStefano Zampini #undef __FUNCT__ 18*3fd1c9e7SStefano Zampini #define __FUNCT__ "MatDiagonalSet_IS" 19*3fd1c9e7SStefano Zampini PetscErrorCode MatDiagonalSet_IS(Mat A,Vec D,InsertMode insmode) 20*3fd1c9e7SStefano Zampini { 21*3fd1c9e7SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 22*3fd1c9e7SStefano Zampini PetscErrorCode ierr; 23*3fd1c9e7SStefano Zampini 24*3fd1c9e7SStefano Zampini PetscFunctionBegin; 25*3fd1c9e7SStefano Zampini if (!D) { /* this code branch is used by MatShift_IS */ 26*3fd1c9e7SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 27*3fd1c9e7SStefano Zampini } else { 28*3fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 29*3fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 30*3fd1c9e7SStefano Zampini } 31*3fd1c9e7SStefano Zampini ierr = VecPointwiseDivide(is->y,is->y,is->counter);CHKERRQ(ierr); 32*3fd1c9e7SStefano Zampini ierr = MatDiagonalSet(is->A,is->y,insmode);CHKERRQ(ierr); 33*3fd1c9e7SStefano Zampini PetscFunctionReturn(0); 34*3fd1c9e7SStefano Zampini } 35*3fd1c9e7SStefano Zampini 36*3fd1c9e7SStefano Zampini #undef __FUNCT__ 37*3fd1c9e7SStefano Zampini #define __FUNCT__ "MatShift_IS" 38*3fd1c9e7SStefano Zampini PetscErrorCode MatShift_IS(Mat A,PetscScalar a) 39*3fd1c9e7SStefano Zampini { 40*3fd1c9e7SStefano Zampini PetscErrorCode ierr; 41*3fd1c9e7SStefano Zampini 42*3fd1c9e7SStefano Zampini PetscFunctionBegin; 43*3fd1c9e7SStefano Zampini ierr = MatDiagonalSet_IS(A,NULL,ADD_VALUES);CHKERRQ(ierr); 44*3fd1c9e7SStefano Zampini PetscFunctionReturn(0); 45*3fd1c9e7SStefano Zampini } 46*3fd1c9e7SStefano Zampini 47*3fd1c9e7SStefano Zampini #undef __FUNCT__ 48f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesLocal_SubMat_IS" 49f26d0771SStefano Zampini static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 50f26d0771SStefano Zampini { 51f26d0771SStefano Zampini PetscErrorCode ierr; 52f26d0771SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 53f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 54f26d0771SStefano Zampini 55f26d0771SStefano Zampini PetscFunctionBegin; 56f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 57f26d0771SStefano 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); 58f26d0771SStefano Zampini #endif 59f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 60f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 61f26d0771SStefano Zampini ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 62f26d0771SStefano Zampini PetscFunctionReturn(0); 63f26d0771SStefano Zampini } 64f26d0771SStefano Zampini 65f26d0771SStefano Zampini #undef __FUNCT__ 66f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesBlockedLocal_SubMat_IS" 67f26d0771SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 68f26d0771SStefano Zampini { 69f26d0771SStefano Zampini PetscErrorCode ierr; 70f26d0771SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 71f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 72f26d0771SStefano Zampini 73f26d0771SStefano Zampini PetscFunctionBegin; 74f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 75f26d0771SStefano 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); 76f26d0771SStefano Zampini #endif 77f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 78f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 79f26d0771SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 80f26d0771SStefano Zampini PetscFunctionReturn(0); 81f26d0771SStefano Zampini } 82f26d0771SStefano Zampini 83f26d0771SStefano Zampini #undef __FUNCT__ 84a8116848SStefano Zampini #define __FUNCT__ "PetscLayoutMapLocal_Private" 85a8116848SStefano Zampini static PetscErrorCode PetscLayoutMapLocal_Private(PetscLayout map,PetscInt N,const PetscInt idxs[],const PetscInt gidxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs) 86a8116848SStefano Zampini { 87a8116848SStefano Zampini PetscInt *owners = map->range; 88a8116848SStefano Zampini PetscInt n = map->n; 89a8116848SStefano Zampini PetscSF sf; 90a8116848SStefano Zampini PetscInt *lidxs,*work = NULL; 91a8116848SStefano Zampini PetscSFNode *ridxs; 92a8116848SStefano Zampini PetscMPIInt rank; 93a8116848SStefano Zampini PetscInt r, p = 0, len = 0; 94a8116848SStefano Zampini PetscErrorCode ierr; 95a8116848SStefano Zampini 96a8116848SStefano Zampini PetscFunctionBegin; 97a8116848SStefano Zampini /* Create SF where leaves are input idxs and roots are owned idxs (code adapted from MatZeroRowsMapLocal_Private) */ 98a8116848SStefano Zampini ierr = MPI_Comm_rank(map->comm,&rank);CHKERRQ(ierr); 99a8116848SStefano Zampini ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr); 100a8116848SStefano Zampini for (r = 0; r < n; ++r) lidxs[r] = -1; 101a8116848SStefano Zampini ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr); 102a8116848SStefano Zampini for (r = 0; r < N; ++r) { 103a8116848SStefano Zampini const PetscInt idx = idxs[r]; 104a8116848SStefano 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); 105a8116848SStefano Zampini if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 106a8116848SStefano Zampini ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr); 107a8116848SStefano Zampini } 108a8116848SStefano Zampini ridxs[r].rank = p; 109a8116848SStefano Zampini ridxs[r].index = idxs[r] - owners[p]; 110a8116848SStefano Zampini } 111a8116848SStefano Zampini ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr); 112a8116848SStefano Zampini ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr); 113a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 114a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 115a8116848SStefano Zampini if ((gidxs && ogidxs) || ogidxs) { 116a8116848SStefano Zampini ierr = PetscMalloc1(n,&work);CHKERRQ(ierr); 117a8116848SStefano Zampini if (!gidxs) { /* if gidxs is not provided, global indices are inferred */ 118a8116848SStefano Zampini PetscInt cum = 0,start,*work2; 119a8116848SStefano Zampini ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr); 120a8116848SStefano Zampini for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++; 121a8116848SStefano Zampini ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRQ(ierr); 122a8116848SStefano Zampini start -= cum; 123a8116848SStefano Zampini cum = 0; 124a8116848SStefano Zampini for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++; 125a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 126a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 127a8116848SStefano Zampini ierr = PetscFree(work2);CHKERRQ(ierr); 128a8116848SStefano Zampini } else { 129a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)gidxs,work,MPIU_REPLACE);CHKERRQ(ierr); 130a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)gidxs,work,MPIU_REPLACE);CHKERRQ(ierr); 131a8116848SStefano Zampini } 132a8116848SStefano Zampini } 133a8116848SStefano Zampini ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 134a8116848SStefano Zampini /* Compress and put in indices */ 135a8116848SStefano Zampini for (r = 0; r < n; ++r) 136a8116848SStefano Zampini if (lidxs[r] >= 0) { 137a8116848SStefano Zampini if (work) work[len] = work[r]; 138a8116848SStefano Zampini lidxs[len++] = r; 139a8116848SStefano Zampini } 140a8116848SStefano Zampini if (on) *on = len; 141a8116848SStefano Zampini if (oidxs) *oidxs = lidxs; 142a8116848SStefano Zampini if (ogidxs) *ogidxs = work; 143a8116848SStefano Zampini PetscFunctionReturn(0); 144a8116848SStefano Zampini } 145a8116848SStefano Zampini 146a8116848SStefano Zampini #undef __FUNCT__ 147a8116848SStefano Zampini #define __FUNCT__ "MatGetSubMatrix_IS" 148a8116848SStefano Zampini static PetscErrorCode MatGetSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat) 149a8116848SStefano Zampini { 150a8116848SStefano Zampini Mat locmat,newlocmat; 151a8116848SStefano Zampini Mat_IS *newmatis; 152a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG) 153a8116848SStefano Zampini Vec rtest,ltest; 154a8116848SStefano Zampini const PetscScalar *array; 155a8116848SStefano Zampini #endif 156a8116848SStefano Zampini const PetscInt *idxs; 157a8116848SStefano Zampini PetscInt i,m,n; 158a8116848SStefano Zampini PetscErrorCode ierr; 159a8116848SStefano Zampini 160a8116848SStefano Zampini PetscFunctionBegin; 161a8116848SStefano Zampini if (scall == MAT_REUSE_MATRIX) { 162a8116848SStefano Zampini PetscBool ismatis; 163a8116848SStefano Zampini 164a8116848SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr); 165a8116848SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type"); 166a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 167a8116848SStefano Zampini if (!newmatis->getsub_ris) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS"); 168a8116848SStefano Zampini if (!newmatis->getsub_cis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS"); 169a8116848SStefano Zampini } 170a8116848SStefano Zampini /* irow and icol may not have duplicate entries */ 171a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG) 172a8116848SStefano Zampini ierr = MatCreateVecs(mat,<est,&rtest);CHKERRQ(ierr); 173a8116848SStefano Zampini ierr = ISGetLocalSize(irow,&n);CHKERRQ(ierr); 174a8116848SStefano Zampini ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 175a8116848SStefano Zampini for (i=0;i<n;i++) { 176a8116848SStefano Zampini ierr = VecSetValue(rtest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 177a8116848SStefano Zampini } 178a8116848SStefano Zampini ierr = VecAssemblyBegin(rtest);CHKERRQ(ierr); 179a8116848SStefano Zampini ierr = VecAssemblyEnd(rtest);CHKERRQ(ierr); 180a8116848SStefano Zampini ierr = VecGetLocalSize(rtest,&n);CHKERRQ(ierr); 181a8116848SStefano Zampini ierr = VecGetOwnershipRange(rtest,&m,NULL);CHKERRQ(ierr); 182a8116848SStefano Zampini ierr = VecGetArrayRead(rtest,&array);CHKERRQ(ierr); 183fd479f66SMatthew 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])); 184a8116848SStefano Zampini ierr = VecRestoreArrayRead(rtest,&array);CHKERRQ(ierr); 185a8116848SStefano Zampini ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 186a8116848SStefano Zampini ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 187a8116848SStefano Zampini ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 188a8116848SStefano Zampini for (i=0;i<n;i++) { 189a8116848SStefano Zampini ierr = VecSetValue(ltest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 190a8116848SStefano Zampini } 191a8116848SStefano Zampini ierr = VecAssemblyBegin(ltest);CHKERRQ(ierr); 192a8116848SStefano Zampini ierr = VecAssemblyEnd(ltest);CHKERRQ(ierr); 193a8116848SStefano Zampini ierr = VecGetLocalSize(ltest,&n);CHKERRQ(ierr); 194a8116848SStefano Zampini ierr = VecGetOwnershipRange(ltest,&m,NULL);CHKERRQ(ierr); 195a8116848SStefano Zampini ierr = VecGetArrayRead(ltest,&array);CHKERRQ(ierr); 196fd479f66SMatthew 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])); 197a8116848SStefano Zampini ierr = VecRestoreArrayRead(ltest,&array);CHKERRQ(ierr); 198a8116848SStefano Zampini ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 199a8116848SStefano Zampini ierr = VecDestroy(&rtest);CHKERRQ(ierr); 200a8116848SStefano Zampini ierr = VecDestroy(<est);CHKERRQ(ierr); 201a8116848SStefano Zampini #endif 202a8116848SStefano Zampini if (scall == MAT_INITIAL_MATRIX) { 203a8116848SStefano Zampini Mat_IS *matis = (Mat_IS*)mat->data; 204a8116848SStefano Zampini ISLocalToGlobalMapping rl2g; 205a8116848SStefano Zampini IS is; 206a8116848SStefano Zampini PetscInt *lidxs,*lgidxs,*newgidxs; 2076dd40735SStefano Zampini PetscInt ll,newloc; 208a8116848SStefano Zampini MPI_Comm comm; 209a8116848SStefano Zampini 210a8116848SStefano Zampini ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 211a8116848SStefano Zampini ierr = ISGetLocalSize(irow,&m);CHKERRQ(ierr); 212a8116848SStefano Zampini ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 213a8116848SStefano Zampini ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 214a8116848SStefano Zampini ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr); 215a8116848SStefano Zampini ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 216a8116848SStefano Zampini /* communicate irow to their owners in the layout */ 217a8116848SStefano Zampini ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 218a8116848SStefano Zampini ierr = PetscLayoutMapLocal_Private(mat->rmap,m,idxs,NULL,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 219a8116848SStefano Zampini ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 220a8116848SStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 221a8116848SStefano Zampini ierr = MatISComputeSF_Private(mat);CHKERRQ(ierr); 222a8116848SStefano Zampini } 223a8116848SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,matis->sf_nroots*sizeof(PetscInt));CHKERRQ(ierr); 224a8116848SStefano Zampini for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1; 225a8116848SStefano Zampini ierr = PetscFree(lidxs);CHKERRQ(ierr); 226a8116848SStefano Zampini ierr = PetscFree(lgidxs);CHKERRQ(ierr); 227a8116848SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 228a8116848SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 229a8116848SStefano Zampini for (i=0,newloc=0;i<matis->sf_nleaves;i++) if (matis->sf_leafdata[i]) newloc++; 230a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 231a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 232a8116848SStefano Zampini for (i=0,newloc=0;i<matis->sf_nleaves;i++) 233a8116848SStefano Zampini if (matis->sf_leafdata[i]) { 234a8116848SStefano Zampini lidxs[newloc] = i; 235a8116848SStefano Zampini newgidxs[newloc++] = matis->sf_leafdata[i]-1; 236a8116848SStefano Zampini } 237a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 238a8116848SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 239a8116848SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 240a8116848SStefano Zampini /* local is to extract local submatrix */ 241a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 242a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris);CHKERRQ(ierr); 243a8116848SStefano Zampini if (mat->congruentlayouts == -1) { /* first time we compare rows and cols layouts */ 244a8116848SStefano Zampini PetscBool cong; 245a8116848SStefano Zampini ierr = PetscLayoutCompare(mat->rmap,mat->cmap,&cong);CHKERRQ(ierr); 246a8116848SStefano Zampini if (cong) mat->congruentlayouts = 1; 247a8116848SStefano Zampini else mat->congruentlayouts = 0; 248a8116848SStefano Zampini } 249a8116848SStefano Zampini if (mat->congruentlayouts && irow == icol) { 250a8116848SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g);CHKERRQ(ierr); 251a8116848SStefano Zampini ierr = PetscObjectReference((PetscObject)newmatis->getsub_ris);CHKERRQ(ierr); 252a8116848SStefano Zampini newmatis->getsub_cis = newmatis->getsub_ris; 253a8116848SStefano Zampini } else { 254a8116848SStefano Zampini ISLocalToGlobalMapping cl2g; 255a8116848SStefano Zampini 256a8116848SStefano Zampini /* communicate icol to their owners in the layout */ 257a8116848SStefano Zampini ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 258a8116848SStefano Zampini ierr = PetscLayoutMapLocal_Private(mat->cmap,n,idxs,NULL,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 259a8116848SStefano Zampini ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 260a8116848SStefano Zampini ierr = PetscMemzero(matis->csf_rootdata,matis->csf_nroots*sizeof(PetscInt));CHKERRQ(ierr); 261a8116848SStefano Zampini for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1; 262a8116848SStefano Zampini ierr = PetscFree(lidxs);CHKERRQ(ierr); 263a8116848SStefano Zampini ierr = PetscFree(lgidxs);CHKERRQ(ierr); 264a8116848SStefano Zampini ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 265a8116848SStefano Zampini ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 266a8116848SStefano Zampini for (i=0,newloc=0;i<matis->csf_nleaves;i++) if (matis->csf_leafdata[i]) newloc++; 267a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 268a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 269a8116848SStefano Zampini for (i=0,newloc=0;i<matis->csf_nleaves;i++) 270a8116848SStefano Zampini if (matis->csf_leafdata[i]) { 271a8116848SStefano Zampini lidxs[newloc] = i; 272a8116848SStefano Zampini newgidxs[newloc++] = matis->csf_leafdata[i]-1; 273a8116848SStefano Zampini } 274a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 275a8116848SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 276a8116848SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 277a8116848SStefano Zampini /* local is to extract local submatrix */ 278a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis);CHKERRQ(ierr); 279a8116848SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr); 280a8116848SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 281a8116848SStefano Zampini } 282a8116848SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 283a8116848SStefano Zampini } else { 284a8116848SStefano Zampini ierr = MatISGetLocalMat(*newmat,&newlocmat);CHKERRQ(ierr); 285a8116848SStefano Zampini } 286a8116848SStefano Zampini ierr = MatISGetLocalMat(mat,&locmat);CHKERRQ(ierr); 287a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 288a8116848SStefano Zampini ierr = MatGetSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat);CHKERRQ(ierr); 289a8116848SStefano Zampini if (scall == MAT_INITIAL_MATRIX) { 290a8116848SStefano Zampini ierr = MatISSetLocalMat(*newmat,newlocmat);CHKERRQ(ierr); 291a8116848SStefano Zampini ierr = MatDestroy(&newlocmat);CHKERRQ(ierr); 292a8116848SStefano Zampini } 293a8116848SStefano Zampini ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 294a8116848SStefano Zampini ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 295a8116848SStefano Zampini PetscFunctionReturn(0); 296a8116848SStefano Zampini } 297a8116848SStefano Zampini 298a8116848SStefano Zampini #undef __FUNCT__ 2992b404112SStefano Zampini #define __FUNCT__ "MatCopy_IS" 300a8116848SStefano Zampini static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str) 3012b404112SStefano Zampini { 3022b404112SStefano Zampini Mat_IS *a = (Mat_IS*)A->data,*b; 3032b404112SStefano Zampini PetscBool ismatis; 3042b404112SStefano Zampini PetscErrorCode ierr; 3052b404112SStefano Zampini 3062b404112SStefano Zampini PetscFunctionBegin; 3072b404112SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr); 3082b404112SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented"); 3092b404112SStefano Zampini b = (Mat_IS*)B->data; 3102b404112SStefano Zampini ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr); 3112b404112SStefano Zampini PetscFunctionReturn(0); 3122b404112SStefano Zampini } 3132b404112SStefano Zampini 3142b404112SStefano Zampini #undef __FUNCT__ 3156bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS" 316a8116848SStefano Zampini static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool *missing,PetscInt *d) 3176bd84002SStefano Zampini { 318527b2640SStefano Zampini Vec v; 319527b2640SStefano Zampini const PetscScalar *array; 320527b2640SStefano Zampini PetscInt i,n; 3216bd84002SStefano Zampini PetscErrorCode ierr; 3226bd84002SStefano Zampini 3236bd84002SStefano Zampini PetscFunctionBegin; 324527b2640SStefano Zampini *missing = PETSC_FALSE; 325527b2640SStefano Zampini ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr); 326527b2640SStefano Zampini ierr = MatGetDiagonal(A,v);CHKERRQ(ierr); 327527b2640SStefano Zampini ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 328527b2640SStefano Zampini ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr); 329527b2640SStefano Zampini for (i=0;i<n;i++) if (array[i] == 0.) break; 330527b2640SStefano Zampini ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr); 331527b2640SStefano Zampini ierr = VecDestroy(&v);CHKERRQ(ierr); 332527b2640SStefano Zampini if (i != n) *missing = PETSC_TRUE; 333527b2640SStefano Zampini if (d) { 334527b2640SStefano Zampini *d = -1; 335527b2640SStefano Zampini if (*missing) { 336527b2640SStefano Zampini PetscInt rstart; 337527b2640SStefano Zampini ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 338527b2640SStefano Zampini *d = i+rstart; 339527b2640SStefano Zampini } 340527b2640SStefano Zampini } 3416bd84002SStefano Zampini PetscFunctionReturn(0); 3426bd84002SStefano Zampini } 3436bd84002SStefano Zampini 3446bd84002SStefano Zampini #undef __FUNCT__ 34528f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private" 346a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B) 34728f4e0baSStefano Zampini { 34828f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 34928f4e0baSStefano Zampini const PetscInt *gidxs; 35028f4e0baSStefano Zampini PetscErrorCode ierr; 35128f4e0baSStefano Zampini 35228f4e0baSStefano Zampini PetscFunctionBegin; 35328f4e0baSStefano Zampini ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr); 35428f4e0baSStefano Zampini ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr); 35528f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 3563bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 35728f4e0baSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 3583bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 35928f4e0baSStefano Zampini ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 360a8116848SStefano Zampini if (B->rmap->mapping != B->cmap->mapping) { /* setup SF for columns */ 361a8116848SStefano Zampini ierr = MatGetSize(matis->A,NULL,&matis->csf_nleaves);CHKERRQ(ierr); 362a8116848SStefano Zampini ierr = MatGetLocalSize(B,NULL,&matis->csf_nroots);CHKERRQ(ierr); 363a8116848SStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf);CHKERRQ(ierr); 364a8116848SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 365a8116848SStefano Zampini ierr = PetscSFSetGraphLayout(matis->csf,B->cmap,matis->csf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 366a8116848SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 367a8116848SStefano Zampini ierr = PetscMalloc2(matis->csf_nroots,&matis->csf_rootdata,matis->csf_nleaves,&matis->csf_leafdata);CHKERRQ(ierr); 368a8116848SStefano Zampini } else { 369a8116848SStefano Zampini matis->csf = matis->sf; 370a8116848SStefano Zampini matis->csf_nleaves = matis->sf_nleaves; 371a8116848SStefano Zampini matis->csf_nroots = matis->sf_nroots; 372a8116848SStefano Zampini matis->csf_leafdata = matis->sf_leafdata; 373a8116848SStefano Zampini matis->csf_rootdata = matis->sf_rootdata; 374a8116848SStefano Zampini } 37528f4e0baSStefano Zampini PetscFunctionReturn(0); 37628f4e0baSStefano Zampini } 3772e1947a5SStefano Zampini 3782e1947a5SStefano Zampini #undef __FUNCT__ 3792e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 380eb82efa4SStefano Zampini /*@ 381a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 382a88811baSStefano Zampini 383a88811baSStefano Zampini Collective on MPI_Comm 384a88811baSStefano Zampini 385a88811baSStefano Zampini Input Parameters: 386a88811baSStefano Zampini + B - the matrix 387a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 388a88811baSStefano Zampini (same value is used for all local rows) 389a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 390a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 391a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 392a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 393a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 394a88811baSStefano Zampini the diagonal entry even if it is zero. 395a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 396a88811baSStefano Zampini submatrix (same value is used for all local rows). 397a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 398a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 399a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 400a88811baSStefano Zampini structure. The size of this array is equal to the number 401a88811baSStefano Zampini of local rows, i.e 'm'. 402a88811baSStefano Zampini 403a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 404a88811baSStefano Zampini 405a88811baSStefano Zampini Level: intermediate 406a88811baSStefano Zampini 407a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 408a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 409a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 410a88811baSStefano Zampini 411a88811baSStefano Zampini .keywords: matrix 412a88811baSStefano Zampini 4133c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS 414a88811baSStefano Zampini @*/ 4152e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 4162e1947a5SStefano Zampini { 4172e1947a5SStefano Zampini PetscErrorCode ierr; 4182e1947a5SStefano Zampini 4192e1947a5SStefano Zampini PetscFunctionBegin; 4202e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 4212e1947a5SStefano Zampini PetscValidType(B,1); 4222e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 4232e1947a5SStefano Zampini PetscFunctionReturn(0); 4242e1947a5SStefano Zampini } 4252e1947a5SStefano Zampini 4262e1947a5SStefano Zampini #undef __FUNCT__ 4272e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 4287230de76SStefano Zampini static PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 4292e1947a5SStefano Zampini { 4302e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 43128f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 4322e1947a5SStefano Zampini PetscErrorCode ierr; 4332e1947a5SStefano Zampini 4342e1947a5SStefano Zampini PetscFunctionBegin; 4356c4ed002SBarry Smith if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 43628f4e0baSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 43728f4e0baSStefano Zampini ierr = MatISComputeSF_Private(B);CHKERRQ(ierr); 43828f4e0baSStefano Zampini } 4392e1947a5SStefano Zampini if (!d_nnz) { 44028f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz; 4412e1947a5SStefano Zampini } else { 44228f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 4432e1947a5SStefano Zampini } 4442e1947a5SStefano Zampini if (!o_nnz) { 44528f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz; 4462e1947a5SStefano Zampini } else { 44728f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 4482e1947a5SStefano Zampini } 44928f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 45028f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 45128f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 45228f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 45328f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves;i++) { 45428f4e0baSStefano Zampini matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 4552e1947a5SStefano Zampini } 45628f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 45728f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 45828f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 4592e1947a5SStefano Zampini } 46028f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 46128f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 46228f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 4632e1947a5SStefano Zampini } 46428f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 4652e1947a5SStefano Zampini PetscFunctionReturn(0); 4662e1947a5SStefano Zampini } 467b4319ba4SBarry Smith 468b4319ba4SBarry Smith #undef __FUNCT__ 4693927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private" 4703927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 4713927de2eSStefano Zampini { 4723927de2eSStefano Zampini Mat_IS *matis = (Mat_IS*)(A->data); 4733927de2eSStefano Zampini PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 474ecf5a873SStefano Zampini const PetscInt *global_indices_r,*global_indices_c; 4753927de2eSStefano Zampini PetscInt i,j,bs,rows,cols; 4763927de2eSStefano Zampini PetscInt lrows,lcols; 4773927de2eSStefano Zampini PetscInt local_rows,local_cols; 4783927de2eSStefano Zampini PetscMPIInt nsubdomains; 4793927de2eSStefano Zampini PetscBool isdense,issbaij; 4803927de2eSStefano Zampini PetscErrorCode ierr; 4813927de2eSStefano Zampini 4823927de2eSStefano Zampini PetscFunctionBegin; 4833927de2eSStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr); 4843927de2eSStefano Zampini ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr); 4853927de2eSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 4863927de2eSStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 4873927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 4883927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 489ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 490ecf5a873SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 4917230de76SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 492ecf5a873SStefano Zampini } else { 493ecf5a873SStefano Zampini global_indices_c = global_indices_r; 494ecf5a873SStefano Zampini } 495ecf5a873SStefano Zampini 4963927de2eSStefano Zampini if (issbaij) { 4973927de2eSStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 4983927de2eSStefano Zampini } 4993927de2eSStefano Zampini /* 500ecf5a873SStefano Zampini An SF reduce is needed to sum up properly on shared rows. 5013927de2eSStefano Zampini Note that generally preallocation is not exact, since it overestimates nonzeros 5023927de2eSStefano Zampini */ 5033927de2eSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 5043927de2eSStefano Zampini ierr = MatISComputeSF_Private(A);CHKERRQ(ierr); 5053927de2eSStefano Zampini } 5063927de2eSStefano Zampini ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr); 5073927de2eSStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr); 5083927de2eSStefano Zampini /* All processes need to compute entire row ownership */ 5093927de2eSStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 5103927de2eSStefano Zampini ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 5113927de2eSStefano Zampini for (i=0;i<nsubdomains;i++) { 5123927de2eSStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 5133927de2eSStefano Zampini row_ownership[j] = i; 5143927de2eSStefano Zampini } 5153927de2eSStefano Zampini } 5167230de76SStefano Zampini ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 5173927de2eSStefano Zampini 5183927de2eSStefano Zampini /* 5193927de2eSStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 5203927de2eSStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 5213927de2eSStefano Zampini */ 5223927de2eSStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 5233927de2eSStefano Zampini /* preallocation as a MATAIJ */ 5243927de2eSStefano Zampini if (isdense) { /* special case for dense local matrices */ 5253927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 526ecf5a873SStefano Zampini PetscInt index_row = global_indices_r[i]; 5273927de2eSStefano Zampini for (j=i;j<local_rows;j++) { 5283927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 529ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[j]; 5303927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 5313927de2eSStefano Zampini my_dnz[i] += 1; 5323927de2eSStefano Zampini } else { /* offdiag block */ 5333927de2eSStefano Zampini my_onz[i] += 1; 5343927de2eSStefano Zampini } 5353927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 5363927de2eSStefano Zampini if (i != j) { 5373927de2eSStefano Zampini owner = row_ownership[index_col]; 5383927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 5393927de2eSStefano Zampini my_dnz[j] += 1; 5403927de2eSStefano Zampini } else { 5413927de2eSStefano Zampini my_onz[j] += 1; 5423927de2eSStefano Zampini } 5433927de2eSStefano Zampini } 5443927de2eSStefano Zampini } 5453927de2eSStefano Zampini } 546bb1015c3SStefano Zampini } else if (matis->A->ops->getrowij) { 547bb1015c3SStefano Zampini const PetscInt *ii,*jj,*jptr; 548bb1015c3SStefano Zampini PetscBool done; 549bb1015c3SStefano Zampini ierr = MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr); 550bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__); 551bb1015c3SStefano Zampini jptr = jj; 552bb1015c3SStefano Zampini for (i=0;i<local_rows;i++) { 553bb1015c3SStefano Zampini PetscInt index_row = global_indices_r[i]; 554bb1015c3SStefano Zampini for (j=0;j<ii[i+1]-ii[i];j++,jptr++) { 555bb1015c3SStefano Zampini PetscInt owner = row_ownership[index_row]; 556bb1015c3SStefano Zampini PetscInt index_col = global_indices_c[*jptr]; 557bb1015c3SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 558bb1015c3SStefano Zampini my_dnz[i] += 1; 559bb1015c3SStefano Zampini } else { /* offdiag block */ 560bb1015c3SStefano Zampini my_onz[i] += 1; 561bb1015c3SStefano Zampini } 562bb1015c3SStefano Zampini /* same as before, interchanging rows and cols */ 563bb1015c3SStefano Zampini if (issbaij && index_col != index_row) { 564bb1015c3SStefano Zampini owner = row_ownership[index_col]; 565bb1015c3SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 566bb1015c3SStefano Zampini my_dnz[*jptr] += 1; 567bb1015c3SStefano Zampini } else { 568bb1015c3SStefano Zampini my_onz[*jptr] += 1; 569bb1015c3SStefano Zampini } 570bb1015c3SStefano Zampini } 571bb1015c3SStefano Zampini } 572bb1015c3SStefano Zampini } 573bb1015c3SStefano Zampini ierr = MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr); 574bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__); 575bb1015c3SStefano Zampini } else { /* loop over rows and use MatGetRow */ 5763927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 5773927de2eSStefano Zampini const PetscInt *cols; 578ecf5a873SStefano Zampini PetscInt ncols,index_row = global_indices_r[i]; 5793927de2eSStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 5803927de2eSStefano Zampini for (j=0;j<ncols;j++) { 5813927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 582ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[cols[j]]; 5833927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 5843927de2eSStefano Zampini my_dnz[i] += 1; 5853927de2eSStefano Zampini } else { /* offdiag block */ 5863927de2eSStefano Zampini my_onz[i] += 1; 5873927de2eSStefano Zampini } 5883927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 589d9a9e74cSStefano Zampini if (issbaij && index_col != index_row) { 5903927de2eSStefano Zampini owner = row_ownership[index_col]; 5913927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 592d9a9e74cSStefano Zampini my_dnz[cols[j]] += 1; 5933927de2eSStefano Zampini } else { 594d9a9e74cSStefano Zampini my_onz[cols[j]] += 1; 5953927de2eSStefano Zampini } 5963927de2eSStefano Zampini } 5973927de2eSStefano Zampini } 5983927de2eSStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 5993927de2eSStefano Zampini } 6003927de2eSStefano Zampini } 601ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 602ecf5a873SStefano Zampini if (global_indices_c != global_indices_r) { 6037230de76SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 604ecf5a873SStefano Zampini } 6053927de2eSStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 606ecf5a873SStefano Zampini 607ecf5a873SStefano Zampini /* Reduce my_dnz and my_onz */ 6083927de2eSStefano Zampini if (maxreduce) { 6093927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 6103927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 611bb1015c3SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 6123927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 6133927de2eSStefano Zampini } else { 6143927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 6153927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 616bb1015c3SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 6173927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 6183927de2eSStefano Zampini } 6193927de2eSStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 6203927de2eSStefano Zampini 6213927de2eSStefano Zampini /* Resize preallocation if overestimated */ 6223927de2eSStefano Zampini for (i=0;i<lrows;i++) { 6233927de2eSStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 6243927de2eSStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 6253927de2eSStefano Zampini } 6263927de2eSStefano Zampini /* set preallocation */ 6273927de2eSStefano Zampini ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 6283927de2eSStefano Zampini for (i=0;i<lrows/bs;i++) { 6293927de2eSStefano Zampini dnz[i] = dnz[i*bs]/bs; 6303927de2eSStefano Zampini onz[i] = onz[i*bs]/bs; 6313927de2eSStefano Zampini } 6323927de2eSStefano Zampini ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 6333927de2eSStefano Zampini ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 6343927de2eSStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 6353927de2eSStefano Zampini if (issbaij) { 6363927de2eSStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 6373927de2eSStefano Zampini } 6383927de2eSStefano Zampini PetscFunctionReturn(0); 6393927de2eSStefano Zampini } 6403927de2eSStefano Zampini 6413927de2eSStefano Zampini #undef __FUNCT__ 642b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 6437230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 644b7ce53b6SStefano Zampini { 645b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 646d9a9e74cSStefano Zampini Mat local_mat; 647b7ce53b6SStefano Zampini /* info on mat */ 6483cfa4ea4SStefano Zampini PetscInt bs,rows,cols,lrows,lcols; 649b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 650686e3a49SStefano Zampini PetscBool isdense,issbaij,isseqaij; 6517c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 652b7ce53b6SStefano Zampini /* values insertion */ 653b7ce53b6SStefano Zampini PetscScalar *array; 654b7ce53b6SStefano Zampini /* work */ 655b7ce53b6SStefano Zampini PetscErrorCode ierr; 656b7ce53b6SStefano Zampini 657b7ce53b6SStefano Zampini PetscFunctionBegin; 658b7ce53b6SStefano Zampini /* get info from mat */ 6597c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 6607c03b4e8SStefano Zampini if (nsubdomains == 1) { 6617c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 6627c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 6637c03b4e8SStefano Zampini } else { 6647c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 6657c03b4e8SStefano Zampini } 6667c03b4e8SStefano Zampini PetscFunctionReturn(0); 6677c03b4e8SStefano Zampini } 668b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 669b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 6703cfa4ea4SStefano Zampini ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 671b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 672b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 673686e3a49SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 674b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 675b7ce53b6SStefano Zampini 676b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 677b7ce53b6SStefano Zampini MatType new_mat_type; 6783927de2eSStefano Zampini PetscBool issbaij_red; 679b7ce53b6SStefano Zampini 680b7ce53b6SStefano Zampini /* determining new matrix type */ 681b2566f29SBarry Smith ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 682b7ce53b6SStefano Zampini if (issbaij_red) { 683b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 684b7ce53b6SStefano Zampini } else { 685b7ce53b6SStefano Zampini if (bs>1) { 686b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 687b7ce53b6SStefano Zampini } else { 688b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 689b7ce53b6SStefano Zampini } 690b7ce53b6SStefano Zampini } 691b7ce53b6SStefano Zampini 6923927de2eSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 6933cfa4ea4SStefano Zampini ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 6943927de2eSStefano Zampini ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 6953927de2eSStefano Zampini ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr); 6963927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 697b7ce53b6SStefano Zampini } else { 6983cfa4ea4SStefano Zampini PetscInt mbs,mrows,mcols,mlrows,mlcols; 699b7ce53b6SStefano Zampini /* some checks */ 700b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 701b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 7023cfa4ea4SStefano Zampini ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 7036c4ed002SBarry Smith if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 7046c4ed002SBarry Smith if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 7056c4ed002SBarry Smith if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 7066c4ed002SBarry Smith if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 7076c4ed002SBarry Smith if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 708b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 709b7ce53b6SStefano Zampini } 710d9a9e74cSStefano Zampini 711d9a9e74cSStefano Zampini if (issbaij) { 712d9a9e74cSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 713d9a9e74cSStefano Zampini } else { 714d9a9e74cSStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 715d9a9e74cSStefano Zampini local_mat = matis->A; 716d9a9e74cSStefano Zampini } 717686e3a49SStefano Zampini 718b7ce53b6SStefano Zampini /* Set values */ 719ecf5a873SStefano Zampini ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 720b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 721ecf5a873SStefano Zampini PetscInt i,*dummy_rows,*dummy_cols; 722ecf5a873SStefano Zampini 723ecf5a873SStefano Zampini if (local_rows != local_cols) { 724ecf5a873SStefano Zampini ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 725ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 726ecf5a873SStefano Zampini for (i=0;i<local_cols;i++) dummy_cols[i] = i; 727ecf5a873SStefano Zampini } else { 728ecf5a873SStefano Zampini ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 729ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 730ecf5a873SStefano Zampini dummy_cols = dummy_rows; 731ecf5a873SStefano Zampini } 732b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 733d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 734ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 735d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 736ecf5a873SStefano Zampini if (dummy_rows != dummy_cols) { 737ecf5a873SStefano Zampini ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 738ecf5a873SStefano Zampini } else { 739ecf5a873SStefano Zampini ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 740ecf5a873SStefano Zampini } 741686e3a49SStefano Zampini } else if (isseqaij) { 742ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 743686e3a49SStefano Zampini PetscBool done; 744686e3a49SStefano Zampini 745d9a9e74cSStefano Zampini ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 746bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__); 747d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 748686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 749ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 750686e3a49SStefano Zampini } 751d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 752bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__); 753d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 754686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 755ecf5a873SStefano Zampini PetscInt i; 756c0962df8SStefano Zampini 757686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 758686e3a49SStefano Zampini PetscInt j; 759ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 760686e3a49SStefano Zampini 761ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 762ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 763ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 764686e3a49SStefano Zampini } 765b7ce53b6SStefano Zampini } 766b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 767d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 768b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 769b7ce53b6SStefano Zampini if (isdense) { 770b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 771b7ce53b6SStefano Zampini } 772b7ce53b6SStefano Zampini PetscFunctionReturn(0); 773b7ce53b6SStefano Zampini } 774b7ce53b6SStefano Zampini 775b7ce53b6SStefano Zampini #undef __FUNCT__ 776b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 777b7ce53b6SStefano Zampini /*@ 778b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 779b7ce53b6SStefano Zampini 780b7ce53b6SStefano Zampini Input Parameter: 781b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 782b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 783b7ce53b6SStefano Zampini 784b7ce53b6SStefano Zampini Output Parameter: 785b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 786b7ce53b6SStefano Zampini 787b7ce53b6SStefano Zampini Level: developer 788b7ce53b6SStefano Zampini 789eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 790b7ce53b6SStefano Zampini 791b7ce53b6SStefano Zampini .seealso: MATIS 792b7ce53b6SStefano Zampini @*/ 793b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 794b7ce53b6SStefano Zampini { 795b7ce53b6SStefano Zampini PetscErrorCode ierr; 796b7ce53b6SStefano Zampini 797b7ce53b6SStefano Zampini PetscFunctionBegin; 798b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 799b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 800b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 801b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 802b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 803b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 8046c4ed002SBarry Smith if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 805b7ce53b6SStefano Zampini } 806b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 807b7ce53b6SStefano Zampini PetscFunctionReturn(0); 808b7ce53b6SStefano Zampini } 809b7ce53b6SStefano Zampini 810b7ce53b6SStefano Zampini #undef __FUNCT__ 811ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 812ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 813ad6194a2SStefano Zampini { 814ad6194a2SStefano Zampini PetscErrorCode ierr; 815ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 816ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 817ad6194a2SStefano Zampini Mat B,localmat; 818ad6194a2SStefano Zampini 819ad6194a2SStefano Zampini PetscFunctionBegin; 820ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 821ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 822ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 823e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 824ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 825ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 826b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 827ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 828ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 829ad6194a2SStefano Zampini *newmat = B; 830ad6194a2SStefano Zampini PetscFunctionReturn(0); 831ad6194a2SStefano Zampini } 832ad6194a2SStefano Zampini 833ad6194a2SStefano Zampini #undef __FUNCT__ 83469796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 835a8116848SStefano Zampini static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 83669796d55SStefano Zampini { 83769796d55SStefano Zampini PetscErrorCode ierr; 83869796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 83969796d55SStefano Zampini PetscBool local_sym; 84069796d55SStefano Zampini 84169796d55SStefano Zampini PetscFunctionBegin; 84269796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 843b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 84469796d55SStefano Zampini PetscFunctionReturn(0); 84569796d55SStefano Zampini } 84669796d55SStefano Zampini 84769796d55SStefano Zampini #undef __FUNCT__ 84869796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 849a8116848SStefano Zampini static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 85069796d55SStefano Zampini { 85169796d55SStefano Zampini PetscErrorCode ierr; 85269796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 85369796d55SStefano Zampini PetscBool local_sym; 85469796d55SStefano Zampini 85569796d55SStefano Zampini PetscFunctionBegin; 85669796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 857b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 85869796d55SStefano Zampini PetscFunctionReturn(0); 85969796d55SStefano Zampini } 86069796d55SStefano Zampini 86169796d55SStefano Zampini #undef __FUNCT__ 862b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 863a8116848SStefano Zampini static PetscErrorCode MatDestroy_IS(Mat A) 864b4319ba4SBarry Smith { 865dfbe8321SBarry Smith PetscErrorCode ierr; 866b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 867b4319ba4SBarry Smith 868b4319ba4SBarry Smith PetscFunctionBegin; 8696bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 870e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 871e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 8726bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 8736bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 874*3fd1c9e7SStefano Zampini ierr = VecDestroy(&b->counter);CHKERRQ(ierr); 875a8116848SStefano Zampini ierr = ISDestroy(&b->getsub_ris);CHKERRQ(ierr); 876a8116848SStefano Zampini ierr = ISDestroy(&b->getsub_cis);CHKERRQ(ierr); 877a8116848SStefano Zampini if (b->sf != b->csf) { 878a8116848SStefano Zampini ierr = PetscSFDestroy(&b->csf);CHKERRQ(ierr); 879a8116848SStefano Zampini ierr = PetscFree2(b->csf_rootdata,b->csf_leafdata);CHKERRQ(ierr); 880a8116848SStefano Zampini } 88128f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 88228f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 883bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 884dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 885bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 886b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 887b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 8882e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 889b4319ba4SBarry Smith PetscFunctionReturn(0); 890b4319ba4SBarry Smith } 891b4319ba4SBarry Smith 892b4319ba4SBarry Smith #undef __FUNCT__ 893b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 894a8116848SStefano Zampini static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 895b4319ba4SBarry Smith { 896dfbe8321SBarry Smith PetscErrorCode ierr; 897b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 898b4319ba4SBarry Smith PetscScalar zero = 0.0; 899b4319ba4SBarry Smith 900b4319ba4SBarry Smith PetscFunctionBegin; 901b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 902e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 903e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 904b4319ba4SBarry Smith 905b4319ba4SBarry Smith /* multiply the local matrix */ 906b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 907b4319ba4SBarry Smith 908b4319ba4SBarry Smith /* scatter product back into global memory */ 9092dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 910e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 911e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 912b4319ba4SBarry Smith PetscFunctionReturn(0); 913b4319ba4SBarry Smith } 914b4319ba4SBarry Smith 915b4319ba4SBarry Smith #undef __FUNCT__ 9162e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 917a8116848SStefano Zampini static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 9182e74eeadSLisandro Dalcin { 919650997f4SStefano Zampini Vec temp_vec; 9202e74eeadSLisandro Dalcin PetscErrorCode ierr; 9212e74eeadSLisandro Dalcin 9222e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 923650997f4SStefano Zampini if (v3 != v2) { 924650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 925650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 926650997f4SStefano Zampini } else { 927650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 928650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 929650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 930650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 931650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 932650997f4SStefano Zampini } 9332e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9342e74eeadSLisandro Dalcin } 9352e74eeadSLisandro Dalcin 9362e74eeadSLisandro Dalcin #undef __FUNCT__ 9372e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 938a8116848SStefano Zampini static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 9392e74eeadSLisandro Dalcin { 9402e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9412e74eeadSLisandro Dalcin PetscErrorCode ierr; 9422e74eeadSLisandro Dalcin 943e176bc59SStefano Zampini PetscFunctionBegin; 9442e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 945e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 946e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9472e74eeadSLisandro Dalcin 9482e74eeadSLisandro Dalcin /* multiply the local matrix */ 949e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 9502e74eeadSLisandro Dalcin 9512e74eeadSLisandro Dalcin /* scatter product back into global vector */ 952e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 953e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 954e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9552e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9562e74eeadSLisandro Dalcin } 9572e74eeadSLisandro Dalcin 9582e74eeadSLisandro Dalcin #undef __FUNCT__ 9592e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 960a8116848SStefano Zampini static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 9612e74eeadSLisandro Dalcin { 962650997f4SStefano Zampini Vec temp_vec; 9632e74eeadSLisandro Dalcin PetscErrorCode ierr; 9642e74eeadSLisandro Dalcin 9652e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 966650997f4SStefano Zampini if (v3 != v2) { 967650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 968650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 969650997f4SStefano Zampini } else { 970650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 971650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 972650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 973650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 974650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 975650997f4SStefano Zampini } 9762e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9772e74eeadSLisandro Dalcin } 9782e74eeadSLisandro Dalcin 9792e74eeadSLisandro Dalcin #undef __FUNCT__ 980b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 981a8116848SStefano Zampini static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 982b4319ba4SBarry Smith { 983b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 984dfbe8321SBarry Smith PetscErrorCode ierr; 985b4319ba4SBarry Smith PetscViewer sviewer; 986b4319ba4SBarry Smith 987b4319ba4SBarry Smith PetscFunctionBegin; 9883f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 989b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 9903f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 9916e520ac8SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 992b4319ba4SBarry Smith PetscFunctionReturn(0); 993b4319ba4SBarry Smith } 994b4319ba4SBarry Smith 995b4319ba4SBarry Smith #undef __FUNCT__ 996b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 997a8116848SStefano Zampini static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 998b4319ba4SBarry Smith { 999dfbe8321SBarry Smith PetscErrorCode ierr; 1000e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 1001b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1002b4319ba4SBarry Smith IS from,to; 1003e176bc59SStefano Zampini Vec cglobal,rglobal; 1004b4319ba4SBarry Smith 1005b4319ba4SBarry Smith PetscFunctionBegin; 1006784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 1007e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 10083bbff08aSStefano Zampini /* Destroy any previous data */ 100970cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 101070cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 1011*3fd1c9e7SStefano Zampini ierr = VecDestroy(&is->counter);CHKERRQ(ierr); 1012e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 1013e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 10141c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 101528f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 101628f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 10173bbff08aSStefano Zampini 10183bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 1019fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1020fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1021fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 1022fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 1023b4319ba4SBarry Smith 1024b4319ba4SBarry Smith /* Create the local matrix A */ 1025e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 1026e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 1027e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 1028e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 1029f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 1030e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 1031e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 1032e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 1033ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 1034ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 1035b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 1036b4319ba4SBarry Smith 1037f26d0771SStefano Zampini if (!is->islocalref) { /* setup scatters and local vectors for MatMult */ 1038b4319ba4SBarry Smith /* Create the local work vectors */ 10393bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 1040b4319ba4SBarry Smith 1041e176bc59SStefano Zampini /* setup the global to local scatters */ 1042e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 1043e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 1044784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 1045e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 1046e176bc59SStefano Zampini if (rmapping != cmapping) { 1047e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 1048e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 1049e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 1050e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 1051e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 1052e176bc59SStefano Zampini } else { 1053e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 1054e176bc59SStefano Zampini is->cctx = is->rctx; 1055e176bc59SStefano Zampini } 1056*3fd1c9e7SStefano Zampini 1057*3fd1c9e7SStefano Zampini /* interface counter vector (local) */ 1058*3fd1c9e7SStefano Zampini ierr = VecDuplicate(is->y,&is->counter);CHKERRQ(ierr); 1059*3fd1c9e7SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 1060*3fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1061*3fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1062*3fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1063*3fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1064*3fd1c9e7SStefano Zampini 1065*3fd1c9e7SStefano Zampini /* free workspace */ 1066e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 1067e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 10686bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 10696bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1070f26d0771SStefano Zampini } 107148ff6bf3SStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 1072b4319ba4SBarry Smith PetscFunctionReturn(0); 1073b4319ba4SBarry Smith } 1074b4319ba4SBarry Smith 10752e74eeadSLisandro Dalcin #undef __FUNCT__ 10762e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 1077a8116848SStefano Zampini static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 10782e74eeadSLisandro Dalcin { 10792e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 10802e74eeadSLisandro Dalcin PetscErrorCode ierr; 108197563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 108297563a80SStefano Zampini PetscInt i,zm,zn; 108397563a80SStefano Zampini #endif 1084f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 10852e74eeadSLisandro Dalcin 10862e74eeadSLisandro Dalcin PetscFunctionBegin; 10872e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 1088f26d0771SStefano 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); 108997563a80SStefano Zampini /* count negative indices */ 109097563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 109197563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 10922e74eeadSLisandro Dalcin #endif 109397563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 109497563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 109597563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 109697563a80SStefano Zampini /* count negative indices (should be the same as before) */ 109797563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 109897563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 109997563a80SStefano 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"); 110097563a80SStefano 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"); 110197563a80SStefano Zampini #endif 11022e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 11032e74eeadSLisandro Dalcin PetscFunctionReturn(0); 11042e74eeadSLisandro Dalcin } 11052e74eeadSLisandro Dalcin 1106b4319ba4SBarry Smith #undef __FUNCT__ 110797563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS" 1108a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 110997563a80SStefano Zampini { 111097563a80SStefano Zampini Mat_IS *is = (Mat_IS*)mat->data; 111197563a80SStefano Zampini PetscErrorCode ierr; 111297563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 111397563a80SStefano Zampini PetscInt i,zm,zn; 111497563a80SStefano Zampini #endif 1115f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 111697563a80SStefano Zampini 111797563a80SStefano Zampini PetscFunctionBegin; 111897563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 1119f26d0771SStefano 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); 112097563a80SStefano Zampini /* count negative indices */ 112197563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 112297563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 112397563a80SStefano Zampini #endif 112497563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 112597563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 112697563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 112797563a80SStefano Zampini /* count negative indices (should be the same as before) */ 112897563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 112997563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 113097563a80SStefano 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"); 113197563a80SStefano 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"); 113297563a80SStefano Zampini #endif 113397563a80SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 113497563a80SStefano Zampini PetscFunctionReturn(0); 113597563a80SStefano Zampini } 113697563a80SStefano Zampini 113797563a80SStefano Zampini #undef __FUNCT__ 1138b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 1139a8116848SStefano Zampini static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1140b4319ba4SBarry Smith { 1141dfbe8321SBarry Smith PetscErrorCode ierr; 1142b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1143b4319ba4SBarry Smith 1144b4319ba4SBarry Smith PetscFunctionBegin; 1145b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1146b4319ba4SBarry Smith PetscFunctionReturn(0); 1147b4319ba4SBarry Smith } 1148b4319ba4SBarry Smith 1149b4319ba4SBarry Smith #undef __FUNCT__ 1150f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 1151a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1152f0006bf2SLisandro Dalcin { 1153f0006bf2SLisandro Dalcin PetscErrorCode ierr; 1154f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 1155f0006bf2SLisandro Dalcin 1156f0006bf2SLisandro Dalcin PetscFunctionBegin; 1157f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1158f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 1159f0006bf2SLisandro Dalcin } 1160f0006bf2SLisandro Dalcin 1161f0006bf2SLisandro Dalcin #undef __FUNCT__ 11622e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 1163a8116848SStefano Zampini static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 11642e74eeadSLisandro Dalcin { 11656e520ac8SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 11666e520ac8SStefano Zampini PetscInt nr,nl,len,i; 11676e520ac8SStefano Zampini PetscInt *lrows; 11682e74eeadSLisandro Dalcin PetscErrorCode ierr; 11692e74eeadSLisandro Dalcin 11702e74eeadSLisandro Dalcin PetscFunctionBegin; 11716e520ac8SStefano Zampini /* get locally owned rows */ 11726e520ac8SStefano Zampini ierr = MatZeroRowsMapLocal_Private(A,n,rows,&len,&lrows);CHKERRQ(ierr); 11736e520ac8SStefano Zampini /* fix right hand side if needed */ 11746e520ac8SStefano Zampini if (x && b) { 11756e520ac8SStefano Zampini const PetscScalar *xx; 11766e520ac8SStefano Zampini PetscScalar *bb; 11776e520ac8SStefano Zampini 11786e520ac8SStefano Zampini ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr); 11796e520ac8SStefano Zampini ierr = VecGetArray(b, &bb);CHKERRQ(ierr); 11806e520ac8SStefano Zampini for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]]; 11816e520ac8SStefano Zampini ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr); 11826e520ac8SStefano Zampini ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr); 11832e74eeadSLisandro Dalcin } 11846e520ac8SStefano Zampini /* get rows associated to the local matrices */ 11856e520ac8SStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 11866e520ac8SStefano Zampini ierr = MatISComputeSF_Private(A);CHKERRQ(ierr); 11876e520ac8SStefano Zampini } 11886e520ac8SStefano Zampini ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr); 11896e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr); 11906e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 11916e520ac8SStefano Zampini for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1; 11926e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 11936e520ac8SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 11946e520ac8SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 11956e520ac8SStefano Zampini ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr); 11966e520ac8SStefano Zampini for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i; 11977230de76SStefano Zampini ierr = MatISZeroRowsLocal_Private(A,nr,lrows,diag);CHKERRQ(ierr); 11986e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 11992e74eeadSLisandro Dalcin PetscFunctionReturn(0); 12002e74eeadSLisandro Dalcin } 12012e74eeadSLisandro Dalcin 12022e74eeadSLisandro Dalcin #undef __FUNCT__ 12037230de76SStefano Zampini #define __FUNCT__ "MatISZeroRowsLocal_Private" 12047230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag) 1205b4319ba4SBarry Smith { 1206b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1207dfbe8321SBarry Smith PetscErrorCode ierr; 1208b4319ba4SBarry Smith 1209b4319ba4SBarry Smith PetscFunctionBegin; 1210958c9bccSBarry Smith if (!n) { 1211b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 1212b4319ba4SBarry Smith } else { 12137230de76SStefano Zampini PetscInt i; 1214b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 12152205254eSKarl Rupp 12162b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 12176e520ac8SStefano Zampini if (diag != 0.) { 12186e520ac8SStefano Zampini const PetscScalar *array; 1219*3fd1c9e7SStefano Zampini ierr = VecGetArrayRead(is->counter,&array);CHKERRQ(ierr); 1220b4319ba4SBarry Smith for (i=0; i<n; i++) { 1221f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 1222b4319ba4SBarry Smith } 1223*3fd1c9e7SStefano Zampini ierr = VecRestoreArrayRead(is->counter,&array);CHKERRQ(ierr); 12246e520ac8SStefano Zampini } 1225b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1226b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1227b4319ba4SBarry Smith } 1228b4319ba4SBarry Smith PetscFunctionReturn(0); 1229b4319ba4SBarry Smith } 1230b4319ba4SBarry Smith 1231b4319ba4SBarry Smith #undef __FUNCT__ 1232b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 1233a8116848SStefano Zampini static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 1234b4319ba4SBarry Smith { 1235b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1236dfbe8321SBarry Smith PetscErrorCode ierr; 1237b4319ba4SBarry Smith 1238b4319ba4SBarry Smith PetscFunctionBegin; 1239b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 1240b4319ba4SBarry Smith PetscFunctionReturn(0); 1241b4319ba4SBarry Smith } 1242b4319ba4SBarry Smith 1243b4319ba4SBarry Smith #undef __FUNCT__ 1244b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 1245a8116848SStefano Zampini static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 1246b4319ba4SBarry Smith { 1247b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1248dfbe8321SBarry Smith PetscErrorCode ierr; 1249b4319ba4SBarry Smith 1250b4319ba4SBarry Smith PetscFunctionBegin; 1251b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 1252b4319ba4SBarry Smith PetscFunctionReturn(0); 1253b4319ba4SBarry Smith } 1254b4319ba4SBarry Smith 1255b4319ba4SBarry Smith #undef __FUNCT__ 1256b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 1257a8116848SStefano Zampini static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 1258b4319ba4SBarry Smith { 1259b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 1260b4319ba4SBarry Smith 1261b4319ba4SBarry Smith PetscFunctionBegin; 1262b4319ba4SBarry Smith *local = is->A; 1263b4319ba4SBarry Smith PetscFunctionReturn(0); 1264b4319ba4SBarry Smith } 1265b4319ba4SBarry Smith 1266b4319ba4SBarry Smith #undef __FUNCT__ 1267b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 1268b4319ba4SBarry Smith /*@ 1269b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 1270b4319ba4SBarry Smith 1271b4319ba4SBarry Smith Input Parameter: 1272b4319ba4SBarry Smith . mat - the matrix 1273b4319ba4SBarry Smith 1274b4319ba4SBarry Smith Output Parameter: 1275eb82efa4SStefano Zampini . local - the local matrix 1276b4319ba4SBarry Smith 1277b4319ba4SBarry Smith Level: advanced 1278b4319ba4SBarry Smith 1279b4319ba4SBarry Smith Notes: 1280b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 1281b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 1282b4319ba4SBarry Smith of the MatSetValues() operation. 1283b4319ba4SBarry Smith 1284b4319ba4SBarry Smith .seealso: MATIS 1285b4319ba4SBarry Smith @*/ 12867087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 1287b4319ba4SBarry Smith { 12884ac538c5SBarry Smith PetscErrorCode ierr; 1289b4319ba4SBarry Smith 1290b4319ba4SBarry Smith PetscFunctionBegin; 12910700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1292b4319ba4SBarry Smith PetscValidPointer(local,2); 12934ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 1294b4319ba4SBarry Smith PetscFunctionReturn(0); 1295b4319ba4SBarry Smith } 1296b4319ba4SBarry Smith 12973b03a366Sstefano_zampini #undef __FUNCT__ 12983b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 1299a8116848SStefano Zampini static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 13003b03a366Sstefano_zampini { 13013b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 13023b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 13033b03a366Sstefano_zampini PetscErrorCode ierr; 13043b03a366Sstefano_zampini 13053b03a366Sstefano_zampini PetscFunctionBegin; 13064e4c7dbeSStefano Zampini if (is->A) { 13073b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 13083b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 1309f23aa3ddSBarry Smith if (orows != nrows || ocols != ncols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local MATIS matrix should be of size %dx%d (you passed a %dx%d matrix)\n",orows,ocols,nrows,ncols); 13104e4c7dbeSStefano Zampini } 13113b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 13123b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 13133b03a366Sstefano_zampini is->A = local; 13143b03a366Sstefano_zampini PetscFunctionReturn(0); 13153b03a366Sstefano_zampini } 13163b03a366Sstefano_zampini 13173b03a366Sstefano_zampini #undef __FUNCT__ 13183b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 13193b03a366Sstefano_zampini /*@ 1320eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 13213b03a366Sstefano_zampini 13223b03a366Sstefano_zampini Input Parameter: 13233b03a366Sstefano_zampini . mat - the matrix 1324eb82efa4SStefano Zampini . local - the local matrix 13253b03a366Sstefano_zampini 13263b03a366Sstefano_zampini Output Parameter: 13273b03a366Sstefano_zampini 13283b03a366Sstefano_zampini Level: advanced 13293b03a366Sstefano_zampini 13303b03a366Sstefano_zampini Notes: 13313b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 13323b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 13333b03a366Sstefano_zampini 13343b03a366Sstefano_zampini .seealso: MATIS 13353b03a366Sstefano_zampini @*/ 13363b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 13373b03a366Sstefano_zampini { 13383b03a366Sstefano_zampini PetscErrorCode ierr; 13393b03a366Sstefano_zampini 13403b03a366Sstefano_zampini PetscFunctionBegin; 13413b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1342b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 13433b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 13443b03a366Sstefano_zampini PetscFunctionReturn(0); 13453b03a366Sstefano_zampini } 13463b03a366Sstefano_zampini 13476726f965SBarry Smith #undef __FUNCT__ 13486726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 1349a8116848SStefano Zampini static PetscErrorCode MatZeroEntries_IS(Mat A) 13506726f965SBarry Smith { 13516726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 13526726f965SBarry Smith PetscErrorCode ierr; 13536726f965SBarry Smith 13546726f965SBarry Smith PetscFunctionBegin; 13556726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 13566726f965SBarry Smith PetscFunctionReturn(0); 13576726f965SBarry Smith } 13586726f965SBarry Smith 13596726f965SBarry Smith #undef __FUNCT__ 13602e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 1361a8116848SStefano Zampini static PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 13622e74eeadSLisandro Dalcin { 13632e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 13642e74eeadSLisandro Dalcin PetscErrorCode ierr; 13652e74eeadSLisandro Dalcin 13662e74eeadSLisandro Dalcin PetscFunctionBegin; 13672e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 13682e74eeadSLisandro Dalcin PetscFunctionReturn(0); 13692e74eeadSLisandro Dalcin } 13702e74eeadSLisandro Dalcin 13712e74eeadSLisandro Dalcin #undef __FUNCT__ 13722e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 1373a8116848SStefano Zampini static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 13742e74eeadSLisandro Dalcin { 13752e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 13762e74eeadSLisandro Dalcin PetscErrorCode ierr; 13772e74eeadSLisandro Dalcin 13782e74eeadSLisandro Dalcin PetscFunctionBegin; 13792e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 1380e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 13812e74eeadSLisandro Dalcin 13822e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 13832e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 1384e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1385e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 13862e74eeadSLisandro Dalcin PetscFunctionReturn(0); 13872e74eeadSLisandro Dalcin } 13882e74eeadSLisandro Dalcin 13892e74eeadSLisandro Dalcin #undef __FUNCT__ 13906726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 1391a8116848SStefano Zampini static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 13926726f965SBarry Smith { 13936726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 13946726f965SBarry Smith PetscErrorCode ierr; 13956726f965SBarry Smith 13966726f965SBarry Smith PetscFunctionBegin; 13974e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 13986726f965SBarry Smith PetscFunctionReturn(0); 13996726f965SBarry Smith } 14006726f965SBarry Smith 1401284134d9SBarry Smith #undef __FUNCT__ 1402f26d0771SStefano Zampini #define __FUNCT__ "MatAXPY_IS" 1403f26d0771SStefano Zampini static PetscErrorCode MatAXPY_IS(Mat Y,PetscScalar a,Mat X,MatStructure str) 1404f26d0771SStefano Zampini { 1405f26d0771SStefano Zampini Mat_IS *y = (Mat_IS*)Y->data; 1406f26d0771SStefano Zampini Mat_IS *x; 1407f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1408f26d0771SStefano Zampini PetscBool ismatis; 1409f26d0771SStefano Zampini #endif 1410f26d0771SStefano Zampini PetscErrorCode ierr; 1411f26d0771SStefano Zampini 1412f26d0771SStefano Zampini PetscFunctionBegin; 1413f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1414f26d0771SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)X,MATIS,&ismatis);CHKERRQ(ierr); 1415f26d0771SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS"); 1416f26d0771SStefano Zampini #endif 1417f26d0771SStefano Zampini x = (Mat_IS*)X->data; 1418f26d0771SStefano Zampini ierr = MatAXPY(y->A,a,x->A,str);CHKERRQ(ierr); 1419f26d0771SStefano Zampini PetscFunctionReturn(0); 1420f26d0771SStefano Zampini } 1421f26d0771SStefano Zampini 1422f26d0771SStefano Zampini #undef __FUNCT__ 1423f26d0771SStefano Zampini #define __FUNCT__ "MatGetLocalSubMatrix_IS" 1424f26d0771SStefano Zampini static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat) 1425f26d0771SStefano Zampini { 1426f26d0771SStefano Zampini Mat lA; 1427f26d0771SStefano Zampini Mat_IS *matis; 1428f26d0771SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 1429f26d0771SStefano Zampini IS is; 1430f26d0771SStefano Zampini const PetscInt *rg,*rl; 1431f26d0771SStefano Zampini PetscInt nrg; 1432f26d0771SStefano Zampini PetscInt N,M,nrl,i,*idxs; 1433f26d0771SStefano Zampini PetscErrorCode ierr; 1434f26d0771SStefano Zampini 1435f26d0771SStefano Zampini PetscFunctionBegin; 1436f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 1437f26d0771SStefano Zampini ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr); 1438f26d0771SStefano Zampini ierr = ISGetIndices(row,&rl);CHKERRQ(ierr); 1439f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr); 1440f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1441f26d0771SStefano Zampini for (i=0;i<nrl;i++) if (rl[i]>=nrg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local row index %d -> %d greater than maximum possible %d\n",i,rl[i],nrg); 1442f26d0771SStefano Zampini #endif 1443f26d0771SStefano Zampini ierr = PetscMalloc1(nrg,&idxs);CHKERRQ(ierr); 1444f26d0771SStefano Zampini /* map from [0,nrl) to row */ 1445f26d0771SStefano Zampini for (i=0;i<nrl;i++) idxs[i] = rl[i]; 1446f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1447f26d0771SStefano Zampini for (i=nrl;i<nrg;i++) idxs[i] = nrg; 1448f26d0771SStefano Zampini #else 1449f26d0771SStefano Zampini for (i=nrl;i<nrg;i++) idxs[i] = -1; 1450f26d0771SStefano Zampini #endif 1451f26d0771SStefano Zampini ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr); 1452f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 1453f26d0771SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 1454f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 1455f26d0771SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 1456f26d0771SStefano Zampini /* compute new l2g map for columns */ 1457f26d0771SStefano Zampini if (col != row || A->rmap->mapping != A->cmap->mapping) { 1458f26d0771SStefano Zampini const PetscInt *cg,*cl; 1459f26d0771SStefano Zampini PetscInt ncg; 1460f26d0771SStefano Zampini PetscInt ncl; 1461f26d0771SStefano Zampini 1462f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 1463f26d0771SStefano Zampini ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr); 1464f26d0771SStefano Zampini ierr = ISGetIndices(col,&cl);CHKERRQ(ierr); 1465f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr); 1466f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1467f26d0771SStefano Zampini for (i=0;i<ncl;i++) if (cl[i]>=ncg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local column index %d -> %d greater than maximum possible %d\n",i,cl[i],ncg); 1468f26d0771SStefano Zampini #endif 1469f26d0771SStefano Zampini ierr = PetscMalloc1(ncg,&idxs);CHKERRQ(ierr); 1470f26d0771SStefano Zampini /* map from [0,ncl) to col */ 1471f26d0771SStefano Zampini for (i=0;i<ncl;i++) idxs[i] = cl[i]; 1472f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1473f26d0771SStefano Zampini for (i=ncl;i<ncg;i++) idxs[i] = ncg; 1474f26d0771SStefano Zampini #else 1475f26d0771SStefano Zampini for (i=ncl;i<ncg;i++) idxs[i] = -1; 1476f26d0771SStefano Zampini #endif 1477f26d0771SStefano Zampini ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr); 1478f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 1479f26d0771SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 1480f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 1481f26d0771SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 1482f26d0771SStefano Zampini } else { 1483f26d0771SStefano Zampini ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr); 1484f26d0771SStefano Zampini cl2g = rl2g; 1485f26d0771SStefano Zampini } 1486f26d0771SStefano Zampini /* create the MATIS submatrix */ 1487f26d0771SStefano Zampini ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 1488f26d0771SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr); 1489f26d0771SStefano Zampini ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1490f26d0771SStefano Zampini ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr); 1491b0aa3428SStefano Zampini matis = (Mat_IS*)((*submat)->data); 1492f26d0771SStefano Zampini matis->islocalref = PETSC_TRUE; 1493f26d0771SStefano Zampini ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr); 1494f26d0771SStefano Zampini ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr); 1495f26d0771SStefano Zampini ierr = MatISSetLocalMat(*submat,lA);CHKERRQ(ierr); 1496f26d0771SStefano Zampini ierr = MatSetUp(*submat);CHKERRQ(ierr); 1497f26d0771SStefano Zampini ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1498f26d0771SStefano Zampini ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1499f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 1500f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 1501f26d0771SStefano Zampini /* remove unsupported ops */ 1502f26d0771SStefano Zampini ierr = PetscMemzero((*submat)->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1503f26d0771SStefano Zampini (*submat)->ops->destroy = MatDestroy_IS; 1504f26d0771SStefano Zampini (*submat)->ops->setvalueslocal = MatSetValuesLocal_SubMat_IS; 1505f26d0771SStefano Zampini (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS; 1506f26d0771SStefano Zampini (*submat)->ops->assemblybegin = MatAssemblyBegin_IS; 1507f26d0771SStefano Zampini (*submat)->ops->assemblyend = MatAssemblyEnd_IS; 1508f26d0771SStefano Zampini PetscFunctionReturn(0); 1509f26d0771SStefano Zampini } 1510f26d0771SStefano Zampini 1511f26d0771SStefano Zampini #undef __FUNCT__ 1512284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 1513284134d9SBarry Smith /*@ 15143c212e90SHong Zhang MatCreateIS - Creates a "process" unassembled matrix, assembled on each 1515284134d9SBarry Smith process but not across processes. 1516284134d9SBarry Smith 1517284134d9SBarry Smith Input Parameters: 1518284134d9SBarry Smith + comm - MPI communicator that will share the matrix 1519e176bc59SStefano Zampini . bs - block size of the matrix 1520df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 1521e176bc59SStefano Zampini . rmap - local to global map for rows 1522e176bc59SStefano Zampini - cmap - local to global map for cols 1523284134d9SBarry Smith 1524284134d9SBarry Smith Output Parameter: 1525284134d9SBarry Smith . A - the resulting matrix 1526284134d9SBarry Smith 15278e6c10adSSatish Balay Level: advanced 15288e6c10adSSatish Balay 15293c212e90SHong Zhang Notes: See MATIS for more details. 15303c212e90SHong Zhang m and n are NOT related to the size of the map; they are the size of the part of the vector owned 1531e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 15323c212e90SHong Zhang If either rmap or cmap are NULL, then the matrix is assumed to be square. 1533284134d9SBarry Smith 1534284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 1535284134d9SBarry Smith @*/ 1536e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 1537284134d9SBarry Smith { 1538284134d9SBarry Smith PetscErrorCode ierr; 1539284134d9SBarry Smith 1540284134d9SBarry Smith PetscFunctionBegin; 1541e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 1542284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1543d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 1544284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1545284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 1546e176bc59SStefano Zampini if (rmap && cmap) { 1547e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 1548e176bc59SStefano Zampini } else if (!rmap) { 1549e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 1550e176bc59SStefano Zampini } else { 1551e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 1552e176bc59SStefano Zampini } 1553284134d9SBarry Smith PetscFunctionReturn(0); 1554284134d9SBarry Smith } 1555284134d9SBarry Smith 1556b4319ba4SBarry Smith /*MC 1557f26d0771SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition methods (e.g. PCBDDC or KSPFETIDP). 1558b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 1559b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 1560b4319ba4SBarry Smith product is handled "implicitly". 1561b4319ba4SBarry Smith 1562b4319ba4SBarry Smith Operations Provided: 15636726f965SBarry Smith + MatMult() 15642e74eeadSLisandro Dalcin . MatMultAdd() 15652e74eeadSLisandro Dalcin . MatMultTranspose() 15662e74eeadSLisandro Dalcin . MatMultTransposeAdd() 15676726f965SBarry Smith . MatZeroEntries() 15686726f965SBarry Smith . MatSetOption() 15692e74eeadSLisandro Dalcin . MatZeroRows() 15702e74eeadSLisandro Dalcin . MatSetValues() 157197563a80SStefano Zampini . MatSetValuesBlocked() 15726726f965SBarry Smith . MatSetValuesLocal() 157397563a80SStefano Zampini . MatSetValuesBlockedLocal() 15742e74eeadSLisandro Dalcin . MatScale() 15752e74eeadSLisandro Dalcin . MatGetDiagonal() 15762b404112SStefano Zampini . MatMissingDiagonal() 15772b404112SStefano Zampini . MatDuplicate() 15782b404112SStefano Zampini . MatCopy() 1579f26d0771SStefano Zampini . MatAXPY() 1580f26d0771SStefano Zampini . MatGetSubMatrix() 1581f26d0771SStefano Zampini . MatGetLocalSubMatrix() 15826726f965SBarry Smith - MatSetLocalToGlobalMapping() 1583b4319ba4SBarry Smith 1584b4319ba4SBarry Smith Options Database Keys: 1585b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1586b4319ba4SBarry Smith 1587b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1588b4319ba4SBarry Smith 1589b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1590b4319ba4SBarry Smith 1591b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1592eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1593b4319ba4SBarry Smith 1594b4319ba4SBarry Smith Level: advanced 1595b4319ba4SBarry Smith 1596f26d0771SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), MatCreateIS(), PCBDDC, KSPFETIDP 1597b4319ba4SBarry Smith 1598b4319ba4SBarry Smith M*/ 1599b4319ba4SBarry Smith 1600b4319ba4SBarry Smith #undef __FUNCT__ 1601b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 16028cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1603b4319ba4SBarry Smith { 1604dfbe8321SBarry Smith PetscErrorCode ierr; 1605b4319ba4SBarry Smith Mat_IS *b; 1606b4319ba4SBarry Smith 1607b4319ba4SBarry Smith PetscFunctionBegin; 1608b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1609b4319ba4SBarry Smith A->data = (void*)b; 1610b4319ba4SBarry Smith 1611e176bc59SStefano Zampini /* matrix ops */ 1612e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1613b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 16142e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 16152e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 16162e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1617b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1618b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 16192e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 162098921651SStefano Zampini A->ops->setvaluesblocked = MatSetValuesBlocked_IS; 1621b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1622f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 16232e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1624b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1625b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1626b4319ba4SBarry Smith A->ops->view = MatView_IS; 16276726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 16282e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 16292e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 16306726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 163169796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 163269796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1633ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 16346bd84002SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_IS; 16352b404112SStefano Zampini A->ops->copy = MatCopy_IS; 1636659959c5SStefano Zampini A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_IS; 1637a8116848SStefano Zampini A->ops->getsubmatrix = MatGetSubMatrix_IS; 1638f26d0771SStefano Zampini A->ops->axpy = MatAXPY_IS; 1639*3fd1c9e7SStefano Zampini A->ops->diagonalset = MatDiagonalSet_IS; 1640*3fd1c9e7SStefano Zampini A->ops->shift = MatShift_IS; 1641b4319ba4SBarry Smith 1642b7ce53b6SStefano Zampini /* special MATIS functions */ 1643bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1644bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1645b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 16462e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 164717667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1648b4319ba4SBarry Smith PetscFunctionReturn(0); 1649b4319ba4SBarry Smith } 1650