1b4319ba4SBarry Smith /* 2b4319ba4SBarry Smith Creates a matrix class for using the Neumann-Neumann type preconditioners. 3b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 4b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 5b4319ba4SBarry Smith product is handled "implicitly". 6b4319ba4SBarry Smith 7b4319ba4SBarry Smith Currently this allows for only one subdomain per processor. 8b4319ba4SBarry Smith */ 9b4319ba4SBarry Smith 10c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/ 11*6989cf23SStefano Zampini #include <../src/mat/impls/aij/mpi/mpiaij.h> 124f2d7cafSStefano Zampini #include <petsc/private/sfimpl.h> 1328f4e0baSStefano Zampini 14f26d0771SStefano Zampini #define MATIS_MAX_ENTRIES_INSERTION 2048 15f26d0771SStefano Zampini 16cf0a3239SStefano Zampini #undef __FUNCT__ 17*6989cf23SStefano Zampini #define __FUNCT__ "MatISContainerDestroy_Private" 18*6989cf23SStefano Zampini static PetscErrorCode MatISContainerDestroy_Private(void *ptr) 19*6989cf23SStefano Zampini { 20*6989cf23SStefano Zampini PetscErrorCode ierr; 21*6989cf23SStefano Zampini 22*6989cf23SStefano Zampini PetscFunctionBeginUser; 23*6989cf23SStefano Zampini ierr = PetscFree(ptr); CHKERRQ(ierr); 24*6989cf23SStefano Zampini PetscFunctionReturn(0); 25*6989cf23SStefano Zampini } 26*6989cf23SStefano Zampini 27*6989cf23SStefano Zampini #undef __FUNCT__ 28*6989cf23SStefano Zampini #define __FUNCT__ "MatConvert_MPIAIJ_IS" 29*6989cf23SStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_IS(Mat A,MatType type,MatReuse reuse,Mat *newmat) 30*6989cf23SStefano Zampini { 31*6989cf23SStefano Zampini Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 32*6989cf23SStefano Zampini Mat_SeqAIJ *diag = (Mat_SeqAIJ*)(aij->A->data); 33*6989cf23SStefano Zampini Mat_SeqAIJ *offd = (Mat_SeqAIJ*)(aij->B->data); 34*6989cf23SStefano Zampini Mat lA; 35*6989cf23SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 36*6989cf23SStefano Zampini IS is; 37*6989cf23SStefano Zampini MPI_Comm comm; 38*6989cf23SStefano Zampini void *ptrs[2]; 39*6989cf23SStefano Zampini const char *names[2] = {"_convert_csr_aux","_convert_csr_data"}; 40*6989cf23SStefano Zampini PetscScalar *dd,*od,*aa,*data; 41*6989cf23SStefano Zampini PetscInt *di,*dj,*oi,*oj; 42*6989cf23SStefano Zampini PetscInt *aux,*ii,*jj; 43*6989cf23SStefano Zampini PetscInt dr,dc,oc,str,stc,nnz,i,jd,jo,cum; 44*6989cf23SStefano Zampini PetscErrorCode ierr; 45*6989cf23SStefano Zampini 46*6989cf23SStefano Zampini PetscFunctionBeginUser; 47*6989cf23SStefano Zampini ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 48*6989cf23SStefano Zampini if (!aij->garray) SETERRQ(comm,PETSC_ERR_SUP,"garray not present"); 49*6989cf23SStefano Zampini 50*6989cf23SStefano Zampini /* access relevant information from MPIAIJ */ 51*6989cf23SStefano Zampini ierr = MatGetOwnershipRange(A,&str,NULL);CHKERRQ(ierr); 52*6989cf23SStefano Zampini ierr = MatGetOwnershipRangeColumn(A,&stc,NULL);CHKERRQ(ierr); 53*6989cf23SStefano Zampini ierr = MatGetLocalSize(A,&dr,&dc);CHKERRQ(ierr); 54*6989cf23SStefano Zampini di = diag->i; 55*6989cf23SStefano Zampini dj = diag->j; 56*6989cf23SStefano Zampini dd = diag->a; 57*6989cf23SStefano Zampini oc = aij->B->cmap->n; 58*6989cf23SStefano Zampini oi = offd->i; 59*6989cf23SStefano Zampini oj = offd->j; 60*6989cf23SStefano Zampini od = offd->a; 61*6989cf23SStefano Zampini nnz = diag->i[dr] + offd->i[dr]; 62*6989cf23SStefano Zampini 63*6989cf23SStefano Zampini /* generate l2g maps for rows and cols */ 64*6989cf23SStefano Zampini ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr); 65*6989cf23SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 66*6989cf23SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 67*6989cf23SStefano Zampini ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr); 68*6989cf23SStefano Zampini for (i=0; i<dc; i++) aux[i] = i+stc; 69*6989cf23SStefano Zampini for (i=0; i<oc; i++) aux[i+dc] = aij->garray[i]; 70*6989cf23SStefano Zampini ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 71*6989cf23SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 72*6989cf23SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 73*6989cf23SStefano Zampini 74*6989cf23SStefano Zampini /* create MATIS object */ 75*6989cf23SStefano Zampini ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 76*6989cf23SStefano Zampini ierr = MatSetSizes(*newmat,dr,dc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 77*6989cf23SStefano Zampini ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr); 78*6989cf23SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr); 79*6989cf23SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 80*6989cf23SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 81*6989cf23SStefano Zampini 82*6989cf23SStefano Zampini /* merge local matrices */ 83*6989cf23SStefano Zampini ierr = PetscMalloc1(nnz+dr+1,&aux);CHKERRQ(ierr); 84*6989cf23SStefano Zampini ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr); 85*6989cf23SStefano Zampini ii = aux; 86*6989cf23SStefano Zampini jj = aux+dr+1; 87*6989cf23SStefano Zampini aa = data; 88*6989cf23SStefano Zampini *ii = *(di++) + *(oi++); 89*6989cf23SStefano Zampini for (jd=0,jo=0,cum=0;*ii<nnz;cum++) 90*6989cf23SStefano Zampini { 91*6989cf23SStefano Zampini for (;jd<*di;jd++) { *jj++ = *dj++; *aa++ = *dd++; } 92*6989cf23SStefano Zampini for (;jo<*oi;jo++) { *jj++ = *oj++ + dc; *aa++ = *od++; } 93*6989cf23SStefano Zampini *(++ii) = *(di++) + *(oi++); 94*6989cf23SStefano Zampini } 95*6989cf23SStefano Zampini for (;cum<dr;cum++) *(++ii) = nnz; 96*6989cf23SStefano Zampini ii = aux; 97*6989cf23SStefano Zampini jj = aux+dr+1; 98*6989cf23SStefano Zampini aa = data; 99*6989cf23SStefano Zampini ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,ii,jj,aa,&lA);CHKERRQ(ierr); 100*6989cf23SStefano Zampini 101*6989cf23SStefano Zampini /* create containers to destroy the data */ 102*6989cf23SStefano Zampini ptrs[0] = aux; 103*6989cf23SStefano Zampini ptrs[1] = data; 104*6989cf23SStefano Zampini for (i=0; i<2; i++) { 105*6989cf23SStefano Zampini PetscContainer c; 106*6989cf23SStefano Zampini 107*6989cf23SStefano Zampini ierr = PetscContainerCreate(PETSC_COMM_SELF,&c);CHKERRQ(ierr); 108*6989cf23SStefano Zampini ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr); 109*6989cf23SStefano Zampini ierr = PetscContainerSetUserDestroy(c,MatISContainerDestroy_Private);CHKERRQ(ierr); 110*6989cf23SStefano Zampini ierr = PetscObjectCompose((PetscObject)lA,names[i],(PetscObject)c);CHKERRQ(ierr); 111*6989cf23SStefano Zampini ierr = PetscContainerDestroy(&c);CHKERRQ(ierr); 112*6989cf23SStefano Zampini } 113*6989cf23SStefano Zampini 114*6989cf23SStefano Zampini /* finalize matrix */ 115*6989cf23SStefano Zampini ierr = MatISSetLocalMat(*newmat,lA);CHKERRQ(ierr); 116*6989cf23SStefano Zampini ierr = MatDestroy(&lA);CHKERRQ(ierr); 117*6989cf23SStefano Zampini ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 118*6989cf23SStefano Zampini ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 119*6989cf23SStefano Zampini PetscFunctionReturn(0); 120*6989cf23SStefano Zampini } 121*6989cf23SStefano Zampini 122*6989cf23SStefano Zampini #undef __FUNCT__ 123cf0a3239SStefano Zampini #define __FUNCT__ "MatISSetUpSF" 124cf0a3239SStefano Zampini /*@ 1253d996552SStefano Zampini MatISSetUpSF - Setup star forest objects used by MatIS. 126cf0a3239SStefano Zampini 127cf0a3239SStefano Zampini Collective on MPI_Comm 128cf0a3239SStefano Zampini 129cf0a3239SStefano Zampini Input Parameters: 130cf0a3239SStefano Zampini + A - the matrix 131cf0a3239SStefano Zampini 132cf0a3239SStefano Zampini Level: advanced 133cf0a3239SStefano Zampini 1343d996552SStefano Zampini Notes: This function does not need to be called by the user. 135cf0a3239SStefano Zampini 136cf0a3239SStefano Zampini .keywords: matrix 137cf0a3239SStefano Zampini 138cf0a3239SStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatISSetPreallocation(), MatISGetLocalMat() 139cf0a3239SStefano Zampini @*/ 140cf0a3239SStefano Zampini PetscErrorCode MatISSetUpSF(Mat A) 141cf0a3239SStefano Zampini { 142cf0a3239SStefano Zampini PetscErrorCode ierr; 143cf0a3239SStefano Zampini 144cf0a3239SStefano Zampini PetscFunctionBegin; 145cf0a3239SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 146cf0a3239SStefano Zampini PetscValidType(A,1); 147cf0a3239SStefano Zampini ierr = PetscTryMethod(A,"MatISSetUpSF_C",(Mat),(A));CHKERRQ(ierr); 148cf0a3239SStefano Zampini PetscFunctionReturn(0); 149cf0a3239SStefano Zampini } 150a72627d2SStefano Zampini 15128f4e0baSStefano Zampini #undef __FUNCT__ 1523fd1c9e7SStefano Zampini #define __FUNCT__ "MatDiagonalSet_IS" 1533fd1c9e7SStefano Zampini PetscErrorCode MatDiagonalSet_IS(Mat A,Vec D,InsertMode insmode) 1543fd1c9e7SStefano Zampini { 1553fd1c9e7SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 1563fd1c9e7SStefano Zampini PetscErrorCode ierr; 1573fd1c9e7SStefano Zampini 1583fd1c9e7SStefano Zampini PetscFunctionBegin; 1593fd1c9e7SStefano Zampini if (!D) { /* this code branch is used by MatShift_IS */ 1603fd1c9e7SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 1613fd1c9e7SStefano Zampini } else { 1623fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1633fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1643fd1c9e7SStefano Zampini } 1653fd1c9e7SStefano Zampini ierr = VecPointwiseDivide(is->y,is->y,is->counter);CHKERRQ(ierr); 1663fd1c9e7SStefano Zampini ierr = MatDiagonalSet(is->A,is->y,insmode);CHKERRQ(ierr); 1673fd1c9e7SStefano Zampini PetscFunctionReturn(0); 1683fd1c9e7SStefano Zampini } 1693fd1c9e7SStefano Zampini 1703fd1c9e7SStefano Zampini #undef __FUNCT__ 1713fd1c9e7SStefano Zampini #define __FUNCT__ "MatShift_IS" 1723fd1c9e7SStefano Zampini PetscErrorCode MatShift_IS(Mat A,PetscScalar a) 1733fd1c9e7SStefano Zampini { 1743fd1c9e7SStefano Zampini PetscErrorCode ierr; 1753fd1c9e7SStefano Zampini 1763fd1c9e7SStefano Zampini PetscFunctionBegin; 1773fd1c9e7SStefano Zampini ierr = MatDiagonalSet_IS(A,NULL,ADD_VALUES);CHKERRQ(ierr); 1783fd1c9e7SStefano Zampini PetscFunctionReturn(0); 1793fd1c9e7SStefano Zampini } 1803fd1c9e7SStefano Zampini 1813fd1c9e7SStefano Zampini #undef __FUNCT__ 182f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesLocal_SubMat_IS" 183f26d0771SStefano Zampini static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 184f26d0771SStefano Zampini { 185f26d0771SStefano Zampini PetscErrorCode ierr; 186f26d0771SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 187f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 188f26d0771SStefano Zampini 189f26d0771SStefano Zampini PetscFunctionBegin; 190f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 191f26d0771SStefano Zampini if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n); 192f26d0771SStefano Zampini #endif 193f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 194f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 195f26d0771SStefano Zampini ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 196f26d0771SStefano Zampini PetscFunctionReturn(0); 197f26d0771SStefano Zampini } 198f26d0771SStefano Zampini 199f26d0771SStefano Zampini #undef __FUNCT__ 200f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesBlockedLocal_SubMat_IS" 201f26d0771SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 202f26d0771SStefano Zampini { 203f26d0771SStefano Zampini PetscErrorCode ierr; 204f26d0771SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 205f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 206f26d0771SStefano Zampini 207f26d0771SStefano Zampini PetscFunctionBegin; 208f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 209f26d0771SStefano Zampini if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column block indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n); 210f26d0771SStefano Zampini #endif 211f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 212f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 213f26d0771SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 214f26d0771SStefano Zampini PetscFunctionReturn(0); 215f26d0771SStefano Zampini } 216f26d0771SStefano Zampini 217f26d0771SStefano Zampini #undef __FUNCT__ 218a8116848SStefano Zampini #define __FUNCT__ "PetscLayoutMapLocal_Private" 219f0ae7da4SStefano Zampini static PetscErrorCode PetscLayoutMapLocal_Private(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs) 220a8116848SStefano Zampini { 221a8116848SStefano Zampini PetscInt *owners = map->range; 222a8116848SStefano Zampini PetscInt n = map->n; 223a8116848SStefano Zampini PetscSF sf; 224a8116848SStefano Zampini PetscInt *lidxs,*work = NULL; 225a8116848SStefano Zampini PetscSFNode *ridxs; 226a8116848SStefano Zampini PetscMPIInt rank; 227a8116848SStefano Zampini PetscInt r, p = 0, len = 0; 228a8116848SStefano Zampini PetscErrorCode ierr; 229a8116848SStefano Zampini 230a8116848SStefano Zampini PetscFunctionBegin; 231a8116848SStefano Zampini /* Create SF where leaves are input idxs and roots are owned idxs (code adapted from MatZeroRowsMapLocal_Private) */ 232a8116848SStefano Zampini ierr = MPI_Comm_rank(map->comm,&rank);CHKERRQ(ierr); 233a8116848SStefano Zampini ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr); 234a8116848SStefano Zampini for (r = 0; r < n; ++r) lidxs[r] = -1; 235a8116848SStefano Zampini ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr); 236a8116848SStefano Zampini for (r = 0; r < N; ++r) { 237a8116848SStefano Zampini const PetscInt idx = idxs[r]; 238a8116848SStefano 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); 239a8116848SStefano Zampini if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 240a8116848SStefano Zampini ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr); 241a8116848SStefano Zampini } 242a8116848SStefano Zampini ridxs[r].rank = p; 243a8116848SStefano Zampini ridxs[r].index = idxs[r] - owners[p]; 244a8116848SStefano Zampini } 245a8116848SStefano Zampini ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr); 246a8116848SStefano Zampini ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr); 247a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 248a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 249f0ae7da4SStefano Zampini if (ogidxs) { /* communicate global idxs */ 250a8116848SStefano Zampini PetscInt cum = 0,start,*work2; 251f0ae7da4SStefano Zampini 252f0ae7da4SStefano Zampini ierr = PetscMalloc1(n,&work);CHKERRQ(ierr); 253a8116848SStefano Zampini ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr); 254a8116848SStefano Zampini for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++; 255a8116848SStefano Zampini ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRQ(ierr); 256a8116848SStefano Zampini start -= cum; 257a8116848SStefano Zampini cum = 0; 258a8116848SStefano Zampini for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++; 259a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 260a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 261a8116848SStefano Zampini ierr = PetscFree(work2);CHKERRQ(ierr); 262a8116848SStefano Zampini } 263a8116848SStefano Zampini ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 264a8116848SStefano Zampini /* Compress and put in indices */ 265a8116848SStefano Zampini for (r = 0; r < n; ++r) 266a8116848SStefano Zampini if (lidxs[r] >= 0) { 267a8116848SStefano Zampini if (work) work[len] = work[r]; 268a8116848SStefano Zampini lidxs[len++] = r; 269a8116848SStefano Zampini } 270a8116848SStefano Zampini if (on) *on = len; 271a8116848SStefano Zampini if (oidxs) *oidxs = lidxs; 272a8116848SStefano Zampini if (ogidxs) *ogidxs = work; 273a8116848SStefano Zampini PetscFunctionReturn(0); 274a8116848SStefano Zampini } 275a8116848SStefano Zampini 276a8116848SStefano Zampini #undef __FUNCT__ 277a8116848SStefano Zampini #define __FUNCT__ "MatGetSubMatrix_IS" 278a8116848SStefano Zampini static PetscErrorCode MatGetSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat) 279a8116848SStefano Zampini { 280a8116848SStefano Zampini Mat locmat,newlocmat; 281a8116848SStefano Zampini Mat_IS *newmatis; 282a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG) 283a8116848SStefano Zampini Vec rtest,ltest; 284a8116848SStefano Zampini const PetscScalar *array; 285a8116848SStefano Zampini #endif 286a8116848SStefano Zampini const PetscInt *idxs; 287a8116848SStefano Zampini PetscInt i,m,n; 288a8116848SStefano Zampini PetscErrorCode ierr; 289a8116848SStefano Zampini 290a8116848SStefano Zampini PetscFunctionBegin; 291a8116848SStefano Zampini if (scall == MAT_REUSE_MATRIX) { 292a8116848SStefano Zampini PetscBool ismatis; 293a8116848SStefano Zampini 294a8116848SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr); 295a8116848SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type"); 296a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 297a8116848SStefano Zampini if (!newmatis->getsub_ris) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS"); 298a8116848SStefano Zampini if (!newmatis->getsub_cis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS"); 299a8116848SStefano Zampini } 300a8116848SStefano Zampini /* irow and icol may not have duplicate entries */ 301a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG) 302a8116848SStefano Zampini ierr = MatCreateVecs(mat,<est,&rtest);CHKERRQ(ierr); 303a8116848SStefano Zampini ierr = ISGetLocalSize(irow,&n);CHKERRQ(ierr); 304a8116848SStefano Zampini ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 305a8116848SStefano Zampini for (i=0;i<n;i++) { 306a8116848SStefano Zampini ierr = VecSetValue(rtest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 307a8116848SStefano Zampini } 308a8116848SStefano Zampini ierr = VecAssemblyBegin(rtest);CHKERRQ(ierr); 309a8116848SStefano Zampini ierr = VecAssemblyEnd(rtest);CHKERRQ(ierr); 310a8116848SStefano Zampini ierr = VecGetLocalSize(rtest,&n);CHKERRQ(ierr); 311a8116848SStefano Zampini ierr = VecGetOwnershipRange(rtest,&m,NULL);CHKERRQ(ierr); 312a8116848SStefano Zampini ierr = VecGetArrayRead(rtest,&array);CHKERRQ(ierr); 313fd479f66SMatthew 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])); 314a8116848SStefano Zampini ierr = VecRestoreArrayRead(rtest,&array);CHKERRQ(ierr); 315a8116848SStefano Zampini ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 316a8116848SStefano Zampini ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 317a8116848SStefano Zampini ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 318a8116848SStefano Zampini for (i=0;i<n;i++) { 319a8116848SStefano Zampini ierr = VecSetValue(ltest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 320a8116848SStefano Zampini } 321a8116848SStefano Zampini ierr = VecAssemblyBegin(ltest);CHKERRQ(ierr); 322a8116848SStefano Zampini ierr = VecAssemblyEnd(ltest);CHKERRQ(ierr); 323a8116848SStefano Zampini ierr = VecGetLocalSize(ltest,&n);CHKERRQ(ierr); 324a8116848SStefano Zampini ierr = VecGetOwnershipRange(ltest,&m,NULL);CHKERRQ(ierr); 325a8116848SStefano Zampini ierr = VecGetArrayRead(ltest,&array);CHKERRQ(ierr); 326fd479f66SMatthew 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])); 327a8116848SStefano Zampini ierr = VecRestoreArrayRead(ltest,&array);CHKERRQ(ierr); 328a8116848SStefano Zampini ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 329a8116848SStefano Zampini ierr = VecDestroy(&rtest);CHKERRQ(ierr); 330a8116848SStefano Zampini ierr = VecDestroy(<est);CHKERRQ(ierr); 331a8116848SStefano Zampini #endif 332a8116848SStefano Zampini if (scall == MAT_INITIAL_MATRIX) { 333a8116848SStefano Zampini Mat_IS *matis = (Mat_IS*)mat->data; 334a8116848SStefano Zampini ISLocalToGlobalMapping rl2g; 335a8116848SStefano Zampini IS is; 336a8116848SStefano Zampini PetscInt *lidxs,*lgidxs,*newgidxs; 3376dd40735SStefano Zampini PetscInt ll,newloc; 338a8116848SStefano Zampini MPI_Comm comm; 339a8116848SStefano Zampini 340a8116848SStefano Zampini ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 341a8116848SStefano Zampini ierr = ISGetLocalSize(irow,&m);CHKERRQ(ierr); 342a8116848SStefano Zampini ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 343a8116848SStefano Zampini ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 344a8116848SStefano Zampini ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr); 345a8116848SStefano Zampini ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 346a8116848SStefano Zampini /* communicate irow to their owners in the layout */ 347a8116848SStefano Zampini ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 348f0ae7da4SStefano Zampini ierr = PetscLayoutMapLocal_Private(mat->rmap,m,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 349a8116848SStefano Zampini ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 3503d996552SStefano Zampini ierr = MatISSetUpSF(mat);CHKERRQ(ierr); 3513d996552SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,matis->sf->nroots*sizeof(PetscInt));CHKERRQ(ierr); 352a8116848SStefano Zampini for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1; 353a8116848SStefano Zampini ierr = PetscFree(lidxs);CHKERRQ(ierr); 354a8116848SStefano Zampini ierr = PetscFree(lgidxs);CHKERRQ(ierr); 355a8116848SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 356a8116848SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 3573d996552SStefano Zampini for (i=0,newloc=0;i<matis->sf->nleaves;i++) if (matis->sf_leafdata[i]) newloc++; 358a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 359a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 3603d996552SStefano Zampini for (i=0,newloc=0;i<matis->sf->nleaves;i++) 361a8116848SStefano Zampini if (matis->sf_leafdata[i]) { 362a8116848SStefano Zampini lidxs[newloc] = i; 363a8116848SStefano Zampini newgidxs[newloc++] = matis->sf_leafdata[i]-1; 364a8116848SStefano Zampini } 365a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 366a8116848SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 367a8116848SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 368a8116848SStefano Zampini /* local is to extract local submatrix */ 369a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 370a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris);CHKERRQ(ierr); 371a8116848SStefano Zampini if (mat->congruentlayouts == -1) { /* first time we compare rows and cols layouts */ 372a8116848SStefano Zampini PetscBool cong; 373a8116848SStefano Zampini ierr = PetscLayoutCompare(mat->rmap,mat->cmap,&cong);CHKERRQ(ierr); 374a8116848SStefano Zampini if (cong) mat->congruentlayouts = 1; 375a8116848SStefano Zampini else mat->congruentlayouts = 0; 376a8116848SStefano Zampini } 377a8116848SStefano Zampini if (mat->congruentlayouts && irow == icol) { 378a8116848SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g);CHKERRQ(ierr); 379a8116848SStefano Zampini ierr = PetscObjectReference((PetscObject)newmatis->getsub_ris);CHKERRQ(ierr); 380a8116848SStefano Zampini newmatis->getsub_cis = newmatis->getsub_ris; 381a8116848SStefano Zampini } else { 382a8116848SStefano Zampini ISLocalToGlobalMapping cl2g; 383a8116848SStefano Zampini 384a8116848SStefano Zampini /* communicate icol to their owners in the layout */ 385a8116848SStefano Zampini ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 386f0ae7da4SStefano Zampini ierr = PetscLayoutMapLocal_Private(mat->cmap,n,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 387a8116848SStefano Zampini ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 3883d996552SStefano Zampini ierr = PetscMemzero(matis->csf_rootdata,matis->csf->nroots*sizeof(PetscInt));CHKERRQ(ierr); 389a8116848SStefano Zampini for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1; 390a8116848SStefano Zampini ierr = PetscFree(lidxs);CHKERRQ(ierr); 391a8116848SStefano Zampini ierr = PetscFree(lgidxs);CHKERRQ(ierr); 392a8116848SStefano Zampini ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 393a8116848SStefano Zampini ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 3943d996552SStefano Zampini for (i=0,newloc=0;i<matis->csf->nleaves;i++) if (matis->csf_leafdata[i]) newloc++; 395a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 396a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 3973d996552SStefano Zampini for (i=0,newloc=0;i<matis->csf->nleaves;i++) 398a8116848SStefano Zampini if (matis->csf_leafdata[i]) { 399a8116848SStefano Zampini lidxs[newloc] = i; 400a8116848SStefano Zampini newgidxs[newloc++] = matis->csf_leafdata[i]-1; 401a8116848SStefano Zampini } 402a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 403a8116848SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 404a8116848SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 405a8116848SStefano Zampini /* local is to extract local submatrix */ 406a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis);CHKERRQ(ierr); 407a8116848SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr); 408a8116848SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 409a8116848SStefano Zampini } 410a8116848SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 411a8116848SStefano Zampini } else { 412a8116848SStefano Zampini ierr = MatISGetLocalMat(*newmat,&newlocmat);CHKERRQ(ierr); 413a8116848SStefano Zampini } 414a8116848SStefano Zampini ierr = MatISGetLocalMat(mat,&locmat);CHKERRQ(ierr); 415a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 416a8116848SStefano Zampini ierr = MatGetSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat);CHKERRQ(ierr); 417a8116848SStefano Zampini if (scall == MAT_INITIAL_MATRIX) { 418a8116848SStefano Zampini ierr = MatISSetLocalMat(*newmat,newlocmat);CHKERRQ(ierr); 419a8116848SStefano Zampini ierr = MatDestroy(&newlocmat);CHKERRQ(ierr); 420a8116848SStefano Zampini } 421a8116848SStefano Zampini ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 422a8116848SStefano Zampini ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 423a8116848SStefano Zampini PetscFunctionReturn(0); 424a8116848SStefano Zampini } 425a8116848SStefano Zampini 426a8116848SStefano Zampini #undef __FUNCT__ 4272b404112SStefano Zampini #define __FUNCT__ "MatCopy_IS" 428a8116848SStefano Zampini static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str) 4292b404112SStefano Zampini { 4302b404112SStefano Zampini Mat_IS *a = (Mat_IS*)A->data,*b; 4312b404112SStefano Zampini PetscBool ismatis; 4322b404112SStefano Zampini PetscErrorCode ierr; 4332b404112SStefano Zampini 4342b404112SStefano Zampini PetscFunctionBegin; 4352b404112SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr); 4362b404112SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented"); 4372b404112SStefano Zampini b = (Mat_IS*)B->data; 4382b404112SStefano Zampini ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr); 4392b404112SStefano Zampini PetscFunctionReturn(0); 4402b404112SStefano Zampini } 4412b404112SStefano Zampini 4422b404112SStefano Zampini #undef __FUNCT__ 4436bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS" 444a8116848SStefano Zampini static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool *missing,PetscInt *d) 4456bd84002SStefano Zampini { 446527b2640SStefano Zampini Vec v; 447527b2640SStefano Zampini const PetscScalar *array; 448527b2640SStefano Zampini PetscInt i,n; 4496bd84002SStefano Zampini PetscErrorCode ierr; 4506bd84002SStefano Zampini 4516bd84002SStefano Zampini PetscFunctionBegin; 452527b2640SStefano Zampini *missing = PETSC_FALSE; 453527b2640SStefano Zampini ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr); 454527b2640SStefano Zampini ierr = MatGetDiagonal(A,v);CHKERRQ(ierr); 455527b2640SStefano Zampini ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 456527b2640SStefano Zampini ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr); 457527b2640SStefano Zampini for (i=0;i<n;i++) if (array[i] == 0.) break; 458527b2640SStefano Zampini ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr); 459527b2640SStefano Zampini ierr = VecDestroy(&v);CHKERRQ(ierr); 460527b2640SStefano Zampini if (i != n) *missing = PETSC_TRUE; 461527b2640SStefano Zampini if (d) { 462527b2640SStefano Zampini *d = -1; 463527b2640SStefano Zampini if (*missing) { 464527b2640SStefano Zampini PetscInt rstart; 465527b2640SStefano Zampini ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 466527b2640SStefano Zampini *d = i+rstart; 467527b2640SStefano Zampini } 468527b2640SStefano Zampini } 4696bd84002SStefano Zampini PetscFunctionReturn(0); 4706bd84002SStefano Zampini } 4716bd84002SStefano Zampini 4726bd84002SStefano Zampini #undef __FUNCT__ 473cf0a3239SStefano Zampini #define __FUNCT__ "MatISSetUpSF_IS" 474cf0a3239SStefano Zampini static PetscErrorCode MatISSetUpSF_IS(Mat B) 47528f4e0baSStefano Zampini { 47628f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 47728f4e0baSStefano Zampini const PetscInt *gidxs; 4784f2d7cafSStefano Zampini PetscInt nleaves; 47928f4e0baSStefano Zampini PetscErrorCode ierr; 48028f4e0baSStefano Zampini 48128f4e0baSStefano Zampini PetscFunctionBegin; 4824f2d7cafSStefano Zampini if (matis->sf) PetscFunctionReturn(0); 48328f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 4843bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 4854f2d7cafSStefano Zampini ierr = ISLocalToGlobalMappingGetSize(B->rmap->mapping,&nleaves);CHKERRQ(ierr); 4864f2d7cafSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 4873bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 4884f2d7cafSStefano Zampini ierr = PetscMalloc2(matis->sf->nroots,&matis->sf_rootdata,matis->sf->nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 489a8116848SStefano Zampini if (B->rmap->mapping != B->cmap->mapping) { /* setup SF for columns */ 4903d996552SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(B->cmap->mapping,&nleaves);CHKERRQ(ierr); 491a8116848SStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf);CHKERRQ(ierr); 492a8116848SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 4933d996552SStefano Zampini ierr = PetscSFSetGraphLayout(matis->csf,B->cmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 494a8116848SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 4953d996552SStefano Zampini ierr = PetscMalloc2(matis->csf->nroots,&matis->csf_rootdata,matis->csf->nleaves,&matis->csf_leafdata);CHKERRQ(ierr); 496a8116848SStefano Zampini } else { 497a8116848SStefano Zampini matis->csf = matis->sf; 498a8116848SStefano Zampini matis->csf_leafdata = matis->sf_leafdata; 499a8116848SStefano Zampini matis->csf_rootdata = matis->sf_rootdata; 500a8116848SStefano Zampini } 50128f4e0baSStefano Zampini PetscFunctionReturn(0); 50228f4e0baSStefano Zampini } 5032e1947a5SStefano Zampini 5042e1947a5SStefano Zampini #undef __FUNCT__ 5052e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 506eb82efa4SStefano Zampini /*@ 507a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 508a88811baSStefano Zampini 509a88811baSStefano Zampini Collective on MPI_Comm 510a88811baSStefano Zampini 511a88811baSStefano Zampini Input Parameters: 512a88811baSStefano Zampini + B - the matrix 513a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 514a88811baSStefano Zampini (same value is used for all local rows) 515a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 516a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 517a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 518a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 519a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 520a88811baSStefano Zampini the diagonal entry even if it is zero. 521a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 522a88811baSStefano Zampini submatrix (same value is used for all local rows). 523a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 524a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 525a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 526a88811baSStefano Zampini structure. The size of this array is equal to the number 527a88811baSStefano Zampini of local rows, i.e 'm'. 528a88811baSStefano Zampini 529a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 530a88811baSStefano Zampini 531a88811baSStefano Zampini Level: intermediate 532a88811baSStefano Zampini 533a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 534a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 535a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 536a88811baSStefano Zampini 537a88811baSStefano Zampini .keywords: matrix 538a88811baSStefano Zampini 5393c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS 540a88811baSStefano Zampini @*/ 5412e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 5422e1947a5SStefano Zampini { 5432e1947a5SStefano Zampini PetscErrorCode ierr; 5442e1947a5SStefano Zampini 5452e1947a5SStefano Zampini PetscFunctionBegin; 5462e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 5472e1947a5SStefano Zampini PetscValidType(B,1); 5482e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 5492e1947a5SStefano Zampini PetscFunctionReturn(0); 5502e1947a5SStefano Zampini } 5512e1947a5SStefano Zampini 5522e1947a5SStefano Zampini #undef __FUNCT__ 5532e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 5547230de76SStefano Zampini static PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 5552e1947a5SStefano Zampini { 5562e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 55728f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 5582e1947a5SStefano Zampini PetscErrorCode ierr; 5592e1947a5SStefano Zampini 5602e1947a5SStefano Zampini PetscFunctionBegin; 5616c4ed002SBarry Smith if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 562cf0a3239SStefano Zampini ierr = MatISSetUpSF(B);CHKERRQ(ierr); 5634f2d7cafSStefano Zampini 5644f2d7cafSStefano Zampini if (!d_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nz; 5654f2d7cafSStefano Zampini else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 5664f2d7cafSStefano Zampini 5674f2d7cafSStefano Zampini if (!o_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nz; 5684f2d7cafSStefano Zampini else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 5694f2d7cafSStefano Zampini 57028f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 57128f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 57228f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 57328f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 5744f2d7cafSStefano Zampini 5754f2d7cafSStefano Zampini for (i=0;i<matis->sf->nleaves;i++) matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 57628f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 5774f2d7cafSStefano Zampini 5784f2d7cafSStefano Zampini for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 57928f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 5804f2d7cafSStefano Zampini 5814f2d7cafSStefano Zampini for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 58228f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 5832e1947a5SStefano Zampini PetscFunctionReturn(0); 5842e1947a5SStefano Zampini } 585b4319ba4SBarry Smith 586b4319ba4SBarry Smith #undef __FUNCT__ 5873927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private" 5883927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 5893927de2eSStefano Zampini { 5903927de2eSStefano Zampini Mat_IS *matis = (Mat_IS*)(A->data); 5913927de2eSStefano Zampini PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 592ecf5a873SStefano Zampini const PetscInt *global_indices_r,*global_indices_c; 5933927de2eSStefano Zampini PetscInt i,j,bs,rows,cols; 5943927de2eSStefano Zampini PetscInt lrows,lcols; 5953927de2eSStefano Zampini PetscInt local_rows,local_cols; 5963927de2eSStefano Zampini PetscMPIInt nsubdomains; 5973927de2eSStefano Zampini PetscBool isdense,issbaij; 5983927de2eSStefano Zampini PetscErrorCode ierr; 5993927de2eSStefano Zampini 6003927de2eSStefano Zampini PetscFunctionBegin; 6013927de2eSStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr); 6023927de2eSStefano Zampini ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr); 6033927de2eSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 6043927de2eSStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 6053927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 6063927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 607ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 608ecf5a873SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 6097230de76SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 610ecf5a873SStefano Zampini } else { 611ecf5a873SStefano Zampini global_indices_c = global_indices_r; 612ecf5a873SStefano Zampini } 613ecf5a873SStefano Zampini 6143927de2eSStefano Zampini if (issbaij) { 6153927de2eSStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 6163927de2eSStefano Zampini } 6173927de2eSStefano Zampini /* 618ecf5a873SStefano Zampini An SF reduce is needed to sum up properly on shared rows. 6193927de2eSStefano Zampini Note that generally preallocation is not exact, since it overestimates nonzeros 6203927de2eSStefano Zampini */ 621cf0a3239SStefano Zampini ierr = MatISSetUpSF(A);CHKERRQ(ierr); 6223927de2eSStefano Zampini ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr); 6233927de2eSStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr); 6243927de2eSStefano Zampini /* All processes need to compute entire row ownership */ 6253927de2eSStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 6263927de2eSStefano Zampini ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 6273927de2eSStefano Zampini for (i=0;i<nsubdomains;i++) { 6283927de2eSStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 6293927de2eSStefano Zampini row_ownership[j] = i; 6303927de2eSStefano Zampini } 6313927de2eSStefano Zampini } 6327230de76SStefano Zampini ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 6333927de2eSStefano Zampini 6343927de2eSStefano Zampini /* 6353927de2eSStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 6363927de2eSStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 6373927de2eSStefano Zampini */ 6383927de2eSStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 6393927de2eSStefano Zampini /* preallocation as a MATAIJ */ 6403927de2eSStefano Zampini if (isdense) { /* special case for dense local matrices */ 6413927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 642ecf5a873SStefano Zampini PetscInt index_row = global_indices_r[i]; 6433927de2eSStefano Zampini for (j=i;j<local_rows;j++) { 6443927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 645ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[j]; 6463927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 6473927de2eSStefano Zampini my_dnz[i] += 1; 6483927de2eSStefano Zampini } else { /* offdiag block */ 6493927de2eSStefano Zampini my_onz[i] += 1; 6503927de2eSStefano Zampini } 6513927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 6523927de2eSStefano Zampini if (i != j) { 6533927de2eSStefano Zampini owner = row_ownership[index_col]; 6543927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 6553927de2eSStefano Zampini my_dnz[j] += 1; 6563927de2eSStefano Zampini } else { 6573927de2eSStefano Zampini my_onz[j] += 1; 6583927de2eSStefano Zampini } 6593927de2eSStefano Zampini } 6603927de2eSStefano Zampini } 6613927de2eSStefano Zampini } 662bb1015c3SStefano Zampini } else if (matis->A->ops->getrowij) { 663bb1015c3SStefano Zampini const PetscInt *ii,*jj,*jptr; 664bb1015c3SStefano Zampini PetscBool done; 665bb1015c3SStefano Zampini ierr = MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr); 666bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__); 667bb1015c3SStefano Zampini jptr = jj; 668bb1015c3SStefano Zampini for (i=0;i<local_rows;i++) { 669bb1015c3SStefano Zampini PetscInt index_row = global_indices_r[i]; 670bb1015c3SStefano Zampini for (j=0;j<ii[i+1]-ii[i];j++,jptr++) { 671bb1015c3SStefano Zampini PetscInt owner = row_ownership[index_row]; 672bb1015c3SStefano Zampini PetscInt index_col = global_indices_c[*jptr]; 673bb1015c3SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 674bb1015c3SStefano Zampini my_dnz[i] += 1; 675bb1015c3SStefano Zampini } else { /* offdiag block */ 676bb1015c3SStefano Zampini my_onz[i] += 1; 677bb1015c3SStefano Zampini } 678bb1015c3SStefano Zampini /* same as before, interchanging rows and cols */ 679bb1015c3SStefano Zampini if (issbaij && index_col != index_row) { 680bb1015c3SStefano Zampini owner = row_ownership[index_col]; 681bb1015c3SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 682bb1015c3SStefano Zampini my_dnz[*jptr] += 1; 683bb1015c3SStefano Zampini } else { 684bb1015c3SStefano Zampini my_onz[*jptr] += 1; 685bb1015c3SStefano Zampini } 686bb1015c3SStefano Zampini } 687bb1015c3SStefano Zampini } 688bb1015c3SStefano Zampini } 689bb1015c3SStefano Zampini ierr = MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr); 690bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__); 691bb1015c3SStefano Zampini } else { /* loop over rows and use MatGetRow */ 6923927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 6933927de2eSStefano Zampini const PetscInt *cols; 694ecf5a873SStefano Zampini PetscInt ncols,index_row = global_indices_r[i]; 6953927de2eSStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 6963927de2eSStefano Zampini for (j=0;j<ncols;j++) { 6973927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 698ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[cols[j]]; 6993927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 7003927de2eSStefano Zampini my_dnz[i] += 1; 7013927de2eSStefano Zampini } else { /* offdiag block */ 7023927de2eSStefano Zampini my_onz[i] += 1; 7033927de2eSStefano Zampini } 7043927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 705d9a9e74cSStefano Zampini if (issbaij && index_col != index_row) { 7063927de2eSStefano Zampini owner = row_ownership[index_col]; 7073927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 708d9a9e74cSStefano Zampini my_dnz[cols[j]] += 1; 7093927de2eSStefano Zampini } else { 710d9a9e74cSStefano Zampini my_onz[cols[j]] += 1; 7113927de2eSStefano Zampini } 7123927de2eSStefano Zampini } 7133927de2eSStefano Zampini } 7143927de2eSStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 7153927de2eSStefano Zampini } 7163927de2eSStefano Zampini } 717ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 718ecf5a873SStefano Zampini if (global_indices_c != global_indices_r) { 7197230de76SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 720ecf5a873SStefano Zampini } 7213927de2eSStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 722ecf5a873SStefano Zampini 723ecf5a873SStefano Zampini /* Reduce my_dnz and my_onz */ 7243927de2eSStefano Zampini if (maxreduce) { 7253927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 7263927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 727bb1015c3SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 7283927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 7293927de2eSStefano Zampini } else { 7303927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 7313927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 732bb1015c3SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 7333927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 7343927de2eSStefano Zampini } 7353927de2eSStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 7363927de2eSStefano Zampini 7373927de2eSStefano Zampini /* Resize preallocation if overestimated */ 7383927de2eSStefano Zampini for (i=0;i<lrows;i++) { 7393927de2eSStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 7403927de2eSStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 7413927de2eSStefano Zampini } 7423927de2eSStefano Zampini /* set preallocation */ 7433927de2eSStefano Zampini ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 7443927de2eSStefano Zampini for (i=0;i<lrows/bs;i++) { 7453927de2eSStefano Zampini dnz[i] = dnz[i*bs]/bs; 7463927de2eSStefano Zampini onz[i] = onz[i*bs]/bs; 7473927de2eSStefano Zampini } 7483927de2eSStefano Zampini ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 7493927de2eSStefano Zampini ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 7503927de2eSStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 7513927de2eSStefano Zampini if (issbaij) { 7523927de2eSStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 7533927de2eSStefano Zampini } 7543927de2eSStefano Zampini PetscFunctionReturn(0); 7553927de2eSStefano Zampini } 7563927de2eSStefano Zampini 7573927de2eSStefano Zampini #undef __FUNCT__ 758b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 7597230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 760b7ce53b6SStefano Zampini { 761b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 762d9a9e74cSStefano Zampini Mat local_mat; 763b7ce53b6SStefano Zampini /* info on mat */ 7643cfa4ea4SStefano Zampini PetscInt bs,rows,cols,lrows,lcols; 765b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 766b9ed4604SStefano Zampini PetscBool isseqdense,isseqsbaij,isseqaij,isseqbaij; 767b9ed4604SStefano Zampini #if defined (PETSC_USE_DEBUG) 768b9ed4604SStefano Zampini PetscBool lb[4],bb[4]; 769b9ed4604SStefano Zampini #endif 7707c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 771b7ce53b6SStefano Zampini /* values insertion */ 772b7ce53b6SStefano Zampini PetscScalar *array; 773b7ce53b6SStefano Zampini /* work */ 774b7ce53b6SStefano Zampini PetscErrorCode ierr; 775b7ce53b6SStefano Zampini 776b7ce53b6SStefano Zampini PetscFunctionBegin; 777b7ce53b6SStefano Zampini /* get info from mat */ 7787c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 7797c03b4e8SStefano Zampini if (nsubdomains == 1) { 7807c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 7817c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 7827c03b4e8SStefano Zampini } else { 7837c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 7847c03b4e8SStefano Zampini } 7857c03b4e8SStefano Zampini PetscFunctionReturn(0); 7867c03b4e8SStefano Zampini } 787b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 788b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 7893cfa4ea4SStefano Zampini ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 790b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 791b9ed4604SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isseqdense);CHKERRQ(ierr); 792686e3a49SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 793b9ed4604SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQBAIJ,&isseqbaij);CHKERRQ(ierr); 794b9ed4604SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr); 795b9ed4604SStefano Zampini if (!isseqdense && !isseqaij && !isseqbaij && !isseqsbaij) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)(matis->A))->type_name); 796b9ed4604SStefano Zampini #if defined (PETSC_USE_DEBUG) 797b9ed4604SStefano Zampini lb[0] = isseqdense; 798b9ed4604SStefano Zampini lb[1] = isseqaij; 799b9ed4604SStefano Zampini lb[2] = isseqbaij; 800b9ed4604SStefano Zampini lb[3] = isseqsbaij; 801b9ed4604SStefano Zampini ierr = MPIU_Allreduce(lb,bb,4,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 802b9ed4604SStefano Zampini if (!bb[0] && !bb[1] && !bb[2] && !bb[3]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local matrices must have the same type"); 803b9ed4604SStefano Zampini #endif 804b7ce53b6SStefano Zampini 805b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 8063927de2eSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 8073cfa4ea4SStefano Zampini ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 8083927de2eSStefano Zampini ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 809b9ed4604SStefano Zampini if (!isseqsbaij) { 810b9ed4604SStefano Zampini ierr = MatSetType(*M,MATAIJ);CHKERRQ(ierr); 811b9ed4604SStefano Zampini } else { 812b9ed4604SStefano Zampini ierr = MatSetType(*M,MATSBAIJ);CHKERRQ(ierr); 813b9ed4604SStefano Zampini } 8143927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 815b7ce53b6SStefano Zampini } else { 8163cfa4ea4SStefano Zampini PetscInt mbs,mrows,mcols,mlrows,mlcols; 817b7ce53b6SStefano Zampini /* some checks */ 818b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 819b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 8203cfa4ea4SStefano Zampini ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 8216c4ed002SBarry Smith if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 8226c4ed002SBarry Smith if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 8236c4ed002SBarry Smith if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 8246c4ed002SBarry Smith if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 8256c4ed002SBarry Smith if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 826b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 827b7ce53b6SStefano Zampini } 828d9a9e74cSStefano Zampini 829b9ed4604SStefano Zampini if (isseqsbaij) { 830d9a9e74cSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 831d9a9e74cSStefano Zampini } else { 832d9a9e74cSStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 833d9a9e74cSStefano Zampini local_mat = matis->A; 834d9a9e74cSStefano Zampini } 835686e3a49SStefano Zampini 836b7ce53b6SStefano Zampini /* Set values */ 837ecf5a873SStefano Zampini ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 838b9ed4604SStefano Zampini if (isseqdense) { /* special case for dense local matrices */ 839ecf5a873SStefano Zampini PetscInt i,*dummy_rows,*dummy_cols; 840ecf5a873SStefano Zampini 841ecf5a873SStefano Zampini if (local_rows != local_cols) { 842ecf5a873SStefano Zampini ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 843ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 844ecf5a873SStefano Zampini for (i=0;i<local_cols;i++) dummy_cols[i] = i; 845ecf5a873SStefano Zampini } else { 846ecf5a873SStefano Zampini ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 847ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 848ecf5a873SStefano Zampini dummy_cols = dummy_rows; 849ecf5a873SStefano Zampini } 850b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 851d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 852ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 853d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 854ecf5a873SStefano Zampini if (dummy_rows != dummy_cols) { 855ecf5a873SStefano Zampini ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 856ecf5a873SStefano Zampini } else { 857ecf5a873SStefano Zampini ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 858ecf5a873SStefano Zampini } 859686e3a49SStefano Zampini } else if (isseqaij) { 860ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 861686e3a49SStefano Zampini PetscBool done; 862686e3a49SStefano Zampini 863d9a9e74cSStefano Zampini ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 864bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__); 865d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 866686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 867ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 868686e3a49SStefano Zampini } 869d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 870bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__); 871d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 872686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 873ecf5a873SStefano Zampini PetscInt i; 874c0962df8SStefano Zampini 875686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 876686e3a49SStefano Zampini PetscInt j; 877ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 878686e3a49SStefano Zampini 879ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 880ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 881ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 882686e3a49SStefano Zampini } 883b7ce53b6SStefano Zampini } 884b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 885d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 886b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 887b9ed4604SStefano Zampini if (isseqdense) { 888b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 889b7ce53b6SStefano Zampini } 890b7ce53b6SStefano Zampini PetscFunctionReturn(0); 891b7ce53b6SStefano Zampini } 892b7ce53b6SStefano Zampini 893b7ce53b6SStefano Zampini #undef __FUNCT__ 894b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 895b7ce53b6SStefano Zampini /*@ 896b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 897b7ce53b6SStefano Zampini 898b7ce53b6SStefano Zampini Input Parameter: 899b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 900b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 901b7ce53b6SStefano Zampini 902b7ce53b6SStefano Zampini Output Parameter: 903b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 904b7ce53b6SStefano Zampini 905b7ce53b6SStefano Zampini Level: developer 906b7ce53b6SStefano Zampini 907eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 908b7ce53b6SStefano Zampini 909b7ce53b6SStefano Zampini .seealso: MATIS 910b7ce53b6SStefano Zampini @*/ 911b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 912b7ce53b6SStefano Zampini { 913b7ce53b6SStefano Zampini PetscErrorCode ierr; 914b7ce53b6SStefano Zampini 915b7ce53b6SStefano Zampini PetscFunctionBegin; 916b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 917b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 918b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 919b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 920b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 921b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 9226c4ed002SBarry Smith if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 923b7ce53b6SStefano Zampini } 924b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 925b7ce53b6SStefano Zampini PetscFunctionReturn(0); 926b7ce53b6SStefano Zampini } 927b7ce53b6SStefano Zampini 928b7ce53b6SStefano Zampini #undef __FUNCT__ 929ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 930ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 931ad6194a2SStefano Zampini { 932ad6194a2SStefano Zampini PetscErrorCode ierr; 933ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 934ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 935ad6194a2SStefano Zampini Mat B,localmat; 936ad6194a2SStefano Zampini 937ad6194a2SStefano Zampini PetscFunctionBegin; 938ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 939ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 940ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 941e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 942ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 943ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 944b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 945ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 946ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 947ad6194a2SStefano Zampini *newmat = B; 948ad6194a2SStefano Zampini PetscFunctionReturn(0); 949ad6194a2SStefano Zampini } 950ad6194a2SStefano Zampini 951ad6194a2SStefano Zampini #undef __FUNCT__ 95269796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 953a8116848SStefano Zampini static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 95469796d55SStefano Zampini { 95569796d55SStefano Zampini PetscErrorCode ierr; 95669796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 95769796d55SStefano Zampini PetscBool local_sym; 95869796d55SStefano Zampini 95969796d55SStefano Zampini PetscFunctionBegin; 96069796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 961b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 96269796d55SStefano Zampini PetscFunctionReturn(0); 96369796d55SStefano Zampini } 96469796d55SStefano Zampini 96569796d55SStefano Zampini #undef __FUNCT__ 96669796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 967a8116848SStefano Zampini static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 96869796d55SStefano Zampini { 96969796d55SStefano Zampini PetscErrorCode ierr; 97069796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 97169796d55SStefano Zampini PetscBool local_sym; 97269796d55SStefano Zampini 97369796d55SStefano Zampini PetscFunctionBegin; 97469796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 975b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 97669796d55SStefano Zampini PetscFunctionReturn(0); 97769796d55SStefano Zampini } 97869796d55SStefano Zampini 97969796d55SStefano Zampini #undef __FUNCT__ 98045471136SStefano Zampini #define __FUNCT__ "MatIsStructurallySymmetric_IS" 98145471136SStefano Zampini static PetscErrorCode MatIsStructurallySymmetric_IS(Mat A,PetscBool *flg) 98245471136SStefano Zampini { 98345471136SStefano Zampini PetscErrorCode ierr; 98445471136SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 98545471136SStefano Zampini PetscBool local_sym; 98645471136SStefano Zampini 98745471136SStefano Zampini PetscFunctionBegin; 98845471136SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 98945471136SStefano Zampini *flg = PETSC_FALSE; 99045471136SStefano Zampini PetscFunctionReturn(0); 99145471136SStefano Zampini } 99245471136SStefano Zampini ierr = MatIsStructurallySymmetric(matis->A,&local_sym);CHKERRQ(ierr); 99345471136SStefano Zampini ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 99445471136SStefano Zampini PetscFunctionReturn(0); 99545471136SStefano Zampini } 99645471136SStefano Zampini 99745471136SStefano Zampini #undef __FUNCT__ 998b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 999a8116848SStefano Zampini static PetscErrorCode MatDestroy_IS(Mat A) 1000b4319ba4SBarry Smith { 1001dfbe8321SBarry Smith PetscErrorCode ierr; 1002b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 1003b4319ba4SBarry Smith 1004b4319ba4SBarry Smith PetscFunctionBegin; 10056bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 1006e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 1007e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 10086bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 10096bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 10103fd1c9e7SStefano Zampini ierr = VecDestroy(&b->counter);CHKERRQ(ierr); 1011a8116848SStefano Zampini ierr = ISDestroy(&b->getsub_ris);CHKERRQ(ierr); 1012a8116848SStefano Zampini ierr = ISDestroy(&b->getsub_cis);CHKERRQ(ierr); 1013a8116848SStefano Zampini if (b->sf != b->csf) { 1014a8116848SStefano Zampini ierr = PetscSFDestroy(&b->csf);CHKERRQ(ierr); 1015a8116848SStefano Zampini ierr = PetscFree2(b->csf_rootdata,b->csf_leafdata);CHKERRQ(ierr); 1016a8116848SStefano Zampini } 101728f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 101828f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 1019bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 1020dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 1021bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 1022b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 1023b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 10242e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 1025cf0a3239SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",NULL);CHKERRQ(ierr); 1026b4319ba4SBarry Smith PetscFunctionReturn(0); 1027b4319ba4SBarry Smith } 1028b4319ba4SBarry Smith 1029b4319ba4SBarry Smith #undef __FUNCT__ 1030b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 1031a8116848SStefano Zampini static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 1032b4319ba4SBarry Smith { 1033dfbe8321SBarry Smith PetscErrorCode ierr; 1034b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1035b4319ba4SBarry Smith PetscScalar zero = 0.0; 1036b4319ba4SBarry Smith 1037b4319ba4SBarry Smith PetscFunctionBegin; 1038b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 1039e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1040e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1041b4319ba4SBarry Smith 1042b4319ba4SBarry Smith /* multiply the local matrix */ 1043b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 1044b4319ba4SBarry Smith 1045b4319ba4SBarry Smith /* scatter product back into global memory */ 10462dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 1047e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1048e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1049b4319ba4SBarry Smith PetscFunctionReturn(0); 1050b4319ba4SBarry Smith } 1051b4319ba4SBarry Smith 1052b4319ba4SBarry Smith #undef __FUNCT__ 10532e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 1054a8116848SStefano Zampini static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 10552e74eeadSLisandro Dalcin { 1056650997f4SStefano Zampini Vec temp_vec; 10572e74eeadSLisandro Dalcin PetscErrorCode ierr; 10582e74eeadSLisandro Dalcin 10592e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 1060650997f4SStefano Zampini if (v3 != v2) { 1061650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 1062650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 1063650997f4SStefano Zampini } else { 1064650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 1065650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 1066650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 1067650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 1068650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 1069650997f4SStefano Zampini } 10702e74eeadSLisandro Dalcin PetscFunctionReturn(0); 10712e74eeadSLisandro Dalcin } 10722e74eeadSLisandro Dalcin 10732e74eeadSLisandro Dalcin #undef __FUNCT__ 10742e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 1075a8116848SStefano Zampini static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 10762e74eeadSLisandro Dalcin { 10772e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 10782e74eeadSLisandro Dalcin PetscErrorCode ierr; 10792e74eeadSLisandro Dalcin 1080e176bc59SStefano Zampini PetscFunctionBegin; 10812e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 1082e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1083e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10842e74eeadSLisandro Dalcin 10852e74eeadSLisandro Dalcin /* multiply the local matrix */ 1086e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 10872e74eeadSLisandro Dalcin 10882e74eeadSLisandro Dalcin /* scatter product back into global vector */ 1089e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 1090e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1091e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10922e74eeadSLisandro Dalcin PetscFunctionReturn(0); 10932e74eeadSLisandro Dalcin } 10942e74eeadSLisandro Dalcin 10952e74eeadSLisandro Dalcin #undef __FUNCT__ 10962e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 1097a8116848SStefano Zampini static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 10982e74eeadSLisandro Dalcin { 1099650997f4SStefano Zampini Vec temp_vec; 11002e74eeadSLisandro Dalcin PetscErrorCode ierr; 11012e74eeadSLisandro Dalcin 11022e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 1103650997f4SStefano Zampini if (v3 != v2) { 1104650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 1105650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 1106650997f4SStefano Zampini } else { 1107650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 1108650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 1109650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 1110650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 1111650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 1112650997f4SStefano Zampini } 11132e74eeadSLisandro Dalcin PetscFunctionReturn(0); 11142e74eeadSLisandro Dalcin } 11152e74eeadSLisandro Dalcin 11162e74eeadSLisandro Dalcin #undef __FUNCT__ 1117b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 1118a8116848SStefano Zampini static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 1119b4319ba4SBarry Smith { 1120b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 1121dfbe8321SBarry Smith PetscErrorCode ierr; 1122b4319ba4SBarry Smith PetscViewer sviewer; 1123ee2491ecSStefano Zampini PetscBool isascii,view = PETSC_TRUE; 1124b4319ba4SBarry Smith 1125b4319ba4SBarry Smith PetscFunctionBegin; 1126ee2491ecSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 1127ee2491ecSStefano Zampini if (isascii) { 1128ee2491ecSStefano Zampini PetscViewerFormat format; 1129ee2491ecSStefano Zampini 1130ee2491ecSStefano Zampini ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1131ee2491ecSStefano Zampini if (format == PETSC_VIEWER_ASCII_INFO) view = PETSC_FALSE; 1132ee2491ecSStefano Zampini } 1133ee2491ecSStefano Zampini if (!view) PetscFunctionReturn(0); 11343f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 1135b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 11363f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 11376e520ac8SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1138b4319ba4SBarry Smith PetscFunctionReturn(0); 1139b4319ba4SBarry Smith } 1140b4319ba4SBarry Smith 1141b4319ba4SBarry Smith #undef __FUNCT__ 1142b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 1143a8116848SStefano Zampini static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 1144b4319ba4SBarry Smith { 1145dfbe8321SBarry Smith PetscErrorCode ierr; 1146e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 1147b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1148b4319ba4SBarry Smith IS from,to; 1149e176bc59SStefano Zampini Vec cglobal,rglobal; 1150b4319ba4SBarry Smith 1151b4319ba4SBarry Smith PetscFunctionBegin; 1152784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 1153e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 11543bbff08aSStefano Zampini /* Destroy any previous data */ 115570cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 115670cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 11573fd1c9e7SStefano Zampini ierr = VecDestroy(&is->counter);CHKERRQ(ierr); 1158e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 1159e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 11601c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 116128f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 116228f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 11633bbff08aSStefano Zampini 11643bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 1165fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1166fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1167fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 1168fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 1169b4319ba4SBarry Smith 1170b4319ba4SBarry Smith /* Create the local matrix A */ 1171e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 1172e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 1173e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 1174e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 1175f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 1176e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 1177e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 1178e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 1179ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 1180ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 1181b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 1182b4319ba4SBarry Smith 1183f26d0771SStefano Zampini if (!is->islocalref) { /* setup scatters and local vectors for MatMult */ 1184b4319ba4SBarry Smith /* Create the local work vectors */ 11853bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 1186b4319ba4SBarry Smith 1187e176bc59SStefano Zampini /* setup the global to local scatters */ 1188e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 1189e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 1190784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 1191e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 1192e176bc59SStefano Zampini if (rmapping != cmapping) { 1193e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 1194e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 1195e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 1196e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 1197e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 1198e176bc59SStefano Zampini } else { 1199e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 1200e176bc59SStefano Zampini is->cctx = is->rctx; 1201e176bc59SStefano Zampini } 12023fd1c9e7SStefano Zampini 12033fd1c9e7SStefano Zampini /* interface counter vector (local) */ 12043fd1c9e7SStefano Zampini ierr = VecDuplicate(is->y,&is->counter);CHKERRQ(ierr); 12053fd1c9e7SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 12063fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12073fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12083fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12093fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12103fd1c9e7SStefano Zampini 12113fd1c9e7SStefano Zampini /* free workspace */ 1212e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 1213e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 12146bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 12156bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1216f26d0771SStefano Zampini } 12179c0b00dbSStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 1218b4319ba4SBarry Smith PetscFunctionReturn(0); 1219b4319ba4SBarry Smith } 1220b4319ba4SBarry Smith 12212e74eeadSLisandro Dalcin #undef __FUNCT__ 12222e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 1223a8116848SStefano Zampini static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 12242e74eeadSLisandro Dalcin { 12252e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 12262e74eeadSLisandro Dalcin PetscErrorCode ierr; 122797563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 122897563a80SStefano Zampini PetscInt i,zm,zn; 122997563a80SStefano Zampini #endif 1230f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 12312e74eeadSLisandro Dalcin 12322e74eeadSLisandro Dalcin PetscFunctionBegin; 12332e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 1234f26d0771SStefano Zampini if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n); 123597563a80SStefano Zampini /* count negative indices */ 123697563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 123797563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 12382e74eeadSLisandro Dalcin #endif 123997563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 124097563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 124197563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 124297563a80SStefano Zampini /* count negative indices (should be the same as before) */ 124397563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 124497563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 124597563a80SStefano 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"); 124697563a80SStefano 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"); 124797563a80SStefano Zampini #endif 12482e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 12492e74eeadSLisandro Dalcin PetscFunctionReturn(0); 12502e74eeadSLisandro Dalcin } 12512e74eeadSLisandro Dalcin 1252b4319ba4SBarry Smith #undef __FUNCT__ 125397563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS" 1254a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 125597563a80SStefano Zampini { 125697563a80SStefano Zampini Mat_IS *is = (Mat_IS*)mat->data; 125797563a80SStefano Zampini PetscErrorCode ierr; 125897563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 125997563a80SStefano Zampini PetscInt i,zm,zn; 126097563a80SStefano Zampini #endif 1261f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 126297563a80SStefano Zampini 126397563a80SStefano Zampini PetscFunctionBegin; 126497563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 1265f26d0771SStefano Zampini if (m > MATIS_MAX_ENTRIES_INSERTION || n > MATIS_MAX_ENTRIES_INSERTION) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column block indices must be <= %D: they are %D %D",MATIS_MAX_ENTRIES_INSERTION,m,n); 126697563a80SStefano Zampini /* count negative indices */ 126797563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 126897563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 126997563a80SStefano Zampini #endif 127097563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 127197563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 127297563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 127397563a80SStefano Zampini /* count negative indices (should be the same as before) */ 127497563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 127597563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 127697563a80SStefano 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"); 127797563a80SStefano 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"); 127897563a80SStefano Zampini #endif 127997563a80SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 128097563a80SStefano Zampini PetscFunctionReturn(0); 128197563a80SStefano Zampini } 128297563a80SStefano Zampini 128397563a80SStefano Zampini #undef __FUNCT__ 1284b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 1285a8116848SStefano Zampini static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1286b4319ba4SBarry Smith { 1287dfbe8321SBarry Smith PetscErrorCode ierr; 1288b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1289b4319ba4SBarry Smith 1290b4319ba4SBarry Smith PetscFunctionBegin; 1291b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1292b4319ba4SBarry Smith PetscFunctionReturn(0); 1293b4319ba4SBarry Smith } 1294b4319ba4SBarry Smith 1295b4319ba4SBarry Smith #undef __FUNCT__ 1296f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 1297a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1298f0006bf2SLisandro Dalcin { 1299f0006bf2SLisandro Dalcin PetscErrorCode ierr; 1300f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 1301f0006bf2SLisandro Dalcin 1302f0006bf2SLisandro Dalcin PetscFunctionBegin; 1303f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1304f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 1305f0006bf2SLisandro Dalcin } 1306f0006bf2SLisandro Dalcin 1307f0006bf2SLisandro Dalcin #undef __FUNCT__ 1308f0ae7da4SStefano Zampini #define __FUNCT__ "MatISZeroRowsColumnsLocal_Private" 1309f0ae7da4SStefano Zampini static PetscErrorCode MatISZeroRowsColumnsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,PetscBool columns) 13102e74eeadSLisandro Dalcin { 1311f0ae7da4SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 13122e74eeadSLisandro Dalcin PetscErrorCode ierr; 13132e74eeadSLisandro Dalcin 13142e74eeadSLisandro Dalcin PetscFunctionBegin; 1315f0ae7da4SStefano Zampini if (!n) { 1316f0ae7da4SStefano Zampini is->pure_neumann = PETSC_TRUE; 1317f0ae7da4SStefano Zampini } else { 1318f0ae7da4SStefano Zampini PetscInt i; 1319f0ae7da4SStefano Zampini is->pure_neumann = PETSC_FALSE; 1320f0ae7da4SStefano Zampini 1321f0ae7da4SStefano Zampini if (columns) { 1322f0ae7da4SStefano Zampini ierr = MatZeroRowsColumns(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 1323f0ae7da4SStefano Zampini } else { 1324f0ae7da4SStefano Zampini ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 13252e74eeadSLisandro Dalcin } 1326f0ae7da4SStefano Zampini if (diag != 0.) { 1327f0ae7da4SStefano Zampini const PetscScalar *array; 1328f0ae7da4SStefano Zampini ierr = VecGetArrayRead(is->counter,&array);CHKERRQ(ierr); 1329f0ae7da4SStefano Zampini for (i=0; i<n; i++) { 1330f0ae7da4SStefano Zampini ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 1331f0ae7da4SStefano Zampini } 1332f0ae7da4SStefano Zampini ierr = VecRestoreArrayRead(is->counter,&array);CHKERRQ(ierr); 1333f0ae7da4SStefano Zampini } 1334f0ae7da4SStefano Zampini ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1335f0ae7da4SStefano Zampini ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1336f0ae7da4SStefano Zampini } 13372e74eeadSLisandro Dalcin PetscFunctionReturn(0); 13382e74eeadSLisandro Dalcin } 13392e74eeadSLisandro Dalcin 13402e74eeadSLisandro Dalcin #undef __FUNCT__ 1341f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_Private_IS" 1342f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Private_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b,PetscBool columns) 1343b4319ba4SBarry Smith { 13446e520ac8SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 13456e520ac8SStefano Zampini PetscInt nr,nl,len,i; 13466e520ac8SStefano Zampini PetscInt *lrows; 1347dfbe8321SBarry Smith PetscErrorCode ierr; 1348b4319ba4SBarry Smith 1349b4319ba4SBarry Smith PetscFunctionBegin; 1350f0ae7da4SStefano Zampini #if defined(PETSC_USE_DEBUG) 1351f0ae7da4SStefano Zampini if (columns || diag != 0. || (x && b)) { 1352f0ae7da4SStefano Zampini PetscBool cong; 1353f0ae7da4SStefano Zampini ierr = PetscLayoutCompare(A->rmap,A->cmap,&cong);CHKERRQ(ierr); 1354f0ae7da4SStefano Zampini if (!cong && columns) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Columns can be zeroed if and only if A->rmap and A->cmap are congruent for MATIS"); 1355f0ae7da4SStefano Zampini if (!cong && diag != 0.) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Nonzero diagonal value supported if and only if A->rmap and A->cmap are congruent for MATIS"); 1356f0ae7da4SStefano Zampini if (!cong && x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"A->rmap and A->cmap need to be congruent"); 1357b4319ba4SBarry Smith } 1358f0ae7da4SStefano Zampini #endif 13596e520ac8SStefano Zampini /* get locally owned rows */ 1360f0ae7da4SStefano Zampini ierr = PetscLayoutMapLocal_Private(A->rmap,n,rows,&len,&lrows,NULL);CHKERRQ(ierr); 13616e520ac8SStefano Zampini /* fix right hand side if needed */ 13626e520ac8SStefano Zampini if (x && b) { 13636e520ac8SStefano Zampini const PetscScalar *xx; 13646e520ac8SStefano Zampini PetscScalar *bb; 13652205254eSKarl Rupp 13666e520ac8SStefano Zampini ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr); 13676e520ac8SStefano Zampini ierr = VecGetArray(b, &bb);CHKERRQ(ierr); 13686e520ac8SStefano Zampini for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]]; 13696e520ac8SStefano Zampini ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr); 13706e520ac8SStefano Zampini ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr); 1371b4319ba4SBarry Smith } 13726e520ac8SStefano Zampini /* get rows associated to the local matrices */ 13733d996552SStefano Zampini ierr = MatISSetUpSF(A);CHKERRQ(ierr); 13746e520ac8SStefano Zampini ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr); 13756e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr); 13766e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 13776e520ac8SStefano Zampini for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1; 13786e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 13796e520ac8SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 13806e520ac8SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 13816e520ac8SStefano Zampini ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr); 13826e520ac8SStefano Zampini for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i; 1383f0ae7da4SStefano Zampini ierr = MatISZeroRowsColumnsLocal_Private(A,nr,lrows,diag,columns);CHKERRQ(ierr); 13846e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 1385b4319ba4SBarry Smith PetscFunctionReturn(0); 1386b4319ba4SBarry Smith } 1387b4319ba4SBarry Smith 1388b4319ba4SBarry Smith #undef __FUNCT__ 1389f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRows_IS" 1390f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 1391b4319ba4SBarry Smith { 1392b4319ba4SBarry Smith PetscErrorCode ierr; 1393b4319ba4SBarry Smith 1394b4319ba4SBarry Smith PetscFunctionBegin; 1395f0ae7da4SStefano Zampini ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_FALSE);CHKERRQ(ierr); 1396f0ae7da4SStefano Zampini PetscFunctionReturn(0); 1397f0ae7da4SStefano Zampini } 1398b4319ba4SBarry Smith 1399f0ae7da4SStefano Zampini #undef __FUNCT__ 1400f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_IS" 1401f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 1402f0ae7da4SStefano Zampini { 1403f0ae7da4SStefano Zampini PetscErrorCode ierr; 1404f0ae7da4SStefano Zampini 1405f0ae7da4SStefano Zampini PetscFunctionBegin; 1406f0ae7da4SStefano Zampini ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_TRUE);CHKERRQ(ierr); 1407b4319ba4SBarry Smith PetscFunctionReturn(0); 1408b4319ba4SBarry Smith } 1409b4319ba4SBarry Smith 1410b4319ba4SBarry Smith #undef __FUNCT__ 1411b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 1412a8116848SStefano Zampini static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 1413b4319ba4SBarry Smith { 1414b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1415dfbe8321SBarry Smith PetscErrorCode ierr; 1416b4319ba4SBarry Smith 1417b4319ba4SBarry Smith PetscFunctionBegin; 1418b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 1419b4319ba4SBarry Smith PetscFunctionReturn(0); 1420b4319ba4SBarry Smith } 1421b4319ba4SBarry Smith 1422b4319ba4SBarry Smith #undef __FUNCT__ 1423b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 1424a8116848SStefano Zampini static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 1425b4319ba4SBarry Smith { 1426b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1427dfbe8321SBarry Smith PetscErrorCode ierr; 1428b4319ba4SBarry Smith 1429b4319ba4SBarry Smith PetscFunctionBegin; 1430b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 1431b4319ba4SBarry Smith PetscFunctionReturn(0); 1432b4319ba4SBarry Smith } 1433b4319ba4SBarry Smith 1434b4319ba4SBarry Smith #undef __FUNCT__ 1435b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 1436a8116848SStefano Zampini static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 1437b4319ba4SBarry Smith { 1438b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 1439b4319ba4SBarry Smith 1440b4319ba4SBarry Smith PetscFunctionBegin; 1441b4319ba4SBarry Smith *local = is->A; 1442b4319ba4SBarry Smith PetscFunctionReturn(0); 1443b4319ba4SBarry Smith } 1444b4319ba4SBarry Smith 1445b4319ba4SBarry Smith #undef __FUNCT__ 1446b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 1447b4319ba4SBarry Smith /*@ 1448b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 1449b4319ba4SBarry Smith 1450b4319ba4SBarry Smith Input Parameter: 1451b4319ba4SBarry Smith . mat - the matrix 1452b4319ba4SBarry Smith 1453b4319ba4SBarry Smith Output Parameter: 1454eb82efa4SStefano Zampini . local - the local matrix 1455b4319ba4SBarry Smith 1456b4319ba4SBarry Smith Level: advanced 1457b4319ba4SBarry Smith 1458b4319ba4SBarry Smith Notes: 1459b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 1460b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 1461b4319ba4SBarry Smith of the MatSetValues() operation. 1462b4319ba4SBarry Smith 1463b4319ba4SBarry Smith .seealso: MATIS 1464b4319ba4SBarry Smith @*/ 14657087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 1466b4319ba4SBarry Smith { 14674ac538c5SBarry Smith PetscErrorCode ierr; 1468b4319ba4SBarry Smith 1469b4319ba4SBarry Smith PetscFunctionBegin; 14700700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1471b4319ba4SBarry Smith PetscValidPointer(local,2); 14724ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 1473b4319ba4SBarry Smith PetscFunctionReturn(0); 1474b4319ba4SBarry Smith } 1475b4319ba4SBarry Smith 14763b03a366Sstefano_zampini #undef __FUNCT__ 14773b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 1478a8116848SStefano Zampini static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 14793b03a366Sstefano_zampini { 14803b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 14813b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 14823b03a366Sstefano_zampini PetscErrorCode ierr; 14833b03a366Sstefano_zampini 14843b03a366Sstefano_zampini PetscFunctionBegin; 14854e4c7dbeSStefano Zampini if (is->A) { 14863b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 14873b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 1488f0ae7da4SStefano Zampini if (orows != nrows || ocols != ncols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local MATIS matrix should be of size %Dx%D (you passed a %Dx%D matrix)",orows,ocols,nrows,ncols); 14894e4c7dbeSStefano Zampini } 14903b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 14913b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 14923b03a366Sstefano_zampini is->A = local; 14933b03a366Sstefano_zampini PetscFunctionReturn(0); 14943b03a366Sstefano_zampini } 14953b03a366Sstefano_zampini 14963b03a366Sstefano_zampini #undef __FUNCT__ 14973b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 14983b03a366Sstefano_zampini /*@ 1499eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 15003b03a366Sstefano_zampini 15013b03a366Sstefano_zampini Input Parameter: 15023b03a366Sstefano_zampini . mat - the matrix 1503eb82efa4SStefano Zampini . local - the local matrix 15043b03a366Sstefano_zampini 15053b03a366Sstefano_zampini Output Parameter: 15063b03a366Sstefano_zampini 15073b03a366Sstefano_zampini Level: advanced 15083b03a366Sstefano_zampini 15093b03a366Sstefano_zampini Notes: 15103b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 15113b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 15123b03a366Sstefano_zampini 15133b03a366Sstefano_zampini .seealso: MATIS 15143b03a366Sstefano_zampini @*/ 15153b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 15163b03a366Sstefano_zampini { 15173b03a366Sstefano_zampini PetscErrorCode ierr; 15183b03a366Sstefano_zampini 15193b03a366Sstefano_zampini PetscFunctionBegin; 15203b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1521b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 15223b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 15233b03a366Sstefano_zampini PetscFunctionReturn(0); 15243b03a366Sstefano_zampini } 15253b03a366Sstefano_zampini 15266726f965SBarry Smith #undef __FUNCT__ 15276726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 1528a8116848SStefano Zampini static PetscErrorCode MatZeroEntries_IS(Mat A) 15296726f965SBarry Smith { 15306726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 15316726f965SBarry Smith PetscErrorCode ierr; 15326726f965SBarry Smith 15336726f965SBarry Smith PetscFunctionBegin; 15346726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 15356726f965SBarry Smith PetscFunctionReturn(0); 15366726f965SBarry Smith } 15376726f965SBarry Smith 15386726f965SBarry Smith #undef __FUNCT__ 15392e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 1540a8116848SStefano Zampini static PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 15412e74eeadSLisandro Dalcin { 15422e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 15432e74eeadSLisandro Dalcin PetscErrorCode ierr; 15442e74eeadSLisandro Dalcin 15452e74eeadSLisandro Dalcin PetscFunctionBegin; 15462e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 15472e74eeadSLisandro Dalcin PetscFunctionReturn(0); 15482e74eeadSLisandro Dalcin } 15492e74eeadSLisandro Dalcin 15502e74eeadSLisandro Dalcin #undef __FUNCT__ 15512e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 1552a8116848SStefano Zampini static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 15532e74eeadSLisandro Dalcin { 15542e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 15552e74eeadSLisandro Dalcin PetscErrorCode ierr; 15562e74eeadSLisandro Dalcin 15572e74eeadSLisandro Dalcin PetscFunctionBegin; 15582e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 1559e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 15602e74eeadSLisandro Dalcin 15612e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 15622e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 1563e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1564e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15652e74eeadSLisandro Dalcin PetscFunctionReturn(0); 15662e74eeadSLisandro Dalcin } 15672e74eeadSLisandro Dalcin 15682e74eeadSLisandro Dalcin #undef __FUNCT__ 15696726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 1570a8116848SStefano Zampini static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 15716726f965SBarry Smith { 15726726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 15736726f965SBarry Smith PetscErrorCode ierr; 15746726f965SBarry Smith 15756726f965SBarry Smith PetscFunctionBegin; 15764e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 15776726f965SBarry Smith PetscFunctionReturn(0); 15786726f965SBarry Smith } 15796726f965SBarry Smith 1580284134d9SBarry Smith #undef __FUNCT__ 1581f26d0771SStefano Zampini #define __FUNCT__ "MatAXPY_IS" 1582f26d0771SStefano Zampini static PetscErrorCode MatAXPY_IS(Mat Y,PetscScalar a,Mat X,MatStructure str) 1583f26d0771SStefano Zampini { 1584f26d0771SStefano Zampini Mat_IS *y = (Mat_IS*)Y->data; 1585f26d0771SStefano Zampini Mat_IS *x; 1586f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1587f26d0771SStefano Zampini PetscBool ismatis; 1588f26d0771SStefano Zampini #endif 1589f26d0771SStefano Zampini PetscErrorCode ierr; 1590f26d0771SStefano Zampini 1591f26d0771SStefano Zampini PetscFunctionBegin; 1592f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1593f26d0771SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)X,MATIS,&ismatis);CHKERRQ(ierr); 1594f26d0771SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS"); 1595f26d0771SStefano Zampini #endif 1596f26d0771SStefano Zampini x = (Mat_IS*)X->data; 1597f26d0771SStefano Zampini ierr = MatAXPY(y->A,a,x->A,str);CHKERRQ(ierr); 1598f26d0771SStefano Zampini PetscFunctionReturn(0); 1599f26d0771SStefano Zampini } 1600f26d0771SStefano Zampini 1601f26d0771SStefano Zampini #undef __FUNCT__ 1602f26d0771SStefano Zampini #define __FUNCT__ "MatGetLocalSubMatrix_IS" 1603f26d0771SStefano Zampini static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat) 1604f26d0771SStefano Zampini { 1605f26d0771SStefano Zampini Mat lA; 1606f26d0771SStefano Zampini Mat_IS *matis; 1607f26d0771SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 1608f26d0771SStefano Zampini IS is; 1609f26d0771SStefano Zampini const PetscInt *rg,*rl; 1610f26d0771SStefano Zampini PetscInt nrg; 1611f26d0771SStefano Zampini PetscInt N,M,nrl,i,*idxs; 1612f26d0771SStefano Zampini PetscErrorCode ierr; 1613f26d0771SStefano Zampini 1614f26d0771SStefano Zampini PetscFunctionBegin; 1615f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 1616f26d0771SStefano Zampini ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr); 1617f26d0771SStefano Zampini ierr = ISGetIndices(row,&rl);CHKERRQ(ierr); 1618f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr); 1619f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1620f0ae7da4SStefano Zampini for (i=0;i<nrl;i++) if (rl[i]>=nrg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local row index %D -> %D greater than maximum possible %D",i,rl[i],nrg); 1621f26d0771SStefano Zampini #endif 1622f26d0771SStefano Zampini ierr = PetscMalloc1(nrg,&idxs);CHKERRQ(ierr); 1623f26d0771SStefano Zampini /* map from [0,nrl) to row */ 1624f26d0771SStefano Zampini for (i=0;i<nrl;i++) idxs[i] = rl[i]; 1625f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1626f26d0771SStefano Zampini for (i=nrl;i<nrg;i++) idxs[i] = nrg; 1627f26d0771SStefano Zampini #else 1628f26d0771SStefano Zampini for (i=nrl;i<nrg;i++) idxs[i] = -1; 1629f26d0771SStefano Zampini #endif 1630f26d0771SStefano Zampini ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr); 1631f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 1632f26d0771SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 1633f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 1634f26d0771SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 1635f26d0771SStefano Zampini /* compute new l2g map for columns */ 1636f26d0771SStefano Zampini if (col != row || A->rmap->mapping != A->cmap->mapping) { 1637f26d0771SStefano Zampini const PetscInt *cg,*cl; 1638f26d0771SStefano Zampini PetscInt ncg; 1639f26d0771SStefano Zampini PetscInt ncl; 1640f26d0771SStefano Zampini 1641f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 1642f26d0771SStefano Zampini ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr); 1643f26d0771SStefano Zampini ierr = ISGetIndices(col,&cl);CHKERRQ(ierr); 1644f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr); 1645f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1646f0ae7da4SStefano Zampini for (i=0;i<ncl;i++) if (cl[i]>=ncg) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local column index %D -> %D greater than maximum possible %D",i,cl[i],ncg); 1647f26d0771SStefano Zampini #endif 1648f26d0771SStefano Zampini ierr = PetscMalloc1(ncg,&idxs);CHKERRQ(ierr); 1649f26d0771SStefano Zampini /* map from [0,ncl) to col */ 1650f26d0771SStefano Zampini for (i=0;i<ncl;i++) idxs[i] = cl[i]; 1651f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1652f26d0771SStefano Zampini for (i=ncl;i<ncg;i++) idxs[i] = ncg; 1653f26d0771SStefano Zampini #else 1654f26d0771SStefano Zampini for (i=ncl;i<ncg;i++) idxs[i] = -1; 1655f26d0771SStefano Zampini #endif 1656f26d0771SStefano Zampini ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr); 1657f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 1658f26d0771SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 1659f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 1660f26d0771SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 1661f26d0771SStefano Zampini } else { 1662f26d0771SStefano Zampini ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr); 1663f26d0771SStefano Zampini cl2g = rl2g; 1664f26d0771SStefano Zampini } 1665f26d0771SStefano Zampini /* create the MATIS submatrix */ 1666f26d0771SStefano Zampini ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 1667f26d0771SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr); 1668f26d0771SStefano Zampini ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1669f26d0771SStefano Zampini ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr); 1670b0aa3428SStefano Zampini matis = (Mat_IS*)((*submat)->data); 1671f26d0771SStefano Zampini matis->islocalref = PETSC_TRUE; 1672f26d0771SStefano Zampini ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr); 1673f26d0771SStefano Zampini ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr); 1674f26d0771SStefano Zampini ierr = MatISSetLocalMat(*submat,lA);CHKERRQ(ierr); 1675f26d0771SStefano Zampini ierr = MatSetUp(*submat);CHKERRQ(ierr); 1676f26d0771SStefano Zampini ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1677f26d0771SStefano Zampini ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1678f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 1679f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 1680f26d0771SStefano Zampini /* remove unsupported ops */ 1681f26d0771SStefano Zampini ierr = PetscMemzero((*submat)->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1682f26d0771SStefano Zampini (*submat)->ops->destroy = MatDestroy_IS; 1683f26d0771SStefano Zampini (*submat)->ops->setvalueslocal = MatSetValuesLocal_SubMat_IS; 1684f26d0771SStefano Zampini (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS; 1685f26d0771SStefano Zampini (*submat)->ops->assemblybegin = MatAssemblyBegin_IS; 1686f26d0771SStefano Zampini (*submat)->ops->assemblyend = MatAssemblyEnd_IS; 1687f26d0771SStefano Zampini PetscFunctionReturn(0); 1688f26d0771SStefano Zampini } 1689f26d0771SStefano Zampini 1690f26d0771SStefano Zampini #undef __FUNCT__ 1691284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 1692284134d9SBarry Smith /*@ 16933c212e90SHong Zhang MatCreateIS - Creates a "process" unassembled matrix, assembled on each 1694284134d9SBarry Smith process but not across processes. 1695284134d9SBarry Smith 1696284134d9SBarry Smith Input Parameters: 1697284134d9SBarry Smith + comm - MPI communicator that will share the matrix 1698e176bc59SStefano Zampini . bs - block size of the matrix 1699df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 1700e176bc59SStefano Zampini . rmap - local to global map for rows 1701e176bc59SStefano Zampini - cmap - local to global map for cols 1702284134d9SBarry Smith 1703284134d9SBarry Smith Output Parameter: 1704284134d9SBarry Smith . A - the resulting matrix 1705284134d9SBarry Smith 17068e6c10adSSatish Balay Level: advanced 17078e6c10adSSatish Balay 17083c212e90SHong Zhang Notes: See MATIS for more details. 17093c212e90SHong Zhang m and n are NOT related to the size of the map; they are the size of the part of the vector owned 1710e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 17113c212e90SHong Zhang If either rmap or cmap are NULL, then the matrix is assumed to be square. 1712284134d9SBarry Smith 1713284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 1714284134d9SBarry Smith @*/ 1715e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 1716284134d9SBarry Smith { 1717284134d9SBarry Smith PetscErrorCode ierr; 1718284134d9SBarry Smith 1719284134d9SBarry Smith PetscFunctionBegin; 1720e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 1721284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1722d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 1723284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1724284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 1725e176bc59SStefano Zampini if (rmap && cmap) { 1726e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 1727e176bc59SStefano Zampini } else if (!rmap) { 1728e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 1729e176bc59SStefano Zampini } else { 1730e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 1731e176bc59SStefano Zampini } 1732284134d9SBarry Smith PetscFunctionReturn(0); 1733284134d9SBarry Smith } 1734284134d9SBarry Smith 1735b4319ba4SBarry Smith /*MC 1736f26d0771SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition methods (e.g. PCBDDC or KSPFETIDP). 1737b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 1738b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 1739b4319ba4SBarry Smith product is handled "implicitly". 1740b4319ba4SBarry Smith 1741b4319ba4SBarry Smith Operations Provided: 17426726f965SBarry Smith + MatMult() 17432e74eeadSLisandro Dalcin . MatMultAdd() 17442e74eeadSLisandro Dalcin . MatMultTranspose() 17452e74eeadSLisandro Dalcin . MatMultTransposeAdd() 17466726f965SBarry Smith . MatZeroEntries() 17476726f965SBarry Smith . MatSetOption() 17482e74eeadSLisandro Dalcin . MatZeroRows() 17492e74eeadSLisandro Dalcin . MatSetValues() 175097563a80SStefano Zampini . MatSetValuesBlocked() 17516726f965SBarry Smith . MatSetValuesLocal() 175297563a80SStefano Zampini . MatSetValuesBlockedLocal() 17532e74eeadSLisandro Dalcin . MatScale() 17542e74eeadSLisandro Dalcin . MatGetDiagonal() 17552b404112SStefano Zampini . MatMissingDiagonal() 17562b404112SStefano Zampini . MatDuplicate() 17572b404112SStefano Zampini . MatCopy() 1758f26d0771SStefano Zampini . MatAXPY() 1759f26d0771SStefano Zampini . MatGetSubMatrix() 1760f26d0771SStefano Zampini . MatGetLocalSubMatrix() 17616726f965SBarry Smith - MatSetLocalToGlobalMapping() 1762b4319ba4SBarry Smith 1763b4319ba4SBarry Smith Options Database Keys: 1764b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1765b4319ba4SBarry Smith 1766b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1767b4319ba4SBarry Smith 1768b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1769b4319ba4SBarry Smith 1770b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1771eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1772b4319ba4SBarry Smith 1773b4319ba4SBarry Smith Level: advanced 1774b4319ba4SBarry Smith 1775f26d0771SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), MatCreateIS(), PCBDDC, KSPFETIDP 1776b4319ba4SBarry Smith 1777b4319ba4SBarry Smith M*/ 1778b4319ba4SBarry Smith 1779b4319ba4SBarry Smith #undef __FUNCT__ 1780b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 17818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1782b4319ba4SBarry Smith { 1783dfbe8321SBarry Smith PetscErrorCode ierr; 1784b4319ba4SBarry Smith Mat_IS *b; 1785b4319ba4SBarry Smith 1786b4319ba4SBarry Smith PetscFunctionBegin; 1787b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1788b4319ba4SBarry Smith A->data = (void*)b; 1789b4319ba4SBarry Smith 1790e176bc59SStefano Zampini /* matrix ops */ 1791e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1792b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 17932e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 17942e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 17952e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1796b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1797b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 17982e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 179998921651SStefano Zampini A->ops->setvaluesblocked = MatSetValuesBlocked_IS; 1800b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1801f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 18022e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1803f0ae7da4SStefano Zampini A->ops->zerorowscolumns = MatZeroRowsColumns_IS; 1804b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1805b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1806b4319ba4SBarry Smith A->ops->view = MatView_IS; 18076726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 18082e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 18092e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 18106726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 181169796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 181269796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 181345471136SStefano Zampini A->ops->isstructurallysymmetric = MatIsStructurallySymmetric_IS; 1814ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 18156bd84002SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_IS; 18162b404112SStefano Zampini A->ops->copy = MatCopy_IS; 1817659959c5SStefano Zampini A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_IS; 1818a8116848SStefano Zampini A->ops->getsubmatrix = MatGetSubMatrix_IS; 1819f26d0771SStefano Zampini A->ops->axpy = MatAXPY_IS; 18203fd1c9e7SStefano Zampini A->ops->diagonalset = MatDiagonalSet_IS; 18213fd1c9e7SStefano Zampini A->ops->shift = MatShift_IS; 1822b4319ba4SBarry Smith 1823b7ce53b6SStefano Zampini /* special MATIS functions */ 1824bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1825bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1826b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 18272e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 1828cf0a3239SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",MatISSetUpSF_IS);CHKERRQ(ierr); 182917667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1830b4319ba4SBarry Smith PetscFunctionReturn(0); 1831b4319ba4SBarry Smith } 1832