1be1d678aSKris Buschelman 2b4319ba4SBarry Smith /* 3b4319ba4SBarry Smith Creates a matrix class for using the Neumann-Neumann type preconditioners. 4b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 5b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 6b4319ba4SBarry Smith product is handled "implicitly". 7b4319ba4SBarry Smith 8b4319ba4SBarry Smith Currently this allows for only one subdomain per processor. 9b4319ba4SBarry Smith */ 10b4319ba4SBarry Smith 11c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/ 1228f4e0baSStefano Zampini 13a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat); 147230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat,PetscInt,const PetscInt[],PetscScalar); 15a72627d2SStefano Zampini 1628f4e0baSStefano Zampini #undef __FUNCT__ 17a8116848SStefano Zampini #define __FUNCT__ "PetscLayoutMapLocal_Private" 18a8116848SStefano Zampini static PetscErrorCode PetscLayoutMapLocal_Private(PetscLayout map,PetscInt N,const PetscInt idxs[],const PetscInt gidxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs) 19a8116848SStefano Zampini { 20a8116848SStefano Zampini PetscInt *owners = map->range; 21a8116848SStefano Zampini PetscInt n = map->n; 22a8116848SStefano Zampini PetscSF sf; 23a8116848SStefano Zampini PetscInt *lidxs,*work = NULL; 24a8116848SStefano Zampini PetscSFNode *ridxs; 25a8116848SStefano Zampini PetscMPIInt rank; 26a8116848SStefano Zampini PetscInt r, p = 0, len = 0; 27a8116848SStefano Zampini PetscErrorCode ierr; 28a8116848SStefano Zampini 29a8116848SStefano Zampini PetscFunctionBegin; 30a8116848SStefano Zampini /* Create SF where leaves are input idxs and roots are owned idxs (code adapted from MatZeroRowsMapLocal_Private) */ 31a8116848SStefano Zampini ierr = MPI_Comm_rank(map->comm,&rank);CHKERRQ(ierr); 32a8116848SStefano Zampini ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr); 33a8116848SStefano Zampini for (r = 0; r < n; ++r) lidxs[r] = -1; 34a8116848SStefano Zampini ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr); 35a8116848SStefano Zampini for (r = 0; r < N; ++r) { 36a8116848SStefano Zampini const PetscInt idx = idxs[r]; 37a8116848SStefano 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); 38a8116848SStefano Zampini if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 39a8116848SStefano Zampini ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr); 40a8116848SStefano Zampini } 41a8116848SStefano Zampini ridxs[r].rank = p; 42a8116848SStefano Zampini ridxs[r].index = idxs[r] - owners[p]; 43a8116848SStefano Zampini } 44a8116848SStefano Zampini ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr); 45a8116848SStefano Zampini ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr); 46a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 47a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 48a8116848SStefano Zampini if ((gidxs && ogidxs) || ogidxs) { 49a8116848SStefano Zampini ierr = PetscMalloc1(n,&work);CHKERRQ(ierr); 50a8116848SStefano Zampini if (!gidxs) { /* if gidxs is not provided, global indices are inferred */ 51a8116848SStefano Zampini PetscInt cum = 0,start,*work2; 52a8116848SStefano Zampini ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr); 53a8116848SStefano Zampini for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++; 54a8116848SStefano Zampini ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRQ(ierr); 55a8116848SStefano Zampini start -= cum; 56a8116848SStefano Zampini cum = 0; 57a8116848SStefano Zampini for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++; 58a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 59a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 60a8116848SStefano Zampini ierr = PetscFree(work2);CHKERRQ(ierr); 61a8116848SStefano Zampini } else { 62a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)gidxs,work,MPIU_REPLACE);CHKERRQ(ierr); 63a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)gidxs,work,MPIU_REPLACE);CHKERRQ(ierr); 64a8116848SStefano Zampini } 65a8116848SStefano Zampini } 66a8116848SStefano Zampini ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 67a8116848SStefano Zampini /* Compress and put in indices */ 68a8116848SStefano Zampini for (r = 0; r < n; ++r) 69a8116848SStefano Zampini if (lidxs[r] >= 0) { 70a8116848SStefano Zampini if (work) work[len] = work[r]; 71a8116848SStefano Zampini lidxs[len++] = r; 72a8116848SStefano Zampini } 73a8116848SStefano Zampini if (on) *on = len; 74a8116848SStefano Zampini if (oidxs) *oidxs = lidxs; 75a8116848SStefano Zampini if (ogidxs) *ogidxs = work; 76a8116848SStefano Zampini PetscFunctionReturn(0); 77a8116848SStefano Zampini } 78a8116848SStefano Zampini 79a8116848SStefano Zampini #undef __FUNCT__ 80a8116848SStefano Zampini #define __FUNCT__ "MatGetSubMatrix_IS" 81a8116848SStefano Zampini static PetscErrorCode MatGetSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat) 82a8116848SStefano Zampini { 83a8116848SStefano Zampini Mat locmat,newlocmat; 84a8116848SStefano Zampini Mat_IS *newmatis; 85a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG) 86a8116848SStefano Zampini Vec rtest,ltest; 87a8116848SStefano Zampini const PetscScalar *array; 88a8116848SStefano Zampini #endif 89a8116848SStefano Zampini const PetscInt *idxs; 90a8116848SStefano Zampini PetscInt i,m,n; 91a8116848SStefano Zampini PetscErrorCode ierr; 92a8116848SStefano Zampini 93a8116848SStefano Zampini PetscFunctionBegin; 94a8116848SStefano Zampini if (scall == MAT_REUSE_MATRIX) { 95a8116848SStefano Zampini PetscBool ismatis; 96a8116848SStefano Zampini 97a8116848SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr); 98a8116848SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type"); 99a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 100a8116848SStefano Zampini if (!newmatis->getsub_ris) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS"); 101a8116848SStefano Zampini if (!newmatis->getsub_cis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS"); 102a8116848SStefano Zampini } 103a8116848SStefano Zampini /* irow and icol may not have duplicate entries */ 104a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG) 105a8116848SStefano Zampini ierr = MatCreateVecs(mat,<est,&rtest);CHKERRQ(ierr); 106a8116848SStefano Zampini ierr = ISGetLocalSize(irow,&n);CHKERRQ(ierr); 107a8116848SStefano Zampini ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 108a8116848SStefano Zampini for (i=0;i<n;i++) { 109a8116848SStefano Zampini ierr = VecSetValue(rtest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 110a8116848SStefano Zampini } 111a8116848SStefano Zampini ierr = VecAssemblyBegin(rtest);CHKERRQ(ierr); 112a8116848SStefano Zampini ierr = VecAssemblyEnd(rtest);CHKERRQ(ierr); 113a8116848SStefano Zampini ierr = VecGetLocalSize(rtest,&n);CHKERRQ(ierr); 114a8116848SStefano Zampini ierr = VecGetOwnershipRange(rtest,&m,NULL);CHKERRQ(ierr); 115a8116848SStefano Zampini ierr = VecGetArrayRead(rtest,&array);CHKERRQ(ierr); 116fd479f66SMatthew 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])); 117a8116848SStefano Zampini ierr = VecRestoreArrayRead(rtest,&array);CHKERRQ(ierr); 118a8116848SStefano Zampini ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 119a8116848SStefano Zampini ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 120a8116848SStefano Zampini ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 121a8116848SStefano Zampini for (i=0;i<n;i++) { 122a8116848SStefano Zampini ierr = VecSetValue(ltest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 123a8116848SStefano Zampini } 124a8116848SStefano Zampini ierr = VecAssemblyBegin(ltest);CHKERRQ(ierr); 125a8116848SStefano Zampini ierr = VecAssemblyEnd(ltest);CHKERRQ(ierr); 126a8116848SStefano Zampini ierr = VecGetLocalSize(ltest,&n);CHKERRQ(ierr); 127a8116848SStefano Zampini ierr = VecGetOwnershipRange(ltest,&m,NULL);CHKERRQ(ierr); 128a8116848SStefano Zampini ierr = VecGetArrayRead(ltest,&array);CHKERRQ(ierr); 129fd479f66SMatthew 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])); 130a8116848SStefano Zampini ierr = VecRestoreArrayRead(ltest,&array);CHKERRQ(ierr); 131a8116848SStefano Zampini ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 132a8116848SStefano Zampini ierr = VecDestroy(&rtest);CHKERRQ(ierr); 133a8116848SStefano Zampini ierr = VecDestroy(<est);CHKERRQ(ierr); 134a8116848SStefano Zampini #endif 135a8116848SStefano Zampini if (scall == MAT_INITIAL_MATRIX) { 136a8116848SStefano Zampini Mat_IS *matis = (Mat_IS*)mat->data; 137a8116848SStefano Zampini ISLocalToGlobalMapping rl2g; 138a8116848SStefano Zampini IS is; 139a8116848SStefano Zampini PetscInt *lidxs,*lgidxs,*newgidxs; 140*6dd40735SStefano Zampini PetscInt ll,newloc; 141a8116848SStefano Zampini MPI_Comm comm; 142a8116848SStefano Zampini 143a8116848SStefano Zampini ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 144a8116848SStefano Zampini ierr = ISGetLocalSize(irow,&m);CHKERRQ(ierr); 145a8116848SStefano Zampini ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 146a8116848SStefano Zampini ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 147a8116848SStefano Zampini ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr); 148a8116848SStefano Zampini ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 149a8116848SStefano Zampini /* communicate irow to their owners in the layout */ 150a8116848SStefano Zampini ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 151a8116848SStefano Zampini ierr = PetscLayoutMapLocal_Private(mat->rmap,m,idxs,NULL,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 152a8116848SStefano Zampini ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 153a8116848SStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 154a8116848SStefano Zampini ierr = MatISComputeSF_Private(mat);CHKERRQ(ierr); 155a8116848SStefano Zampini } 156a8116848SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,matis->sf_nroots*sizeof(PetscInt));CHKERRQ(ierr); 157a8116848SStefano Zampini for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1; 158a8116848SStefano Zampini ierr = PetscFree(lidxs);CHKERRQ(ierr); 159a8116848SStefano Zampini ierr = PetscFree(lgidxs);CHKERRQ(ierr); 160a8116848SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 161a8116848SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 162a8116848SStefano Zampini for (i=0,newloc=0;i<matis->sf_nleaves;i++) if (matis->sf_leafdata[i]) newloc++; 163a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 164a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 165a8116848SStefano Zampini for (i=0,newloc=0;i<matis->sf_nleaves;i++) 166a8116848SStefano Zampini if (matis->sf_leafdata[i]) { 167a8116848SStefano Zampini lidxs[newloc] = i; 168a8116848SStefano Zampini newgidxs[newloc++] = matis->sf_leafdata[i]-1; 169a8116848SStefano Zampini } 170a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 171a8116848SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 172a8116848SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 173a8116848SStefano Zampini /* local is to extract local submatrix */ 174a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 175a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris);CHKERRQ(ierr); 176a8116848SStefano Zampini if (mat->congruentlayouts == -1) { /* first time we compare rows and cols layouts */ 177a8116848SStefano Zampini PetscBool cong; 178a8116848SStefano Zampini ierr = PetscLayoutCompare(mat->rmap,mat->cmap,&cong);CHKERRQ(ierr); 179a8116848SStefano Zampini if (cong) mat->congruentlayouts = 1; 180a8116848SStefano Zampini else mat->congruentlayouts = 0; 181a8116848SStefano Zampini } 182a8116848SStefano Zampini if (mat->congruentlayouts && irow == icol) { 183a8116848SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g);CHKERRQ(ierr); 184a8116848SStefano Zampini ierr = PetscObjectReference((PetscObject)newmatis->getsub_ris);CHKERRQ(ierr); 185a8116848SStefano Zampini newmatis->getsub_cis = newmatis->getsub_ris; 186a8116848SStefano Zampini } else { 187a8116848SStefano Zampini ISLocalToGlobalMapping cl2g; 188a8116848SStefano Zampini 189a8116848SStefano Zampini /* communicate icol to their owners in the layout */ 190a8116848SStefano Zampini ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 191a8116848SStefano Zampini ierr = PetscLayoutMapLocal_Private(mat->cmap,n,idxs,NULL,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 192a8116848SStefano Zampini ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 193a8116848SStefano Zampini ierr = PetscMemzero(matis->csf_rootdata,matis->csf_nroots*sizeof(PetscInt));CHKERRQ(ierr); 194a8116848SStefano Zampini for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1; 195a8116848SStefano Zampini ierr = PetscFree(lidxs);CHKERRQ(ierr); 196a8116848SStefano Zampini ierr = PetscFree(lgidxs);CHKERRQ(ierr); 197a8116848SStefano Zampini ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 198a8116848SStefano Zampini ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 199a8116848SStefano Zampini for (i=0,newloc=0;i<matis->csf_nleaves;i++) if (matis->csf_leafdata[i]) newloc++; 200a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 201a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 202a8116848SStefano Zampini for (i=0,newloc=0;i<matis->csf_nleaves;i++) 203a8116848SStefano Zampini if (matis->csf_leafdata[i]) { 204a8116848SStefano Zampini lidxs[newloc] = i; 205a8116848SStefano Zampini newgidxs[newloc++] = matis->csf_leafdata[i]-1; 206a8116848SStefano Zampini } 207a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 208a8116848SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 209a8116848SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 210a8116848SStefano Zampini /* local is to extract local submatrix */ 211a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis);CHKERRQ(ierr); 212a8116848SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr); 213a8116848SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 214a8116848SStefano Zampini } 215a8116848SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 216a8116848SStefano Zampini } else { 217a8116848SStefano Zampini ierr = MatISGetLocalMat(*newmat,&newlocmat);CHKERRQ(ierr); 218a8116848SStefano Zampini } 219a8116848SStefano Zampini ierr = MatISGetLocalMat(mat,&locmat);CHKERRQ(ierr); 220a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 221a8116848SStefano Zampini ierr = MatGetSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat);CHKERRQ(ierr); 222a8116848SStefano Zampini if (scall == MAT_INITIAL_MATRIX) { 223a8116848SStefano Zampini ierr = MatISSetLocalMat(*newmat,newlocmat);CHKERRQ(ierr); 224a8116848SStefano Zampini ierr = MatDestroy(&newlocmat);CHKERRQ(ierr); 225a8116848SStefano Zampini } 226a8116848SStefano Zampini ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 227a8116848SStefano Zampini ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 228a8116848SStefano Zampini PetscFunctionReturn(0); 229a8116848SStefano Zampini } 230a8116848SStefano Zampini 231a8116848SStefano Zampini #undef __FUNCT__ 232659959c5SStefano Zampini #define __FUNCT__ "MatGetLocalSubMatrix_IS" 233a8116848SStefano Zampini static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat) 234659959c5SStefano Zampini { 235659959c5SStefano Zampini Mat lA,lsubmat; 236659959c5SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 237659959c5SStefano Zampini IS is,isn; 238659959c5SStefano Zampini const PetscInt *rg,*rl; 239659959c5SStefano Zampini #if defined(PETSC_USE_DEBUG) 240659959c5SStefano Zampini PetscInt nrg; 241659959c5SStefano Zampini #endif 242659959c5SStefano Zampini PetscInt N,M,nrl,i,*idxs; 243659959c5SStefano Zampini PetscErrorCode ierr; 244659959c5SStefano Zampini 245659959c5SStefano Zampini PetscFunctionBegin; 246659959c5SStefano Zampini /* extract local submatrix (it will complain if row and col exceed the local sizes) */ 247659959c5SStefano Zampini ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr); 248659959c5SStefano Zampini ierr = MatGetSubMatrix(lA,row,col,MAT_INITIAL_MATRIX,&lsubmat);CHKERRQ(ierr); 249659959c5SStefano Zampini /* compute new l2g map for rows */ 250659959c5SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 251659959c5SStefano Zampini ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr); 252659959c5SStefano Zampini ierr = ISGetIndices(row,&rl);CHKERRQ(ierr); 253659959c5SStefano Zampini #if defined(PETSC_USE_DEBUG) 254659959c5SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr); 255659959c5SStefano 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); 256659959c5SStefano Zampini #endif 257659959c5SStefano Zampini ierr = PetscMalloc1(nrl,&idxs);CHKERRQ(ierr); 258659959c5SStefano Zampini for (i=0;i<nrl;i++) idxs[i] = rg[rl[i]]; 259659959c5SStefano Zampini ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr); 260659959c5SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 261659959c5SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrl,idxs,PETSC_COPY_VALUES,&is);CHKERRQ(ierr); 262659959c5SStefano Zampini ierr = PetscFree(idxs);CHKERRQ(ierr); 263659959c5SStefano Zampini ierr = ISRenumber(is,NULL,&M,&isn);CHKERRQ(ierr); 264659959c5SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 265659959c5SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(isn,&rl2g);CHKERRQ(ierr); 266659959c5SStefano Zampini ierr = ISDestroy(&isn);CHKERRQ(ierr); 267659959c5SStefano Zampini /* compute new l2g map for columns */ 268659959c5SStefano Zampini if (col != row || A->rmap->mapping != A->cmap->mapping) { 269659959c5SStefano Zampini const PetscInt *cg,*cl; 270659959c5SStefano Zampini #if defined(PETSC_USE_DEBUG) 271659959c5SStefano Zampini PetscInt ncg; 272659959c5SStefano Zampini #endif 273659959c5SStefano Zampini PetscInt ncl; 274659959c5SStefano Zampini 275659959c5SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 276659959c5SStefano Zampini ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr); 277659959c5SStefano Zampini ierr = ISGetIndices(col,&cl);CHKERRQ(ierr); 278659959c5SStefano Zampini #if defined(PETSC_USE_DEBUG) 279659959c5SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr); 280659959c5SStefano 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); 281659959c5SStefano Zampini #endif 282659959c5SStefano Zampini ierr = PetscMalloc1(ncl,&idxs);CHKERRQ(ierr); 283659959c5SStefano Zampini for (i=0;i<ncl;i++) idxs[i] = cg[cl[i]]; 284659959c5SStefano Zampini ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr); 285659959c5SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 286659959c5SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncl,idxs,PETSC_COPY_VALUES,&is);CHKERRQ(ierr); 287659959c5SStefano Zampini ierr = PetscFree(idxs);CHKERRQ(ierr); 288659959c5SStefano Zampini ierr = ISRenumber(is,NULL,&N,&isn);CHKERRQ(ierr); 289659959c5SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 290659959c5SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(isn,&cl2g);CHKERRQ(ierr); 291659959c5SStefano Zampini ierr = ISDestroy(&isn);CHKERRQ(ierr); 292659959c5SStefano Zampini } else { 2930a40dfa3SStefano Zampini ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr); 294659959c5SStefano Zampini cl2g = rl2g; 295659959c5SStefano Zampini N = M; 296659959c5SStefano Zampini } 297659959c5SStefano Zampini /* create the MATIS submatrix */ 298659959c5SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr); 299659959c5SStefano Zampini ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 300659959c5SStefano Zampini ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr); 301659959c5SStefano Zampini ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr); 302659959c5SStefano Zampini ierr = MatISSetLocalMat(*submat,lsubmat);CHKERRQ(ierr); 303659959c5SStefano Zampini ierr = MatDestroy(&lsubmat);CHKERRQ(ierr); 304659959c5SStefano Zampini ierr = MatSetUp(*submat);CHKERRQ(ierr); 305659959c5SStefano Zampini ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 306659959c5SStefano Zampini ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 307659959c5SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 308659959c5SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 309659959c5SStefano Zampini PetscFunctionReturn(0); 310659959c5SStefano Zampini } 311659959c5SStefano Zampini 312659959c5SStefano Zampini #undef __FUNCT__ 3132b404112SStefano Zampini #define __FUNCT__ "MatCopy_IS" 314a8116848SStefano Zampini static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str) 3152b404112SStefano Zampini { 3162b404112SStefano Zampini Mat_IS *a = (Mat_IS*)A->data,*b; 3172b404112SStefano Zampini PetscBool ismatis; 3182b404112SStefano Zampini PetscErrorCode ierr; 3192b404112SStefano Zampini 3202b404112SStefano Zampini PetscFunctionBegin; 3212b404112SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr); 3222b404112SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented"); 3232b404112SStefano Zampini b = (Mat_IS*)B->data; 3242b404112SStefano Zampini ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr); 3252b404112SStefano Zampini PetscFunctionReturn(0); 3262b404112SStefano Zampini } 3272b404112SStefano Zampini 3282b404112SStefano Zampini #undef __FUNCT__ 3296bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS" 330a8116848SStefano Zampini static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool *missing,PetscInt *d) 3316bd84002SStefano Zampini { 332527b2640SStefano Zampini Vec v; 333527b2640SStefano Zampini const PetscScalar *array; 334527b2640SStefano Zampini PetscInt i,n; 3356bd84002SStefano Zampini PetscErrorCode ierr; 3366bd84002SStefano Zampini 3376bd84002SStefano Zampini PetscFunctionBegin; 338527b2640SStefano Zampini *missing = PETSC_FALSE; 339527b2640SStefano Zampini ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr); 340527b2640SStefano Zampini ierr = MatGetDiagonal(A,v);CHKERRQ(ierr); 341527b2640SStefano Zampini ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 342527b2640SStefano Zampini ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr); 343527b2640SStefano Zampini for (i=0;i<n;i++) if (array[i] == 0.) break; 344527b2640SStefano Zampini ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr); 345527b2640SStefano Zampini ierr = VecDestroy(&v);CHKERRQ(ierr); 346527b2640SStefano Zampini if (i != n) *missing = PETSC_TRUE; 347527b2640SStefano Zampini if (d) { 348527b2640SStefano Zampini *d = -1; 349527b2640SStefano Zampini if (*missing) { 350527b2640SStefano Zampini PetscInt rstart; 351527b2640SStefano Zampini ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 352527b2640SStefano Zampini *d = i+rstart; 353527b2640SStefano Zampini } 354527b2640SStefano Zampini } 3556bd84002SStefano Zampini PetscFunctionReturn(0); 3566bd84002SStefano Zampini } 3576bd84002SStefano Zampini 3586bd84002SStefano Zampini #undef __FUNCT__ 35928f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private" 360a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B) 36128f4e0baSStefano Zampini { 36228f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 36328f4e0baSStefano Zampini const PetscInt *gidxs; 36428f4e0baSStefano Zampini PetscErrorCode ierr; 36528f4e0baSStefano Zampini 36628f4e0baSStefano Zampini PetscFunctionBegin; 36728f4e0baSStefano Zampini ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr); 36828f4e0baSStefano Zampini ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr); 36928f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 3703bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 37128f4e0baSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 3723bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 37328f4e0baSStefano Zampini ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 374a8116848SStefano Zampini if (B->rmap->mapping != B->cmap->mapping) { /* setup SF for columns */ 375a8116848SStefano Zampini ierr = MatGetSize(matis->A,NULL,&matis->csf_nleaves);CHKERRQ(ierr); 376a8116848SStefano Zampini ierr = MatGetLocalSize(B,NULL,&matis->csf_nroots);CHKERRQ(ierr); 377a8116848SStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf);CHKERRQ(ierr); 378a8116848SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 379a8116848SStefano Zampini ierr = PetscSFSetGraphLayout(matis->csf,B->cmap,matis->csf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 380a8116848SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 381a8116848SStefano Zampini ierr = PetscMalloc2(matis->csf_nroots,&matis->csf_rootdata,matis->csf_nleaves,&matis->csf_leafdata);CHKERRQ(ierr); 382a8116848SStefano Zampini } else { 383a8116848SStefano Zampini matis->csf = matis->sf; 384a8116848SStefano Zampini matis->csf_nleaves = matis->sf_nleaves; 385a8116848SStefano Zampini matis->csf_nroots = matis->sf_nroots; 386a8116848SStefano Zampini matis->csf_leafdata = matis->sf_leafdata; 387a8116848SStefano Zampini matis->csf_rootdata = matis->sf_rootdata; 388a8116848SStefano Zampini } 38928f4e0baSStefano Zampini PetscFunctionReturn(0); 39028f4e0baSStefano Zampini } 3912e1947a5SStefano Zampini 3922e1947a5SStefano Zampini #undef __FUNCT__ 3932e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 394eb82efa4SStefano Zampini /*@ 395a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 396a88811baSStefano Zampini 397a88811baSStefano Zampini Collective on MPI_Comm 398a88811baSStefano Zampini 399a88811baSStefano Zampini Input Parameters: 400a88811baSStefano Zampini + B - the matrix 401a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 402a88811baSStefano Zampini (same value is used for all local rows) 403a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 404a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 405a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 406a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 407a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 408a88811baSStefano Zampini the diagonal entry even if it is zero. 409a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 410a88811baSStefano Zampini submatrix (same value is used for all local rows). 411a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 412a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 413a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 414a88811baSStefano Zampini structure. The size of this array is equal to the number 415a88811baSStefano Zampini of local rows, i.e 'm'. 416a88811baSStefano Zampini 417a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 418a88811baSStefano Zampini 419a88811baSStefano Zampini Level: intermediate 420a88811baSStefano Zampini 421a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 422a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 423a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 424a88811baSStefano Zampini 425a88811baSStefano Zampini .keywords: matrix 426a88811baSStefano Zampini 4273c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS 428a88811baSStefano Zampini @*/ 4292e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 4302e1947a5SStefano Zampini { 4312e1947a5SStefano Zampini PetscErrorCode ierr; 4322e1947a5SStefano Zampini 4332e1947a5SStefano Zampini PetscFunctionBegin; 4342e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 4352e1947a5SStefano Zampini PetscValidType(B,1); 4362e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 4372e1947a5SStefano Zampini PetscFunctionReturn(0); 4382e1947a5SStefano Zampini } 4392e1947a5SStefano Zampini 4402e1947a5SStefano Zampini #undef __FUNCT__ 4412e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 4427230de76SStefano Zampini static PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 4432e1947a5SStefano Zampini { 4442e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 44528f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 4462e1947a5SStefano Zampini PetscErrorCode ierr; 4472e1947a5SStefano Zampini 4482e1947a5SStefano Zampini PetscFunctionBegin; 4496c4ed002SBarry Smith if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 45028f4e0baSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 45128f4e0baSStefano Zampini ierr = MatISComputeSF_Private(B);CHKERRQ(ierr); 45228f4e0baSStefano Zampini } 4532e1947a5SStefano Zampini if (!d_nnz) { 45428f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz; 4552e1947a5SStefano Zampini } else { 45628f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 4572e1947a5SStefano Zampini } 4582e1947a5SStefano Zampini if (!o_nnz) { 45928f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz; 4602e1947a5SStefano Zampini } else { 46128f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 4622e1947a5SStefano Zampini } 46328f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 46428f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 46528f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 46628f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 46728f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves;i++) { 46828f4e0baSStefano Zampini matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 4692e1947a5SStefano Zampini } 47028f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 47128f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 47228f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 4732e1947a5SStefano Zampini } 47428f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 47528f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 47628f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 4772e1947a5SStefano Zampini } 47828f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 4792e1947a5SStefano Zampini PetscFunctionReturn(0); 4802e1947a5SStefano Zampini } 481b4319ba4SBarry Smith 482b4319ba4SBarry Smith #undef __FUNCT__ 4833927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private" 4843927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 4853927de2eSStefano Zampini { 4863927de2eSStefano Zampini Mat_IS *matis = (Mat_IS*)(A->data); 4873927de2eSStefano Zampini PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 488ecf5a873SStefano Zampini const PetscInt *global_indices_r,*global_indices_c; 4893927de2eSStefano Zampini PetscInt i,j,bs,rows,cols; 4903927de2eSStefano Zampini PetscInt lrows,lcols; 4913927de2eSStefano Zampini PetscInt local_rows,local_cols; 4923927de2eSStefano Zampini PetscMPIInt nsubdomains; 4933927de2eSStefano Zampini PetscBool isdense,issbaij; 4943927de2eSStefano Zampini PetscErrorCode ierr; 4953927de2eSStefano Zampini 4963927de2eSStefano Zampini PetscFunctionBegin; 4973927de2eSStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr); 4983927de2eSStefano Zampini ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr); 4993927de2eSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 5003927de2eSStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 5013927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 5023927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 503ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 504ecf5a873SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 5057230de76SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 506ecf5a873SStefano Zampini } else { 507ecf5a873SStefano Zampini global_indices_c = global_indices_r; 508ecf5a873SStefano Zampini } 509ecf5a873SStefano Zampini 5103927de2eSStefano Zampini if (issbaij) { 5113927de2eSStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 5123927de2eSStefano Zampini } 5133927de2eSStefano Zampini /* 514ecf5a873SStefano Zampini An SF reduce is needed to sum up properly on shared rows. 5153927de2eSStefano Zampini Note that generally preallocation is not exact, since it overestimates nonzeros 5163927de2eSStefano Zampini */ 5173927de2eSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 5183927de2eSStefano Zampini ierr = MatISComputeSF_Private(A);CHKERRQ(ierr); 5193927de2eSStefano Zampini } 5203927de2eSStefano Zampini ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr); 5213927de2eSStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr); 5223927de2eSStefano Zampini /* All processes need to compute entire row ownership */ 5233927de2eSStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 5243927de2eSStefano Zampini ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 5253927de2eSStefano Zampini for (i=0;i<nsubdomains;i++) { 5263927de2eSStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 5273927de2eSStefano Zampini row_ownership[j] = i; 5283927de2eSStefano Zampini } 5293927de2eSStefano Zampini } 5307230de76SStefano Zampini ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 5313927de2eSStefano Zampini 5323927de2eSStefano Zampini /* 5333927de2eSStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 5343927de2eSStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 5353927de2eSStefano Zampini */ 5363927de2eSStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 5373927de2eSStefano Zampini /* preallocation as a MATAIJ */ 5383927de2eSStefano Zampini if (isdense) { /* special case for dense local matrices */ 5393927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 540ecf5a873SStefano Zampini PetscInt index_row = global_indices_r[i]; 5413927de2eSStefano Zampini for (j=i;j<local_rows;j++) { 5423927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 543ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[j]; 5443927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 5453927de2eSStefano Zampini my_dnz[i] += 1; 5463927de2eSStefano Zampini } else { /* offdiag block */ 5473927de2eSStefano Zampini my_onz[i] += 1; 5483927de2eSStefano Zampini } 5493927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 5503927de2eSStefano Zampini if (i != j) { 5513927de2eSStefano Zampini owner = row_ownership[index_col]; 5523927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 5533927de2eSStefano Zampini my_dnz[j] += 1; 5543927de2eSStefano Zampini } else { 5553927de2eSStefano Zampini my_onz[j] += 1; 5563927de2eSStefano Zampini } 5573927de2eSStefano Zampini } 5583927de2eSStefano Zampini } 5593927de2eSStefano Zampini } 560ecf5a873SStefano Zampini } else { /* TODO: this could be optimized using MatGetRowIJ */ 5613927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 5623927de2eSStefano Zampini const PetscInt *cols; 563ecf5a873SStefano Zampini PetscInt ncols,index_row = global_indices_r[i]; 5643927de2eSStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 5653927de2eSStefano Zampini for (j=0;j<ncols;j++) { 5663927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 567ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[cols[j]]; 5683927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 5693927de2eSStefano Zampini my_dnz[i] += 1; 5703927de2eSStefano Zampini } else { /* offdiag block */ 5713927de2eSStefano Zampini my_onz[i] += 1; 5723927de2eSStefano Zampini } 5733927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 574d9a9e74cSStefano Zampini if (issbaij && index_col != index_row) { 5753927de2eSStefano Zampini owner = row_ownership[index_col]; 5763927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 577d9a9e74cSStefano Zampini my_dnz[cols[j]] += 1; 5783927de2eSStefano Zampini } else { 579d9a9e74cSStefano Zampini my_onz[cols[j]] += 1; 5803927de2eSStefano Zampini } 5813927de2eSStefano Zampini } 5823927de2eSStefano Zampini } 5833927de2eSStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 5843927de2eSStefano Zampini } 5853927de2eSStefano Zampini } 586ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 587ecf5a873SStefano Zampini if (global_indices_c != global_indices_r) { 5887230de76SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 589ecf5a873SStefano Zampini } 5903927de2eSStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 591ecf5a873SStefano Zampini 592ecf5a873SStefano Zampini /* Reduce my_dnz and my_onz */ 5933927de2eSStefano Zampini if (maxreduce) { 5943927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 5953927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 5963927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 5973927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 5983927de2eSStefano Zampini } else { 5993927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 6003927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 6013927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 6023927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 6033927de2eSStefano Zampini } 6043927de2eSStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 6053927de2eSStefano Zampini 6063927de2eSStefano Zampini /* Resize preallocation if overestimated */ 6073927de2eSStefano Zampini for (i=0;i<lrows;i++) { 6083927de2eSStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 6093927de2eSStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 6103927de2eSStefano Zampini } 6113927de2eSStefano Zampini /* set preallocation */ 6123927de2eSStefano Zampini ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 6133927de2eSStefano Zampini for (i=0;i<lrows/bs;i++) { 6143927de2eSStefano Zampini dnz[i] = dnz[i*bs]/bs; 6153927de2eSStefano Zampini onz[i] = onz[i*bs]/bs; 6163927de2eSStefano Zampini } 6173927de2eSStefano Zampini ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 6183927de2eSStefano Zampini ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 6193927de2eSStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 6203927de2eSStefano Zampini if (issbaij) { 6213927de2eSStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 6223927de2eSStefano Zampini } 6233927de2eSStefano Zampini PetscFunctionReturn(0); 6243927de2eSStefano Zampini } 6253927de2eSStefano Zampini 6263927de2eSStefano Zampini #undef __FUNCT__ 627b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 6287230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 629b7ce53b6SStefano Zampini { 630b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 631d9a9e74cSStefano Zampini Mat local_mat; 632b7ce53b6SStefano Zampini /* info on mat */ 6333cfa4ea4SStefano Zampini PetscInt bs,rows,cols,lrows,lcols; 634b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 635686e3a49SStefano Zampini PetscBool isdense,issbaij,isseqaij; 6367c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 637b7ce53b6SStefano Zampini /* values insertion */ 638b7ce53b6SStefano Zampini PetscScalar *array; 639b7ce53b6SStefano Zampini /* work */ 640b7ce53b6SStefano Zampini PetscErrorCode ierr; 641b7ce53b6SStefano Zampini 642b7ce53b6SStefano Zampini PetscFunctionBegin; 643b7ce53b6SStefano Zampini /* get info from mat */ 6447c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 6457c03b4e8SStefano Zampini if (nsubdomains == 1) { 6467c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 6477c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 6487c03b4e8SStefano Zampini } else { 6497c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 6507c03b4e8SStefano Zampini } 6517c03b4e8SStefano Zampini PetscFunctionReturn(0); 6527c03b4e8SStefano Zampini } 653b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 654b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 6553cfa4ea4SStefano Zampini ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 656b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 657b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 658686e3a49SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 659b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 660b7ce53b6SStefano Zampini 661b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 662b7ce53b6SStefano Zampini MatType new_mat_type; 6633927de2eSStefano Zampini PetscBool issbaij_red; 664b7ce53b6SStefano Zampini 665b7ce53b6SStefano Zampini /* determining new matrix type */ 666b2566f29SBarry Smith ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 667b7ce53b6SStefano Zampini if (issbaij_red) { 668b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 669b7ce53b6SStefano Zampini } else { 670b7ce53b6SStefano Zampini if (bs>1) { 671b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 672b7ce53b6SStefano Zampini } else { 673b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 674b7ce53b6SStefano Zampini } 675b7ce53b6SStefano Zampini } 676b7ce53b6SStefano Zampini 6773927de2eSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 6783cfa4ea4SStefano Zampini ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 6793927de2eSStefano Zampini ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 6803927de2eSStefano Zampini ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr); 6813927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 682b7ce53b6SStefano Zampini } else { 6833cfa4ea4SStefano Zampini PetscInt mbs,mrows,mcols,mlrows,mlcols; 684b7ce53b6SStefano Zampini /* some checks */ 685b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 686b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 6873cfa4ea4SStefano Zampini ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 6886c4ed002SBarry Smith if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 6896c4ed002SBarry Smith if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 6906c4ed002SBarry Smith if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 6916c4ed002SBarry Smith if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 6926c4ed002SBarry Smith if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 693b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 694b7ce53b6SStefano Zampini } 695d9a9e74cSStefano Zampini 696d9a9e74cSStefano Zampini if (issbaij) { 697d9a9e74cSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 698d9a9e74cSStefano Zampini } else { 699d9a9e74cSStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 700d9a9e74cSStefano Zampini local_mat = matis->A; 701d9a9e74cSStefano Zampini } 702686e3a49SStefano Zampini 703b7ce53b6SStefano Zampini /* Set values */ 704ecf5a873SStefano Zampini ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 705b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 706ecf5a873SStefano Zampini PetscInt i,*dummy_rows,*dummy_cols; 707ecf5a873SStefano Zampini 708ecf5a873SStefano Zampini if (local_rows != local_cols) { 709ecf5a873SStefano Zampini ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 710ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 711ecf5a873SStefano Zampini for (i=0;i<local_cols;i++) dummy_cols[i] = i; 712ecf5a873SStefano Zampini } else { 713ecf5a873SStefano Zampini ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 714ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 715ecf5a873SStefano Zampini dummy_cols = dummy_rows; 716ecf5a873SStefano Zampini } 717b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 718d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 719ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 720d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 721ecf5a873SStefano Zampini if (dummy_rows != dummy_cols) { 722ecf5a873SStefano Zampini ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 723ecf5a873SStefano Zampini } else { 724ecf5a873SStefano Zampini ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 725ecf5a873SStefano Zampini } 726686e3a49SStefano Zampini } else if (isseqaij) { 727ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 728686e3a49SStefano Zampini PetscBool done; 729686e3a49SStefano Zampini 730d9a9e74cSStefano Zampini ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 7316c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 732d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 733686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 734ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 735686e3a49SStefano Zampini } 736d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 7376c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 738d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 739686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 740ecf5a873SStefano Zampini PetscInt i; 741c0962df8SStefano Zampini 742686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 743686e3a49SStefano Zampini PetscInt j; 744ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 745686e3a49SStefano Zampini 746ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 747ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 748ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 749686e3a49SStefano Zampini } 750b7ce53b6SStefano Zampini } 751b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 752d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 753b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 754b7ce53b6SStefano Zampini if (isdense) { 755b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 756b7ce53b6SStefano Zampini } 757b7ce53b6SStefano Zampini PetscFunctionReturn(0); 758b7ce53b6SStefano Zampini } 759b7ce53b6SStefano Zampini 760b7ce53b6SStefano Zampini #undef __FUNCT__ 761b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 762b7ce53b6SStefano Zampini /*@ 763b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 764b7ce53b6SStefano Zampini 765b7ce53b6SStefano Zampini Input Parameter: 766b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 767b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 768b7ce53b6SStefano Zampini 769b7ce53b6SStefano Zampini Output Parameter: 770b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 771b7ce53b6SStefano Zampini 772b7ce53b6SStefano Zampini Level: developer 773b7ce53b6SStefano Zampini 774eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 775b7ce53b6SStefano Zampini 776b7ce53b6SStefano Zampini .seealso: MATIS 777b7ce53b6SStefano Zampini @*/ 778b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 779b7ce53b6SStefano Zampini { 780b7ce53b6SStefano Zampini PetscErrorCode ierr; 781b7ce53b6SStefano Zampini 782b7ce53b6SStefano Zampini PetscFunctionBegin; 783b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 784b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 785b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 786b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 787b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 788b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 7896c4ed002SBarry Smith if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 790b7ce53b6SStefano Zampini } 791b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 792b7ce53b6SStefano Zampini PetscFunctionReturn(0); 793b7ce53b6SStefano Zampini } 794b7ce53b6SStefano Zampini 795b7ce53b6SStefano Zampini #undef __FUNCT__ 796ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 797ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 798ad6194a2SStefano Zampini { 799ad6194a2SStefano Zampini PetscErrorCode ierr; 800ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 801ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 802ad6194a2SStefano Zampini Mat B,localmat; 803ad6194a2SStefano Zampini 804ad6194a2SStefano Zampini PetscFunctionBegin; 805ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 806ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 807ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 808e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 809ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 810ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 811b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 812ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 813ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 814ad6194a2SStefano Zampini *newmat = B; 815ad6194a2SStefano Zampini PetscFunctionReturn(0); 816ad6194a2SStefano Zampini } 817ad6194a2SStefano Zampini 818ad6194a2SStefano Zampini #undef __FUNCT__ 81969796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 820a8116848SStefano Zampini static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 82169796d55SStefano Zampini { 82269796d55SStefano Zampini PetscErrorCode ierr; 82369796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 82469796d55SStefano Zampini PetscBool local_sym; 82569796d55SStefano Zampini 82669796d55SStefano Zampini PetscFunctionBegin; 82769796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 828b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 82969796d55SStefano Zampini PetscFunctionReturn(0); 83069796d55SStefano Zampini } 83169796d55SStefano Zampini 83269796d55SStefano Zampini #undef __FUNCT__ 83369796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 834a8116848SStefano Zampini static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 83569796d55SStefano Zampini { 83669796d55SStefano Zampini PetscErrorCode ierr; 83769796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 83869796d55SStefano Zampini PetscBool local_sym; 83969796d55SStefano Zampini 84069796d55SStefano Zampini PetscFunctionBegin; 84169796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 842b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 84369796d55SStefano Zampini PetscFunctionReturn(0); 84469796d55SStefano Zampini } 84569796d55SStefano Zampini 84669796d55SStefano Zampini #undef __FUNCT__ 847b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 848a8116848SStefano Zampini static PetscErrorCode MatDestroy_IS(Mat A) 849b4319ba4SBarry Smith { 850dfbe8321SBarry Smith PetscErrorCode ierr; 851b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 852b4319ba4SBarry Smith 853b4319ba4SBarry Smith PetscFunctionBegin; 8546bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 855e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 856e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 8576bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 8586bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 859a8116848SStefano Zampini ierr = ISDestroy(&b->getsub_ris);CHKERRQ(ierr); 860a8116848SStefano Zampini ierr = ISDestroy(&b->getsub_cis);CHKERRQ(ierr); 861a8116848SStefano Zampini if (b->sf != b->csf) { 862a8116848SStefano Zampini ierr = PetscSFDestroy(&b->csf);CHKERRQ(ierr); 863a8116848SStefano Zampini ierr = PetscFree2(b->csf_rootdata,b->csf_leafdata);CHKERRQ(ierr); 864a8116848SStefano Zampini } 86528f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 86628f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 867bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 868dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 869bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 870b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 871b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 8722e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 873b4319ba4SBarry Smith PetscFunctionReturn(0); 874b4319ba4SBarry Smith } 875b4319ba4SBarry Smith 876b4319ba4SBarry Smith #undef __FUNCT__ 877b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 878a8116848SStefano Zampini static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 879b4319ba4SBarry Smith { 880dfbe8321SBarry Smith PetscErrorCode ierr; 881b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 882b4319ba4SBarry Smith PetscScalar zero = 0.0; 883b4319ba4SBarry Smith 884b4319ba4SBarry Smith PetscFunctionBegin; 885b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 886e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 887e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 888b4319ba4SBarry Smith 889b4319ba4SBarry Smith /* multiply the local matrix */ 890b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 891b4319ba4SBarry Smith 892b4319ba4SBarry Smith /* scatter product back into global memory */ 8932dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 894e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 895e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 896b4319ba4SBarry Smith PetscFunctionReturn(0); 897b4319ba4SBarry Smith } 898b4319ba4SBarry Smith 899b4319ba4SBarry Smith #undef __FUNCT__ 9002e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 901a8116848SStefano Zampini static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 9022e74eeadSLisandro Dalcin { 903650997f4SStefano Zampini Vec temp_vec; 9042e74eeadSLisandro Dalcin PetscErrorCode ierr; 9052e74eeadSLisandro Dalcin 9062e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 907650997f4SStefano Zampini if (v3 != v2) { 908650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 909650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 910650997f4SStefano Zampini } else { 911650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 912650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 913650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 914650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 915650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 916650997f4SStefano Zampini } 9172e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9182e74eeadSLisandro Dalcin } 9192e74eeadSLisandro Dalcin 9202e74eeadSLisandro Dalcin #undef __FUNCT__ 9212e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 922a8116848SStefano Zampini static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 9232e74eeadSLisandro Dalcin { 9242e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9252e74eeadSLisandro Dalcin PetscErrorCode ierr; 9262e74eeadSLisandro Dalcin 927e176bc59SStefano Zampini PetscFunctionBegin; 9282e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 929e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 930e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 9312e74eeadSLisandro Dalcin 9322e74eeadSLisandro Dalcin /* multiply the local matrix */ 933e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 9342e74eeadSLisandro Dalcin 9352e74eeadSLisandro Dalcin /* scatter product back into global vector */ 936e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 937e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 938e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9392e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9402e74eeadSLisandro Dalcin } 9412e74eeadSLisandro Dalcin 9422e74eeadSLisandro Dalcin #undef __FUNCT__ 9432e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 944a8116848SStefano Zampini static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 9452e74eeadSLisandro Dalcin { 946650997f4SStefano Zampini Vec temp_vec; 9472e74eeadSLisandro Dalcin PetscErrorCode ierr; 9482e74eeadSLisandro Dalcin 9492e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 950650997f4SStefano Zampini if (v3 != v2) { 951650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 952650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 953650997f4SStefano Zampini } else { 954650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 955650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 956650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 957650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 958650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 959650997f4SStefano Zampini } 9602e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9612e74eeadSLisandro Dalcin } 9622e74eeadSLisandro Dalcin 9632e74eeadSLisandro Dalcin #undef __FUNCT__ 964b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 965a8116848SStefano Zampini static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 966b4319ba4SBarry Smith { 967b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 968dfbe8321SBarry Smith PetscErrorCode ierr; 969b4319ba4SBarry Smith PetscViewer sviewer; 970b4319ba4SBarry Smith 971b4319ba4SBarry Smith PetscFunctionBegin; 9723f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 973b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 9743f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 9756e520ac8SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 976b4319ba4SBarry Smith PetscFunctionReturn(0); 977b4319ba4SBarry Smith } 978b4319ba4SBarry Smith 979b4319ba4SBarry Smith #undef __FUNCT__ 980b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 981a8116848SStefano Zampini static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 982b4319ba4SBarry Smith { 983dfbe8321SBarry Smith PetscErrorCode ierr; 984e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 985b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 986b4319ba4SBarry Smith IS from,to; 987e176bc59SStefano Zampini Vec cglobal,rglobal; 988b4319ba4SBarry Smith 989b4319ba4SBarry Smith PetscFunctionBegin; 990784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 991e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 9923bbff08aSStefano Zampini /* Destroy any previous data */ 99370cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 99470cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 995e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 996e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 9971c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 99828f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 99928f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 10003bbff08aSStefano Zampini 10013bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 1002fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1003fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1004fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 1005fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 1006b4319ba4SBarry Smith 1007b4319ba4SBarry Smith /* Create the local matrix A */ 1008e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 1009e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 1010e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 1011e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 1012f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 1013e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 1014e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 1015e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 1016ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 1017ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 1018b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 1019b4319ba4SBarry Smith 1020b4319ba4SBarry Smith /* Create the local work vectors */ 10213bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 1022b4319ba4SBarry Smith 1023e176bc59SStefano Zampini /* setup the global to local scatters */ 1024e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 1025e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 1026784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 1027e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 1028e176bc59SStefano Zampini if (rmapping != cmapping) { 1029e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 1030e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 1031e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 1032e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 1033e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 1034e176bc59SStefano Zampini } else { 1035e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 1036e176bc59SStefano Zampini is->cctx = is->rctx; 1037e176bc59SStefano Zampini } 1038e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 1039e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 10406bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 10416bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 104248ff6bf3SStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 1043b4319ba4SBarry Smith PetscFunctionReturn(0); 1044b4319ba4SBarry Smith } 1045b4319ba4SBarry Smith 10462e74eeadSLisandro Dalcin #undef __FUNCT__ 10472e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 1048a8116848SStefano Zampini static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 10492e74eeadSLisandro Dalcin { 10502e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 10512e74eeadSLisandro Dalcin PetscErrorCode ierr; 105297563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 105397563a80SStefano Zampini PetscInt i,zm,zn; 105497563a80SStefano Zampini #endif 105597563a80SStefano Zampini PetscInt rows_l[2048],cols_l[2048]; 10562e74eeadSLisandro Dalcin 10572e74eeadSLisandro Dalcin PetscFunctionBegin; 10582e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 1059f23aa3ddSBarry Smith if (m > 2048 || n > 2048) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n); 106097563a80SStefano Zampini /* count negative indices */ 106197563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 106297563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 10632e74eeadSLisandro Dalcin #endif 106497563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 106597563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 106697563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 106797563a80SStefano Zampini /* count negative indices (should be the same as before) */ 106897563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 106997563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 107097563a80SStefano 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"); 107197563a80SStefano 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"); 107297563a80SStefano Zampini #endif 10732e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 10742e74eeadSLisandro Dalcin PetscFunctionReturn(0); 10752e74eeadSLisandro Dalcin } 10762e74eeadSLisandro Dalcin 1077b4319ba4SBarry Smith #undef __FUNCT__ 107897563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS" 1079a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 108097563a80SStefano Zampini { 108197563a80SStefano Zampini Mat_IS *is = (Mat_IS*)mat->data; 108297563a80SStefano Zampini PetscErrorCode ierr; 108397563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 108497563a80SStefano Zampini PetscInt i,zm,zn; 108597563a80SStefano Zampini #endif 108697563a80SStefano Zampini PetscInt rows_l[2048],cols_l[2048]; 108797563a80SStefano Zampini 108897563a80SStefano Zampini PetscFunctionBegin; 108997563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 109097563a80SStefano Zampini if (m > 2048 || n > 2048) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n); 109197563a80SStefano Zampini /* count negative indices */ 109297563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 109397563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 109497563a80SStefano Zampini #endif 109597563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 109697563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 109797563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 109897563a80SStefano Zampini /* count negative indices (should be the same as before) */ 109997563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 110097563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 110197563a80SStefano 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"); 110297563a80SStefano 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"); 110397563a80SStefano Zampini #endif 110497563a80SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 110597563a80SStefano Zampini PetscFunctionReturn(0); 110697563a80SStefano Zampini } 110797563a80SStefano Zampini 110897563a80SStefano Zampini #undef __FUNCT__ 1109b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 1110a8116848SStefano Zampini static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1111b4319ba4SBarry Smith { 1112dfbe8321SBarry Smith PetscErrorCode ierr; 1113b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1114b4319ba4SBarry Smith 1115b4319ba4SBarry Smith PetscFunctionBegin; 1116b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1117b4319ba4SBarry Smith PetscFunctionReturn(0); 1118b4319ba4SBarry Smith } 1119b4319ba4SBarry Smith 1120b4319ba4SBarry Smith #undef __FUNCT__ 1121f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 1122a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1123f0006bf2SLisandro Dalcin { 1124f0006bf2SLisandro Dalcin PetscErrorCode ierr; 1125f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 1126f0006bf2SLisandro Dalcin 1127f0006bf2SLisandro Dalcin PetscFunctionBegin; 1128f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1129f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 1130f0006bf2SLisandro Dalcin } 1131f0006bf2SLisandro Dalcin 1132f0006bf2SLisandro Dalcin #undef __FUNCT__ 11332e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 1134a8116848SStefano Zampini static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 11352e74eeadSLisandro Dalcin { 11366e520ac8SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 11376e520ac8SStefano Zampini PetscInt nr,nl,len,i; 11386e520ac8SStefano Zampini PetscInt *lrows; 11392e74eeadSLisandro Dalcin PetscErrorCode ierr; 11402e74eeadSLisandro Dalcin 11412e74eeadSLisandro Dalcin PetscFunctionBegin; 11426e520ac8SStefano Zampini /* get locally owned rows */ 11436e520ac8SStefano Zampini ierr = MatZeroRowsMapLocal_Private(A,n,rows,&len,&lrows);CHKERRQ(ierr); 11446e520ac8SStefano Zampini /* fix right hand side if needed */ 11456e520ac8SStefano Zampini if (x && b) { 11466e520ac8SStefano Zampini const PetscScalar *xx; 11476e520ac8SStefano Zampini PetscScalar *bb; 11486e520ac8SStefano Zampini 11496e520ac8SStefano Zampini ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr); 11506e520ac8SStefano Zampini ierr = VecGetArray(b, &bb);CHKERRQ(ierr); 11516e520ac8SStefano Zampini for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]]; 11526e520ac8SStefano Zampini ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr); 11536e520ac8SStefano Zampini ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr); 11542e74eeadSLisandro Dalcin } 11556e520ac8SStefano Zampini /* get rows associated to the local matrices */ 11566e520ac8SStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 11576e520ac8SStefano Zampini ierr = MatISComputeSF_Private(A);CHKERRQ(ierr); 11586e520ac8SStefano Zampini } 11596e520ac8SStefano Zampini ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr); 11606e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr); 11616e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 11626e520ac8SStefano Zampini for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1; 11636e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 11646e520ac8SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 11656e520ac8SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 11666e520ac8SStefano Zampini ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr); 11676e520ac8SStefano Zampini for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i; 11687230de76SStefano Zampini ierr = MatISZeroRowsLocal_Private(A,nr,lrows,diag);CHKERRQ(ierr); 11696e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 11702e74eeadSLisandro Dalcin PetscFunctionReturn(0); 11712e74eeadSLisandro Dalcin } 11722e74eeadSLisandro Dalcin 11732e74eeadSLisandro Dalcin #undef __FUNCT__ 11747230de76SStefano Zampini #define __FUNCT__ "MatISZeroRowsLocal_Private" 11757230de76SStefano Zampini static PetscErrorCode MatISZeroRowsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag) 1176b4319ba4SBarry Smith { 1177b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1178dfbe8321SBarry Smith PetscErrorCode ierr; 1179b4319ba4SBarry Smith 1180b4319ba4SBarry Smith PetscFunctionBegin; 11816e520ac8SStefano Zampini if (diag != 0.) { 1182b4319ba4SBarry Smith /* 1183b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 1184b4319ba4SBarry Smith work properly in the interface nodes. 1185b4319ba4SBarry Smith */ 1186b4319ba4SBarry Smith Vec counter; 1187e176bc59SStefano Zampini ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr); 1188e176bc59SStefano Zampini ierr = VecSet(counter,0.);CHKERRQ(ierr); 1189e176bc59SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 1190e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1191e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1192e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1193e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 11946bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 1195b4319ba4SBarry Smith } 11962b404112SStefano Zampini 1197958c9bccSBarry Smith if (!n) { 1198b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 1199b4319ba4SBarry Smith } else { 12007230de76SStefano Zampini PetscInt i; 1201b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 12022205254eSKarl Rupp 12032b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 12046e520ac8SStefano Zampini if (diag != 0.) { 12056e520ac8SStefano Zampini const PetscScalar *array; 12066e520ac8SStefano Zampini ierr = VecGetArrayRead(is->y,&array);CHKERRQ(ierr); 1207b4319ba4SBarry Smith for (i=0; i<n; i++) { 1208f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 1209b4319ba4SBarry Smith } 12106e520ac8SStefano Zampini ierr = VecRestoreArrayRead(is->y,&array);CHKERRQ(ierr); 12116e520ac8SStefano Zampini } 1212b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1213b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1214b4319ba4SBarry Smith } 1215b4319ba4SBarry Smith PetscFunctionReturn(0); 1216b4319ba4SBarry Smith } 1217b4319ba4SBarry Smith 1218b4319ba4SBarry Smith #undef __FUNCT__ 1219b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 1220a8116848SStefano Zampini static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 1221b4319ba4SBarry Smith { 1222b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1223dfbe8321SBarry Smith PetscErrorCode ierr; 1224b4319ba4SBarry Smith 1225b4319ba4SBarry Smith PetscFunctionBegin; 1226b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 1227b4319ba4SBarry Smith PetscFunctionReturn(0); 1228b4319ba4SBarry Smith } 1229b4319ba4SBarry Smith 1230b4319ba4SBarry Smith #undef __FUNCT__ 1231b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 1232a8116848SStefano Zampini static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 1233b4319ba4SBarry Smith { 1234b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1235dfbe8321SBarry Smith PetscErrorCode ierr; 1236b4319ba4SBarry Smith 1237b4319ba4SBarry Smith PetscFunctionBegin; 1238b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 1239b4319ba4SBarry Smith PetscFunctionReturn(0); 1240b4319ba4SBarry Smith } 1241b4319ba4SBarry Smith 1242b4319ba4SBarry Smith #undef __FUNCT__ 1243b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 1244a8116848SStefano Zampini static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 1245b4319ba4SBarry Smith { 1246b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 1247b4319ba4SBarry Smith 1248b4319ba4SBarry Smith PetscFunctionBegin; 1249b4319ba4SBarry Smith *local = is->A; 1250b4319ba4SBarry Smith PetscFunctionReturn(0); 1251b4319ba4SBarry Smith } 1252b4319ba4SBarry Smith 1253b4319ba4SBarry Smith #undef __FUNCT__ 1254b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 1255b4319ba4SBarry Smith /*@ 1256b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 1257b4319ba4SBarry Smith 1258b4319ba4SBarry Smith Input Parameter: 1259b4319ba4SBarry Smith . mat - the matrix 1260b4319ba4SBarry Smith 1261b4319ba4SBarry Smith Output Parameter: 1262eb82efa4SStefano Zampini . local - the local matrix 1263b4319ba4SBarry Smith 1264b4319ba4SBarry Smith Level: advanced 1265b4319ba4SBarry Smith 1266b4319ba4SBarry Smith Notes: 1267b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 1268b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 1269b4319ba4SBarry Smith of the MatSetValues() operation. 1270b4319ba4SBarry Smith 1271b4319ba4SBarry Smith .seealso: MATIS 1272b4319ba4SBarry Smith @*/ 12737087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 1274b4319ba4SBarry Smith { 12754ac538c5SBarry Smith PetscErrorCode ierr; 1276b4319ba4SBarry Smith 1277b4319ba4SBarry Smith PetscFunctionBegin; 12780700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1279b4319ba4SBarry Smith PetscValidPointer(local,2); 12804ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 1281b4319ba4SBarry Smith PetscFunctionReturn(0); 1282b4319ba4SBarry Smith } 1283b4319ba4SBarry Smith 12843b03a366Sstefano_zampini #undef __FUNCT__ 12853b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 1286a8116848SStefano Zampini static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 12873b03a366Sstefano_zampini { 12883b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 12893b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 12903b03a366Sstefano_zampini PetscErrorCode ierr; 12913b03a366Sstefano_zampini 12923b03a366Sstefano_zampini PetscFunctionBegin; 12934e4c7dbeSStefano Zampini if (is->A) { 12943b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 12953b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 1296f23aa3ddSBarry 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); 12974e4c7dbeSStefano Zampini } 12983b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 12993b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 13003b03a366Sstefano_zampini is->A = local; 13013b03a366Sstefano_zampini PetscFunctionReturn(0); 13023b03a366Sstefano_zampini } 13033b03a366Sstefano_zampini 13043b03a366Sstefano_zampini #undef __FUNCT__ 13053b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 13063b03a366Sstefano_zampini /*@ 1307eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 13083b03a366Sstefano_zampini 13093b03a366Sstefano_zampini Input Parameter: 13103b03a366Sstefano_zampini . mat - the matrix 1311eb82efa4SStefano Zampini . local - the local matrix 13123b03a366Sstefano_zampini 13133b03a366Sstefano_zampini Output Parameter: 13143b03a366Sstefano_zampini 13153b03a366Sstefano_zampini Level: advanced 13163b03a366Sstefano_zampini 13173b03a366Sstefano_zampini Notes: 13183b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 13193b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 13203b03a366Sstefano_zampini 13213b03a366Sstefano_zampini .seealso: MATIS 13223b03a366Sstefano_zampini @*/ 13233b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 13243b03a366Sstefano_zampini { 13253b03a366Sstefano_zampini PetscErrorCode ierr; 13263b03a366Sstefano_zampini 13273b03a366Sstefano_zampini PetscFunctionBegin; 13283b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1329b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 13303b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 13313b03a366Sstefano_zampini PetscFunctionReturn(0); 13323b03a366Sstefano_zampini } 13333b03a366Sstefano_zampini 13346726f965SBarry Smith #undef __FUNCT__ 13356726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 1336a8116848SStefano Zampini static PetscErrorCode MatZeroEntries_IS(Mat A) 13376726f965SBarry Smith { 13386726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 13396726f965SBarry Smith PetscErrorCode ierr; 13406726f965SBarry Smith 13416726f965SBarry Smith PetscFunctionBegin; 13426726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 13436726f965SBarry Smith PetscFunctionReturn(0); 13446726f965SBarry Smith } 13456726f965SBarry Smith 13466726f965SBarry Smith #undef __FUNCT__ 13472e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 1348a8116848SStefano Zampini static PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 13492e74eeadSLisandro Dalcin { 13502e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 13512e74eeadSLisandro Dalcin PetscErrorCode ierr; 13522e74eeadSLisandro Dalcin 13532e74eeadSLisandro Dalcin PetscFunctionBegin; 13542e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 13552e74eeadSLisandro Dalcin PetscFunctionReturn(0); 13562e74eeadSLisandro Dalcin } 13572e74eeadSLisandro Dalcin 13582e74eeadSLisandro Dalcin #undef __FUNCT__ 13592e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 1360a8116848SStefano Zampini static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 13612e74eeadSLisandro Dalcin { 13622e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 13632e74eeadSLisandro Dalcin PetscErrorCode ierr; 13642e74eeadSLisandro Dalcin 13652e74eeadSLisandro Dalcin PetscFunctionBegin; 13662e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 1367e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 13682e74eeadSLisandro Dalcin 13692e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 13702e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 1371e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1372e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 13732e74eeadSLisandro Dalcin PetscFunctionReturn(0); 13742e74eeadSLisandro Dalcin } 13752e74eeadSLisandro Dalcin 13762e74eeadSLisandro Dalcin #undef __FUNCT__ 13776726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 1378a8116848SStefano Zampini static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 13796726f965SBarry Smith { 13806726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 13816726f965SBarry Smith PetscErrorCode ierr; 13826726f965SBarry Smith 13836726f965SBarry Smith PetscFunctionBegin; 13844e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 13856726f965SBarry Smith PetscFunctionReturn(0); 13866726f965SBarry Smith } 13876726f965SBarry Smith 1388284134d9SBarry Smith #undef __FUNCT__ 1389284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 1390284134d9SBarry Smith /*@ 13913c212e90SHong Zhang MatCreateIS - Creates a "process" unassembled matrix, assembled on each 1392284134d9SBarry Smith process but not across processes. 1393284134d9SBarry Smith 1394284134d9SBarry Smith Input Parameters: 1395284134d9SBarry Smith + comm - MPI communicator that will share the matrix 1396e176bc59SStefano Zampini . bs - block size of the matrix 1397df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 1398e176bc59SStefano Zampini . rmap - local to global map for rows 1399e176bc59SStefano Zampini - cmap - local to global map for cols 1400284134d9SBarry Smith 1401284134d9SBarry Smith Output Parameter: 1402284134d9SBarry Smith . A - the resulting matrix 1403284134d9SBarry Smith 14048e6c10adSSatish Balay Level: advanced 14058e6c10adSSatish Balay 14063c212e90SHong Zhang Notes: See MATIS for more details. 14073c212e90SHong Zhang m and n are NOT related to the size of the map; they are the size of the part of the vector owned 1408e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 14093c212e90SHong Zhang If either rmap or cmap are NULL, then the matrix is assumed to be square. 1410284134d9SBarry Smith 1411284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 1412284134d9SBarry Smith @*/ 1413e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 1414284134d9SBarry Smith { 1415284134d9SBarry Smith PetscErrorCode ierr; 1416284134d9SBarry Smith 1417284134d9SBarry Smith PetscFunctionBegin; 1418e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 1419284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1420d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 1421284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1422284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 1423e176bc59SStefano Zampini if (rmap && cmap) { 1424e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 1425e176bc59SStefano Zampini } else if (!rmap) { 1426e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 1427e176bc59SStefano Zampini } else { 1428e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 1429e176bc59SStefano Zampini } 1430284134d9SBarry Smith PetscFunctionReturn(0); 1431284134d9SBarry Smith } 1432284134d9SBarry Smith 1433b4319ba4SBarry Smith /*MC 1434eb82efa4SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC). 1435b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 1436b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 1437b4319ba4SBarry Smith product is handled "implicitly". 1438b4319ba4SBarry Smith 1439b4319ba4SBarry Smith Operations Provided: 14406726f965SBarry Smith + MatMult() 14412e74eeadSLisandro Dalcin . MatMultAdd() 14422e74eeadSLisandro Dalcin . MatMultTranspose() 14432e74eeadSLisandro Dalcin . MatMultTransposeAdd() 14446726f965SBarry Smith . MatZeroEntries() 14456726f965SBarry Smith . MatSetOption() 14462e74eeadSLisandro Dalcin . MatZeroRows() 14472e74eeadSLisandro Dalcin . MatSetValues() 144897563a80SStefano Zampini . MatSetValuesBlocked() 14496726f965SBarry Smith . MatSetValuesLocal() 145097563a80SStefano Zampini . MatSetValuesBlockedLocal() 14512e74eeadSLisandro Dalcin . MatScale() 14522e74eeadSLisandro Dalcin . MatGetDiagonal() 14532b404112SStefano Zampini . MatMissingDiagonal() 14542b404112SStefano Zampini . MatDuplicate() 14552b404112SStefano Zampini . MatCopy() 14566726f965SBarry Smith - MatSetLocalToGlobalMapping() 1457b4319ba4SBarry Smith 1458b4319ba4SBarry Smith Options Database Keys: 1459b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1460b4319ba4SBarry Smith 1461b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1462b4319ba4SBarry Smith 1463b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1464b4319ba4SBarry Smith 1465b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1466eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1467b4319ba4SBarry Smith 1468b4319ba4SBarry Smith Level: advanced 1469b4319ba4SBarry Smith 14703c212e90SHong Zhang .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC, MatCreateIS() 1471b4319ba4SBarry Smith 1472b4319ba4SBarry Smith M*/ 1473b4319ba4SBarry Smith 1474b4319ba4SBarry Smith #undef __FUNCT__ 1475b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 14768cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1477b4319ba4SBarry Smith { 1478dfbe8321SBarry Smith PetscErrorCode ierr; 1479b4319ba4SBarry Smith Mat_IS *b; 1480b4319ba4SBarry Smith 1481b4319ba4SBarry Smith PetscFunctionBegin; 1482b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1483b4319ba4SBarry Smith A->data = (void*)b; 1484b4319ba4SBarry Smith 1485e176bc59SStefano Zampini /* matrix ops */ 1486e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1487b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 14882e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 14892e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 14902e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1491b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1492b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 14932e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 149498921651SStefano Zampini A->ops->setvaluesblocked = MatSetValuesBlocked_IS; 1495b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1496f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 14972e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1498b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1499b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1500b4319ba4SBarry Smith A->ops->view = MatView_IS; 15016726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 15022e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 15032e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 15046726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 150569796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 150669796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1507ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 15086bd84002SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_IS; 15092b404112SStefano Zampini A->ops->copy = MatCopy_IS; 1510659959c5SStefano Zampini A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_IS; 1511a8116848SStefano Zampini A->ops->getsubmatrix = MatGetSubMatrix_IS; 1512b4319ba4SBarry Smith 1513b7ce53b6SStefano Zampini /* special MATIS functions */ 1514bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1515bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1516b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 15172e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 151817667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1519b4319ba4SBarry Smith PetscFunctionReturn(0); 1520b4319ba4SBarry Smith } 1521