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*/ 116989cf23SStefano 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__ 176989cf23SStefano Zampini #define __FUNCT__ "MatISContainerDestroy_Private" 186989cf23SStefano Zampini static PetscErrorCode MatISContainerDestroy_Private(void *ptr) 196989cf23SStefano Zampini { 206989cf23SStefano Zampini PetscErrorCode ierr; 216989cf23SStefano Zampini 226989cf23SStefano Zampini PetscFunctionBeginUser; 236989cf23SStefano Zampini ierr = PetscFree(ptr); CHKERRQ(ierr); 246989cf23SStefano Zampini PetscFunctionReturn(0); 256989cf23SStefano Zampini } 266989cf23SStefano Zampini 276989cf23SStefano Zampini #undef __FUNCT__ 286989cf23SStefano Zampini #define __FUNCT__ "MatConvert_MPIAIJ_IS" 296989cf23SStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_IS(Mat A,MatType type,MatReuse reuse,Mat *newmat) 306989cf23SStefano Zampini { 316989cf23SStefano Zampini Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 326989cf23SStefano Zampini Mat_SeqAIJ *diag = (Mat_SeqAIJ*)(aij->A->data); 336989cf23SStefano Zampini Mat_SeqAIJ *offd = (Mat_SeqAIJ*)(aij->B->data); 346989cf23SStefano Zampini Mat lA; 356989cf23SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 366989cf23SStefano Zampini IS is; 376989cf23SStefano Zampini MPI_Comm comm; 386989cf23SStefano Zampini void *ptrs[2]; 396989cf23SStefano Zampini const char *names[2] = {"_convert_csr_aux","_convert_csr_data"}; 406989cf23SStefano Zampini PetscScalar *dd,*od,*aa,*data; 416989cf23SStefano Zampini PetscInt *di,*dj,*oi,*oj; 426989cf23SStefano Zampini PetscInt *aux,*ii,*jj; 43*e363d98aSStefano Zampini PetscInt lc,dr,dc,oc,str,stc,nnz,i,jd,jo,cum; 446989cf23SStefano Zampini PetscErrorCode ierr; 456989cf23SStefano Zampini 466989cf23SStefano Zampini PetscFunctionBeginUser; 476989cf23SStefano Zampini ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 486989cf23SStefano Zampini if (!aij->garray) SETERRQ(comm,PETSC_ERR_SUP,"garray not present"); 496989cf23SStefano Zampini 506989cf23SStefano Zampini /* access relevant information from MPIAIJ */ 516989cf23SStefano Zampini ierr = MatGetOwnershipRange(A,&str,NULL);CHKERRQ(ierr); 526989cf23SStefano Zampini ierr = MatGetOwnershipRangeColumn(A,&stc,NULL);CHKERRQ(ierr); 536989cf23SStefano Zampini ierr = MatGetLocalSize(A,&dr,&dc);CHKERRQ(ierr); 546989cf23SStefano Zampini di = diag->i; 556989cf23SStefano Zampini dj = diag->j; 566989cf23SStefano Zampini dd = diag->a; 576989cf23SStefano Zampini oc = aij->B->cmap->n; 586989cf23SStefano Zampini oi = offd->i; 596989cf23SStefano Zampini oj = offd->j; 606989cf23SStefano Zampini od = offd->a; 616989cf23SStefano Zampini nnz = diag->i[dr] + offd->i[dr]; 626989cf23SStefano Zampini 636989cf23SStefano Zampini /* generate l2g maps for rows and cols */ 646989cf23SStefano Zampini ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr); 656989cf23SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 666989cf23SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 67*e363d98aSStefano Zampini if (dr) { 686989cf23SStefano Zampini ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr); 696989cf23SStefano Zampini for (i=0; i<dc; i++) aux[i] = i+stc; 706989cf23SStefano Zampini for (i=0; i<oc; i++) aux[i+dc] = aij->garray[i]; 716989cf23SStefano Zampini ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 72*e363d98aSStefano Zampini lc = dc+oc; 73*e363d98aSStefano Zampini } else { 74*e363d98aSStefano Zampini ierr = ISCreateGeneral(comm,0,NULL,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 75*e363d98aSStefano Zampini lc = 0; 76*e363d98aSStefano Zampini } 776989cf23SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 786989cf23SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 796989cf23SStefano Zampini 806989cf23SStefano Zampini /* create MATIS object */ 816989cf23SStefano Zampini ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 826989cf23SStefano Zampini ierr = MatSetSizes(*newmat,dr,dc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 836989cf23SStefano Zampini ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr); 846989cf23SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr); 856989cf23SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 866989cf23SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 876989cf23SStefano Zampini 886989cf23SStefano Zampini /* merge local matrices */ 896989cf23SStefano Zampini ierr = PetscMalloc1(nnz+dr+1,&aux);CHKERRQ(ierr); 906989cf23SStefano Zampini ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr); 916989cf23SStefano Zampini ii = aux; 926989cf23SStefano Zampini jj = aux+dr+1; 936989cf23SStefano Zampini aa = data; 946989cf23SStefano Zampini *ii = *(di++) + *(oi++); 956989cf23SStefano Zampini for (jd=0,jo=0,cum=0;*ii<nnz;cum++) 966989cf23SStefano Zampini { 976989cf23SStefano Zampini for (;jd<*di;jd++) { *jj++ = *dj++; *aa++ = *dd++; } 986989cf23SStefano Zampini for (;jo<*oi;jo++) { *jj++ = *oj++ + dc; *aa++ = *od++; } 996989cf23SStefano Zampini *(++ii) = *(di++) + *(oi++); 1006989cf23SStefano Zampini } 1016989cf23SStefano Zampini for (;cum<dr;cum++) *(++ii) = nnz; 1026989cf23SStefano Zampini ii = aux; 1036989cf23SStefano Zampini jj = aux+dr+1; 1046989cf23SStefano Zampini aa = data; 105*e363d98aSStefano Zampini ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,lc,ii,jj,aa,&lA);CHKERRQ(ierr); 1066989cf23SStefano Zampini 1076989cf23SStefano Zampini /* create containers to destroy the data */ 1086989cf23SStefano Zampini ptrs[0] = aux; 1096989cf23SStefano Zampini ptrs[1] = data; 1106989cf23SStefano Zampini for (i=0; i<2; i++) { 1116989cf23SStefano Zampini PetscContainer c; 1126989cf23SStefano Zampini 1136989cf23SStefano Zampini ierr = PetscContainerCreate(PETSC_COMM_SELF,&c);CHKERRQ(ierr); 1146989cf23SStefano Zampini ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr); 1156989cf23SStefano Zampini ierr = PetscContainerSetUserDestroy(c,MatISContainerDestroy_Private);CHKERRQ(ierr); 1166989cf23SStefano Zampini ierr = PetscObjectCompose((PetscObject)lA,names[i],(PetscObject)c);CHKERRQ(ierr); 1176989cf23SStefano Zampini ierr = PetscContainerDestroy(&c);CHKERRQ(ierr); 1186989cf23SStefano Zampini } 1196989cf23SStefano Zampini 1206989cf23SStefano Zampini /* finalize matrix */ 1216989cf23SStefano Zampini ierr = MatISSetLocalMat(*newmat,lA);CHKERRQ(ierr); 1226989cf23SStefano Zampini ierr = MatDestroy(&lA);CHKERRQ(ierr); 1236989cf23SStefano Zampini ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1246989cf23SStefano Zampini ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1256989cf23SStefano Zampini PetscFunctionReturn(0); 1266989cf23SStefano Zampini } 1276989cf23SStefano Zampini 1286989cf23SStefano Zampini #undef __FUNCT__ 129cf0a3239SStefano Zampini #define __FUNCT__ "MatISSetUpSF" 130cf0a3239SStefano Zampini /*@ 1313d996552SStefano Zampini MatISSetUpSF - Setup star forest objects used by MatIS. 132cf0a3239SStefano Zampini 133cf0a3239SStefano Zampini Collective on MPI_Comm 134cf0a3239SStefano Zampini 135cf0a3239SStefano Zampini Input Parameters: 136cf0a3239SStefano Zampini + A - the matrix 137cf0a3239SStefano Zampini 138cf0a3239SStefano Zampini Level: advanced 139cf0a3239SStefano Zampini 1403d996552SStefano Zampini Notes: This function does not need to be called by the user. 141cf0a3239SStefano Zampini 142cf0a3239SStefano Zampini .keywords: matrix 143cf0a3239SStefano Zampini 144cf0a3239SStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatISSetPreallocation(), MatISGetLocalMat() 145cf0a3239SStefano Zampini @*/ 146cf0a3239SStefano Zampini PetscErrorCode MatISSetUpSF(Mat A) 147cf0a3239SStefano Zampini { 148cf0a3239SStefano Zampini PetscErrorCode ierr; 149cf0a3239SStefano Zampini 150cf0a3239SStefano Zampini PetscFunctionBegin; 151cf0a3239SStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 152cf0a3239SStefano Zampini PetscValidType(A,1); 153cf0a3239SStefano Zampini ierr = PetscTryMethod(A,"MatISSetUpSF_C",(Mat),(A));CHKERRQ(ierr); 154cf0a3239SStefano Zampini PetscFunctionReturn(0); 155cf0a3239SStefano Zampini } 156a72627d2SStefano Zampini 15728f4e0baSStefano Zampini #undef __FUNCT__ 1583fd1c9e7SStefano Zampini #define __FUNCT__ "MatDiagonalSet_IS" 1593fd1c9e7SStefano Zampini PetscErrorCode MatDiagonalSet_IS(Mat A,Vec D,InsertMode insmode) 1603fd1c9e7SStefano Zampini { 1613fd1c9e7SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 1623fd1c9e7SStefano Zampini PetscErrorCode ierr; 1633fd1c9e7SStefano Zampini 1643fd1c9e7SStefano Zampini PetscFunctionBegin; 1653fd1c9e7SStefano Zampini if (!D) { /* this code branch is used by MatShift_IS */ 1663fd1c9e7SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 1673fd1c9e7SStefano Zampini } else { 1683fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1693fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1703fd1c9e7SStefano Zampini } 1713fd1c9e7SStefano Zampini ierr = VecPointwiseDivide(is->y,is->y,is->counter);CHKERRQ(ierr); 1723fd1c9e7SStefano Zampini ierr = MatDiagonalSet(is->A,is->y,insmode);CHKERRQ(ierr); 1733fd1c9e7SStefano Zampini PetscFunctionReturn(0); 1743fd1c9e7SStefano Zampini } 1753fd1c9e7SStefano Zampini 1763fd1c9e7SStefano Zampini #undef __FUNCT__ 1773fd1c9e7SStefano Zampini #define __FUNCT__ "MatShift_IS" 1783fd1c9e7SStefano Zampini PetscErrorCode MatShift_IS(Mat A,PetscScalar a) 1793fd1c9e7SStefano Zampini { 1803fd1c9e7SStefano Zampini PetscErrorCode ierr; 1813fd1c9e7SStefano Zampini 1823fd1c9e7SStefano Zampini PetscFunctionBegin; 1833fd1c9e7SStefano Zampini ierr = MatDiagonalSet_IS(A,NULL,ADD_VALUES);CHKERRQ(ierr); 1843fd1c9e7SStefano Zampini PetscFunctionReturn(0); 1853fd1c9e7SStefano Zampini } 1863fd1c9e7SStefano Zampini 1873fd1c9e7SStefano Zampini #undef __FUNCT__ 188f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesLocal_SubMat_IS" 189f26d0771SStefano Zampini static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 190f26d0771SStefano Zampini { 191f26d0771SStefano Zampini PetscErrorCode ierr; 192f26d0771SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 193f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 194f26d0771SStefano Zampini 195f26d0771SStefano Zampini PetscFunctionBegin; 196f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 197f26d0771SStefano 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); 198f26d0771SStefano Zampini #endif 199f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 200f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 201f26d0771SStefano Zampini ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 202f26d0771SStefano Zampini PetscFunctionReturn(0); 203f26d0771SStefano Zampini } 204f26d0771SStefano Zampini 205f26d0771SStefano Zampini #undef __FUNCT__ 206f26d0771SStefano Zampini #define __FUNCT__ "MatSetValuesBlockedLocal_SubMat_IS" 207f26d0771SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 208f26d0771SStefano Zampini { 209f26d0771SStefano Zampini PetscErrorCode ierr; 210f26d0771SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 211f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 212f26d0771SStefano Zampini 213f26d0771SStefano Zampini PetscFunctionBegin; 214f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 215f26d0771SStefano 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); 216f26d0771SStefano Zampini #endif 217f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 218f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 219f26d0771SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 220f26d0771SStefano Zampini PetscFunctionReturn(0); 221f26d0771SStefano Zampini } 222f26d0771SStefano Zampini 223f26d0771SStefano Zampini #undef __FUNCT__ 224a8116848SStefano Zampini #define __FUNCT__ "PetscLayoutMapLocal_Private" 225f0ae7da4SStefano Zampini static PetscErrorCode PetscLayoutMapLocal_Private(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs) 226a8116848SStefano Zampini { 227a8116848SStefano Zampini PetscInt *owners = map->range; 228a8116848SStefano Zampini PetscInt n = map->n; 229a8116848SStefano Zampini PetscSF sf; 230a8116848SStefano Zampini PetscInt *lidxs,*work = NULL; 231a8116848SStefano Zampini PetscSFNode *ridxs; 232a8116848SStefano Zampini PetscMPIInt rank; 233a8116848SStefano Zampini PetscInt r, p = 0, len = 0; 234a8116848SStefano Zampini PetscErrorCode ierr; 235a8116848SStefano Zampini 236a8116848SStefano Zampini PetscFunctionBegin; 237a8116848SStefano Zampini /* Create SF where leaves are input idxs and roots are owned idxs (code adapted from MatZeroRowsMapLocal_Private) */ 238a8116848SStefano Zampini ierr = MPI_Comm_rank(map->comm,&rank);CHKERRQ(ierr); 239a8116848SStefano Zampini ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr); 240a8116848SStefano Zampini for (r = 0; r < n; ++r) lidxs[r] = -1; 241a8116848SStefano Zampini ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr); 242a8116848SStefano Zampini for (r = 0; r < N; ++r) { 243a8116848SStefano Zampini const PetscInt idx = idxs[r]; 244a8116848SStefano 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); 245a8116848SStefano Zampini if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 246a8116848SStefano Zampini ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr); 247a8116848SStefano Zampini } 248a8116848SStefano Zampini ridxs[r].rank = p; 249a8116848SStefano Zampini ridxs[r].index = idxs[r] - owners[p]; 250a8116848SStefano Zampini } 251a8116848SStefano Zampini ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr); 252a8116848SStefano Zampini ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr); 253a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 254a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 255f0ae7da4SStefano Zampini if (ogidxs) { /* communicate global idxs */ 256a8116848SStefano Zampini PetscInt cum = 0,start,*work2; 257f0ae7da4SStefano Zampini 258f0ae7da4SStefano Zampini ierr = PetscMalloc1(n,&work);CHKERRQ(ierr); 259a8116848SStefano Zampini ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr); 260a8116848SStefano Zampini for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++; 261a8116848SStefano Zampini ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRQ(ierr); 262a8116848SStefano Zampini start -= cum; 263a8116848SStefano Zampini cum = 0; 264a8116848SStefano Zampini for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++; 265a8116848SStefano Zampini ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 266a8116848SStefano Zampini ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 267a8116848SStefano Zampini ierr = PetscFree(work2);CHKERRQ(ierr); 268a8116848SStefano Zampini } 269a8116848SStefano Zampini ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 270a8116848SStefano Zampini /* Compress and put in indices */ 271a8116848SStefano Zampini for (r = 0; r < n; ++r) 272a8116848SStefano Zampini if (lidxs[r] >= 0) { 273a8116848SStefano Zampini if (work) work[len] = work[r]; 274a8116848SStefano Zampini lidxs[len++] = r; 275a8116848SStefano Zampini } 276a8116848SStefano Zampini if (on) *on = len; 277a8116848SStefano Zampini if (oidxs) *oidxs = lidxs; 278a8116848SStefano Zampini if (ogidxs) *ogidxs = work; 279a8116848SStefano Zampini PetscFunctionReturn(0); 280a8116848SStefano Zampini } 281a8116848SStefano Zampini 282a8116848SStefano Zampini #undef __FUNCT__ 283a8116848SStefano Zampini #define __FUNCT__ "MatGetSubMatrix_IS" 284a8116848SStefano Zampini static PetscErrorCode MatGetSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat) 285a8116848SStefano Zampini { 286a8116848SStefano Zampini Mat locmat,newlocmat; 287a8116848SStefano Zampini Mat_IS *newmatis; 288a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG) 289a8116848SStefano Zampini Vec rtest,ltest; 290a8116848SStefano Zampini const PetscScalar *array; 291a8116848SStefano Zampini #endif 292a8116848SStefano Zampini const PetscInt *idxs; 293a8116848SStefano Zampini PetscInt i,m,n; 294a8116848SStefano Zampini PetscErrorCode ierr; 295a8116848SStefano Zampini 296a8116848SStefano Zampini PetscFunctionBegin; 297a8116848SStefano Zampini if (scall == MAT_REUSE_MATRIX) { 298a8116848SStefano Zampini PetscBool ismatis; 299a8116848SStefano Zampini 300a8116848SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr); 301a8116848SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type"); 302a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 303a8116848SStefano Zampini if (!newmatis->getsub_ris) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS"); 304a8116848SStefano Zampini if (!newmatis->getsub_cis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS"); 305a8116848SStefano Zampini } 306a8116848SStefano Zampini /* irow and icol may not have duplicate entries */ 307a8116848SStefano Zampini #if defined(PETSC_USE_DEBUG) 308a8116848SStefano Zampini ierr = MatCreateVecs(mat,<est,&rtest);CHKERRQ(ierr); 309a8116848SStefano Zampini ierr = ISGetLocalSize(irow,&n);CHKERRQ(ierr); 310a8116848SStefano Zampini ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 311a8116848SStefano Zampini for (i=0;i<n;i++) { 312a8116848SStefano Zampini ierr = VecSetValue(rtest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 313a8116848SStefano Zampini } 314a8116848SStefano Zampini ierr = VecAssemblyBegin(rtest);CHKERRQ(ierr); 315a8116848SStefano Zampini ierr = VecAssemblyEnd(rtest);CHKERRQ(ierr); 316a8116848SStefano Zampini ierr = VecGetLocalSize(rtest,&n);CHKERRQ(ierr); 317a8116848SStefano Zampini ierr = VecGetOwnershipRange(rtest,&m,NULL);CHKERRQ(ierr); 318a8116848SStefano Zampini ierr = VecGetArrayRead(rtest,&array);CHKERRQ(ierr); 319fd479f66SMatthew 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])); 320a8116848SStefano Zampini ierr = VecRestoreArrayRead(rtest,&array);CHKERRQ(ierr); 321a8116848SStefano Zampini ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 322a8116848SStefano Zampini ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 323a8116848SStefano Zampini ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 324a8116848SStefano Zampini for (i=0;i<n;i++) { 325a8116848SStefano Zampini ierr = VecSetValue(ltest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 326a8116848SStefano Zampini } 327a8116848SStefano Zampini ierr = VecAssemblyBegin(ltest);CHKERRQ(ierr); 328a8116848SStefano Zampini ierr = VecAssemblyEnd(ltest);CHKERRQ(ierr); 329a8116848SStefano Zampini ierr = VecGetLocalSize(ltest,&n);CHKERRQ(ierr); 330a8116848SStefano Zampini ierr = VecGetOwnershipRange(ltest,&m,NULL);CHKERRQ(ierr); 331a8116848SStefano Zampini ierr = VecGetArrayRead(ltest,&array);CHKERRQ(ierr); 332fd479f66SMatthew 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])); 333a8116848SStefano Zampini ierr = VecRestoreArrayRead(ltest,&array);CHKERRQ(ierr); 334a8116848SStefano Zampini ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 335a8116848SStefano Zampini ierr = VecDestroy(&rtest);CHKERRQ(ierr); 336a8116848SStefano Zampini ierr = VecDestroy(<est);CHKERRQ(ierr); 337a8116848SStefano Zampini #endif 338a8116848SStefano Zampini if (scall == MAT_INITIAL_MATRIX) { 339a8116848SStefano Zampini Mat_IS *matis = (Mat_IS*)mat->data; 340a8116848SStefano Zampini ISLocalToGlobalMapping rl2g; 341a8116848SStefano Zampini IS is; 342a8116848SStefano Zampini PetscInt *lidxs,*lgidxs,*newgidxs; 3436dd40735SStefano Zampini PetscInt ll,newloc; 344a8116848SStefano Zampini MPI_Comm comm; 345a8116848SStefano Zampini 346a8116848SStefano Zampini ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 347a8116848SStefano Zampini ierr = ISGetLocalSize(irow,&m);CHKERRQ(ierr); 348a8116848SStefano Zampini ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 349a8116848SStefano Zampini ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 350a8116848SStefano Zampini ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr); 351a8116848SStefano Zampini ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 352a8116848SStefano Zampini /* communicate irow to their owners in the layout */ 353a8116848SStefano Zampini ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 354f0ae7da4SStefano Zampini ierr = PetscLayoutMapLocal_Private(mat->rmap,m,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 355a8116848SStefano Zampini ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 3563d996552SStefano Zampini ierr = MatISSetUpSF(mat);CHKERRQ(ierr); 3573d996552SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,matis->sf->nroots*sizeof(PetscInt));CHKERRQ(ierr); 358a8116848SStefano Zampini for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1; 359a8116848SStefano Zampini ierr = PetscFree(lidxs);CHKERRQ(ierr); 360a8116848SStefano Zampini ierr = PetscFree(lgidxs);CHKERRQ(ierr); 361a8116848SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 362a8116848SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 3633d996552SStefano Zampini for (i=0,newloc=0;i<matis->sf->nleaves;i++) if (matis->sf_leafdata[i]) newloc++; 364a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 365a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 3663d996552SStefano Zampini for (i=0,newloc=0;i<matis->sf->nleaves;i++) 367a8116848SStefano Zampini if (matis->sf_leafdata[i]) { 368a8116848SStefano Zampini lidxs[newloc] = i; 369a8116848SStefano Zampini newgidxs[newloc++] = matis->sf_leafdata[i]-1; 370a8116848SStefano Zampini } 371a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 372a8116848SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 373a8116848SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 374a8116848SStefano Zampini /* local is to extract local submatrix */ 375a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 376a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris);CHKERRQ(ierr); 377a8116848SStefano Zampini if (mat->congruentlayouts == -1) { /* first time we compare rows and cols layouts */ 378a8116848SStefano Zampini PetscBool cong; 379a8116848SStefano Zampini ierr = PetscLayoutCompare(mat->rmap,mat->cmap,&cong);CHKERRQ(ierr); 380a8116848SStefano Zampini if (cong) mat->congruentlayouts = 1; 381a8116848SStefano Zampini else mat->congruentlayouts = 0; 382a8116848SStefano Zampini } 383a8116848SStefano Zampini if (mat->congruentlayouts && irow == icol) { 384a8116848SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g);CHKERRQ(ierr); 385a8116848SStefano Zampini ierr = PetscObjectReference((PetscObject)newmatis->getsub_ris);CHKERRQ(ierr); 386a8116848SStefano Zampini newmatis->getsub_cis = newmatis->getsub_ris; 387a8116848SStefano Zampini } else { 388a8116848SStefano Zampini ISLocalToGlobalMapping cl2g; 389a8116848SStefano Zampini 390a8116848SStefano Zampini /* communicate icol to their owners in the layout */ 391a8116848SStefano Zampini ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 392f0ae7da4SStefano Zampini ierr = PetscLayoutMapLocal_Private(mat->cmap,n,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 393a8116848SStefano Zampini ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 3943d996552SStefano Zampini ierr = PetscMemzero(matis->csf_rootdata,matis->csf->nroots*sizeof(PetscInt));CHKERRQ(ierr); 395a8116848SStefano Zampini for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1; 396a8116848SStefano Zampini ierr = PetscFree(lidxs);CHKERRQ(ierr); 397a8116848SStefano Zampini ierr = PetscFree(lgidxs);CHKERRQ(ierr); 398a8116848SStefano Zampini ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 399a8116848SStefano Zampini ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 4003d996552SStefano Zampini for (i=0,newloc=0;i<matis->csf->nleaves;i++) if (matis->csf_leafdata[i]) newloc++; 401a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 402a8116848SStefano Zampini ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 4033d996552SStefano Zampini for (i=0,newloc=0;i<matis->csf->nleaves;i++) 404a8116848SStefano Zampini if (matis->csf_leafdata[i]) { 405a8116848SStefano Zampini lidxs[newloc] = i; 406a8116848SStefano Zampini newgidxs[newloc++] = matis->csf_leafdata[i]-1; 407a8116848SStefano Zampini } 408a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 409a8116848SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 410a8116848SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 411a8116848SStefano Zampini /* local is to extract local submatrix */ 412a8116848SStefano Zampini ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis);CHKERRQ(ierr); 413a8116848SStefano Zampini ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr); 414a8116848SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 415a8116848SStefano Zampini } 416a8116848SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 417a8116848SStefano Zampini } else { 418a8116848SStefano Zampini ierr = MatISGetLocalMat(*newmat,&newlocmat);CHKERRQ(ierr); 419a8116848SStefano Zampini } 420a8116848SStefano Zampini ierr = MatISGetLocalMat(mat,&locmat);CHKERRQ(ierr); 421a8116848SStefano Zampini newmatis = (Mat_IS*)(*newmat)->data; 422a8116848SStefano Zampini ierr = MatGetSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat);CHKERRQ(ierr); 423a8116848SStefano Zampini if (scall == MAT_INITIAL_MATRIX) { 424a8116848SStefano Zampini ierr = MatISSetLocalMat(*newmat,newlocmat);CHKERRQ(ierr); 425a8116848SStefano Zampini ierr = MatDestroy(&newlocmat);CHKERRQ(ierr); 426a8116848SStefano Zampini } 427a8116848SStefano Zampini ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 428a8116848SStefano Zampini ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 429a8116848SStefano Zampini PetscFunctionReturn(0); 430a8116848SStefano Zampini } 431a8116848SStefano Zampini 432a8116848SStefano Zampini #undef __FUNCT__ 4332b404112SStefano Zampini #define __FUNCT__ "MatCopy_IS" 434a8116848SStefano Zampini static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str) 4352b404112SStefano Zampini { 4362b404112SStefano Zampini Mat_IS *a = (Mat_IS*)A->data,*b; 4372b404112SStefano Zampini PetscBool ismatis; 4382b404112SStefano Zampini PetscErrorCode ierr; 4392b404112SStefano Zampini 4402b404112SStefano Zampini PetscFunctionBegin; 4412b404112SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr); 4422b404112SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented"); 4432b404112SStefano Zampini b = (Mat_IS*)B->data; 4442b404112SStefano Zampini ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr); 4452b404112SStefano Zampini PetscFunctionReturn(0); 4462b404112SStefano Zampini } 4472b404112SStefano Zampini 4482b404112SStefano Zampini #undef __FUNCT__ 4496bd84002SStefano Zampini #define __FUNCT__ "MatMissingDiagonal_IS" 450a8116848SStefano Zampini static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool *missing,PetscInt *d) 4516bd84002SStefano Zampini { 452527b2640SStefano Zampini Vec v; 453527b2640SStefano Zampini const PetscScalar *array; 454527b2640SStefano Zampini PetscInt i,n; 4556bd84002SStefano Zampini PetscErrorCode ierr; 4566bd84002SStefano Zampini 4576bd84002SStefano Zampini PetscFunctionBegin; 458527b2640SStefano Zampini *missing = PETSC_FALSE; 459527b2640SStefano Zampini ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr); 460527b2640SStefano Zampini ierr = MatGetDiagonal(A,v);CHKERRQ(ierr); 461527b2640SStefano Zampini ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 462527b2640SStefano Zampini ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr); 463527b2640SStefano Zampini for (i=0;i<n;i++) if (array[i] == 0.) break; 464527b2640SStefano Zampini ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr); 465527b2640SStefano Zampini ierr = VecDestroy(&v);CHKERRQ(ierr); 466527b2640SStefano Zampini if (i != n) *missing = PETSC_TRUE; 467527b2640SStefano Zampini if (d) { 468527b2640SStefano Zampini *d = -1; 469527b2640SStefano Zampini if (*missing) { 470527b2640SStefano Zampini PetscInt rstart; 471527b2640SStefano Zampini ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 472527b2640SStefano Zampini *d = i+rstart; 473527b2640SStefano Zampini } 474527b2640SStefano Zampini } 4756bd84002SStefano Zampini PetscFunctionReturn(0); 4766bd84002SStefano Zampini } 4776bd84002SStefano Zampini 4786bd84002SStefano Zampini #undef __FUNCT__ 479cf0a3239SStefano Zampini #define __FUNCT__ "MatISSetUpSF_IS" 480cf0a3239SStefano Zampini static PetscErrorCode MatISSetUpSF_IS(Mat B) 48128f4e0baSStefano Zampini { 48228f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 48328f4e0baSStefano Zampini const PetscInt *gidxs; 4844f2d7cafSStefano Zampini PetscInt nleaves; 48528f4e0baSStefano Zampini PetscErrorCode ierr; 48628f4e0baSStefano Zampini 48728f4e0baSStefano Zampini PetscFunctionBegin; 4884f2d7cafSStefano Zampini if (matis->sf) PetscFunctionReturn(0); 48928f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 4903bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 4914f2d7cafSStefano Zampini ierr = ISLocalToGlobalMappingGetSize(B->rmap->mapping,&nleaves);CHKERRQ(ierr); 4924f2d7cafSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 4933bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 4944f2d7cafSStefano Zampini ierr = PetscMalloc2(matis->sf->nroots,&matis->sf_rootdata,matis->sf->nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 495a8116848SStefano Zampini if (B->rmap->mapping != B->cmap->mapping) { /* setup SF for columns */ 4963d996552SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(B->cmap->mapping,&nleaves);CHKERRQ(ierr); 497a8116848SStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf);CHKERRQ(ierr); 498a8116848SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 4993d996552SStefano Zampini ierr = PetscSFSetGraphLayout(matis->csf,B->cmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 500a8116848SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 5013d996552SStefano Zampini ierr = PetscMalloc2(matis->csf->nroots,&matis->csf_rootdata,matis->csf->nleaves,&matis->csf_leafdata);CHKERRQ(ierr); 502a8116848SStefano Zampini } else { 503a8116848SStefano Zampini matis->csf = matis->sf; 504a8116848SStefano Zampini matis->csf_leafdata = matis->sf_leafdata; 505a8116848SStefano Zampini matis->csf_rootdata = matis->sf_rootdata; 506a8116848SStefano Zampini } 50728f4e0baSStefano Zampini PetscFunctionReturn(0); 50828f4e0baSStefano Zampini } 5092e1947a5SStefano Zampini 5102e1947a5SStefano Zampini #undef __FUNCT__ 5112e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 512eb82efa4SStefano Zampini /*@ 513a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 514a88811baSStefano Zampini 515a88811baSStefano Zampini Collective on MPI_Comm 516a88811baSStefano Zampini 517a88811baSStefano Zampini Input Parameters: 518a88811baSStefano Zampini + B - the matrix 519a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 520a88811baSStefano Zampini (same value is used for all local rows) 521a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 522a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 523a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 524a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 525a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 526a88811baSStefano Zampini the diagonal entry even if it is zero. 527a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 528a88811baSStefano Zampini submatrix (same value is used for all local rows). 529a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 530a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 531a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 532a88811baSStefano Zampini structure. The size of this array is equal to the number 533a88811baSStefano Zampini of local rows, i.e 'm'. 534a88811baSStefano Zampini 535a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 536a88811baSStefano Zampini 537a88811baSStefano Zampini Level: intermediate 538a88811baSStefano Zampini 539a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 540a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 541a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 542a88811baSStefano Zampini 543a88811baSStefano Zampini .keywords: matrix 544a88811baSStefano Zampini 5453c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS 546a88811baSStefano Zampini @*/ 5472e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 5482e1947a5SStefano Zampini { 5492e1947a5SStefano Zampini PetscErrorCode ierr; 5502e1947a5SStefano Zampini 5512e1947a5SStefano Zampini PetscFunctionBegin; 5522e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 5532e1947a5SStefano Zampini PetscValidType(B,1); 5542e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 5552e1947a5SStefano Zampini PetscFunctionReturn(0); 5562e1947a5SStefano Zampini } 5572e1947a5SStefano Zampini 5582e1947a5SStefano Zampini #undef __FUNCT__ 5592e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 5607230de76SStefano Zampini static PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 5612e1947a5SStefano Zampini { 5622e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 56328f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 5642e1947a5SStefano Zampini PetscErrorCode ierr; 5652e1947a5SStefano Zampini 5662e1947a5SStefano Zampini PetscFunctionBegin; 5676c4ed002SBarry Smith if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 568cf0a3239SStefano Zampini ierr = MatISSetUpSF(B);CHKERRQ(ierr); 5694f2d7cafSStefano Zampini 5704f2d7cafSStefano Zampini if (!d_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nz; 5714f2d7cafSStefano Zampini else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 5724f2d7cafSStefano Zampini 5734f2d7cafSStefano Zampini if (!o_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nz; 5744f2d7cafSStefano Zampini else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 5754f2d7cafSStefano Zampini 57628f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 57728f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 57828f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 57928f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 5804f2d7cafSStefano Zampini 5814f2d7cafSStefano Zampini for (i=0;i<matis->sf->nleaves;i++) matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 58228f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 5834f2d7cafSStefano Zampini 5844f2d7cafSStefano Zampini for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 58528f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 5864f2d7cafSStefano Zampini 5874f2d7cafSStefano Zampini for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 58828f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 5892e1947a5SStefano Zampini PetscFunctionReturn(0); 5902e1947a5SStefano Zampini } 591b4319ba4SBarry Smith 592b4319ba4SBarry Smith #undef __FUNCT__ 5933927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private" 5943927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 5953927de2eSStefano Zampini { 5963927de2eSStefano Zampini Mat_IS *matis = (Mat_IS*)(A->data); 5973927de2eSStefano Zampini PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 598ecf5a873SStefano Zampini const PetscInt *global_indices_r,*global_indices_c; 5993927de2eSStefano Zampini PetscInt i,j,bs,rows,cols; 6003927de2eSStefano Zampini PetscInt lrows,lcols; 6013927de2eSStefano Zampini PetscInt local_rows,local_cols; 6023927de2eSStefano Zampini PetscMPIInt nsubdomains; 6033927de2eSStefano Zampini PetscBool isdense,issbaij; 6043927de2eSStefano Zampini PetscErrorCode ierr; 6053927de2eSStefano Zampini 6063927de2eSStefano Zampini PetscFunctionBegin; 6073927de2eSStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr); 6083927de2eSStefano Zampini ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr); 6093927de2eSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 6103927de2eSStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 6113927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 6123927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 613ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 614ecf5a873SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 6157230de76SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 616ecf5a873SStefano Zampini } else { 617ecf5a873SStefano Zampini global_indices_c = global_indices_r; 618ecf5a873SStefano Zampini } 619ecf5a873SStefano Zampini 6203927de2eSStefano Zampini if (issbaij) { 6213927de2eSStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 6223927de2eSStefano Zampini } 6233927de2eSStefano Zampini /* 624ecf5a873SStefano Zampini An SF reduce is needed to sum up properly on shared rows. 6253927de2eSStefano Zampini Note that generally preallocation is not exact, since it overestimates nonzeros 6263927de2eSStefano Zampini */ 627cf0a3239SStefano Zampini ierr = MatISSetUpSF(A);CHKERRQ(ierr); 6283927de2eSStefano Zampini ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr); 6293927de2eSStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr); 6303927de2eSStefano Zampini /* All processes need to compute entire row ownership */ 6313927de2eSStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 6323927de2eSStefano Zampini ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 6333927de2eSStefano Zampini for (i=0;i<nsubdomains;i++) { 6343927de2eSStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 6353927de2eSStefano Zampini row_ownership[j] = i; 6363927de2eSStefano Zampini } 6373927de2eSStefano Zampini } 6387230de76SStefano Zampini ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 6393927de2eSStefano Zampini 6403927de2eSStefano Zampini /* 6413927de2eSStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 6423927de2eSStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 6433927de2eSStefano Zampini */ 6443927de2eSStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 6453927de2eSStefano Zampini /* preallocation as a MATAIJ */ 6463927de2eSStefano Zampini if (isdense) { /* special case for dense local matrices */ 6473927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 648ecf5a873SStefano Zampini PetscInt index_row = global_indices_r[i]; 6493927de2eSStefano Zampini for (j=i;j<local_rows;j++) { 6503927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 651ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[j]; 6523927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 6533927de2eSStefano Zampini my_dnz[i] += 1; 6543927de2eSStefano Zampini } else { /* offdiag block */ 6553927de2eSStefano Zampini my_onz[i] += 1; 6563927de2eSStefano Zampini } 6573927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 6583927de2eSStefano Zampini if (i != j) { 6593927de2eSStefano Zampini owner = row_ownership[index_col]; 6603927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 6613927de2eSStefano Zampini my_dnz[j] += 1; 6623927de2eSStefano Zampini } else { 6633927de2eSStefano Zampini my_onz[j] += 1; 6643927de2eSStefano Zampini } 6653927de2eSStefano Zampini } 6663927de2eSStefano Zampini } 6673927de2eSStefano Zampini } 668bb1015c3SStefano Zampini } else if (matis->A->ops->getrowij) { 669bb1015c3SStefano Zampini const PetscInt *ii,*jj,*jptr; 670bb1015c3SStefano Zampini PetscBool done; 671bb1015c3SStefano Zampini ierr = MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr); 672bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__); 673bb1015c3SStefano Zampini jptr = jj; 674bb1015c3SStefano Zampini for (i=0;i<local_rows;i++) { 675bb1015c3SStefano Zampini PetscInt index_row = global_indices_r[i]; 676bb1015c3SStefano Zampini for (j=0;j<ii[i+1]-ii[i];j++,jptr++) { 677bb1015c3SStefano Zampini PetscInt owner = row_ownership[index_row]; 678bb1015c3SStefano Zampini PetscInt index_col = global_indices_c[*jptr]; 679bb1015c3SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 680bb1015c3SStefano Zampini my_dnz[i] += 1; 681bb1015c3SStefano Zampini } else { /* offdiag block */ 682bb1015c3SStefano Zampini my_onz[i] += 1; 683bb1015c3SStefano Zampini } 684bb1015c3SStefano Zampini /* same as before, interchanging rows and cols */ 685bb1015c3SStefano Zampini if (issbaij && index_col != index_row) { 686bb1015c3SStefano Zampini owner = row_ownership[index_col]; 687bb1015c3SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 688bb1015c3SStefano Zampini my_dnz[*jptr] += 1; 689bb1015c3SStefano Zampini } else { 690bb1015c3SStefano Zampini my_onz[*jptr] += 1; 691bb1015c3SStefano Zampini } 692bb1015c3SStefano Zampini } 693bb1015c3SStefano Zampini } 694bb1015c3SStefano Zampini } 695bb1015c3SStefano Zampini ierr = MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr); 696bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__); 697bb1015c3SStefano Zampini } else { /* loop over rows and use MatGetRow */ 6983927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 6993927de2eSStefano Zampini const PetscInt *cols; 700ecf5a873SStefano Zampini PetscInt ncols,index_row = global_indices_r[i]; 7013927de2eSStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 7023927de2eSStefano Zampini for (j=0;j<ncols;j++) { 7033927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 704ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[cols[j]]; 7053927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 7063927de2eSStefano Zampini my_dnz[i] += 1; 7073927de2eSStefano Zampini } else { /* offdiag block */ 7083927de2eSStefano Zampini my_onz[i] += 1; 7093927de2eSStefano Zampini } 7103927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 711d9a9e74cSStefano Zampini if (issbaij && index_col != index_row) { 7123927de2eSStefano Zampini owner = row_ownership[index_col]; 7133927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 714d9a9e74cSStefano Zampini my_dnz[cols[j]] += 1; 7153927de2eSStefano Zampini } else { 716d9a9e74cSStefano Zampini my_onz[cols[j]] += 1; 7173927de2eSStefano Zampini } 7183927de2eSStefano Zampini } 7193927de2eSStefano Zampini } 7203927de2eSStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 7213927de2eSStefano Zampini } 7223927de2eSStefano Zampini } 723ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 724ecf5a873SStefano Zampini if (global_indices_c != global_indices_r) { 7257230de76SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 726ecf5a873SStefano Zampini } 7273927de2eSStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 728ecf5a873SStefano Zampini 729ecf5a873SStefano Zampini /* Reduce my_dnz and my_onz */ 7303927de2eSStefano Zampini if (maxreduce) { 7313927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 7323927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 733bb1015c3SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 7343927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 7353927de2eSStefano Zampini } else { 7363927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 7373927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 738bb1015c3SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 7393927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 7403927de2eSStefano Zampini } 7413927de2eSStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 7423927de2eSStefano Zampini 7433927de2eSStefano Zampini /* Resize preallocation if overestimated */ 7443927de2eSStefano Zampini for (i=0;i<lrows;i++) { 7453927de2eSStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 7463927de2eSStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 7473927de2eSStefano Zampini } 7483927de2eSStefano Zampini /* set preallocation */ 7493927de2eSStefano Zampini ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 7503927de2eSStefano Zampini for (i=0;i<lrows/bs;i++) { 7513927de2eSStefano Zampini dnz[i] = dnz[i*bs]/bs; 7523927de2eSStefano Zampini onz[i] = onz[i*bs]/bs; 7533927de2eSStefano Zampini } 7543927de2eSStefano Zampini ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 7553927de2eSStefano Zampini ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 7563927de2eSStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 7573927de2eSStefano Zampini if (issbaij) { 7583927de2eSStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 7593927de2eSStefano Zampini } 7603927de2eSStefano Zampini PetscFunctionReturn(0); 7613927de2eSStefano Zampini } 7623927de2eSStefano Zampini 7633927de2eSStefano Zampini #undef __FUNCT__ 764b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 7657230de76SStefano Zampini static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 766b7ce53b6SStefano Zampini { 767b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 768d9a9e74cSStefano Zampini Mat local_mat; 769b7ce53b6SStefano Zampini /* info on mat */ 7703cfa4ea4SStefano Zampini PetscInt bs,rows,cols,lrows,lcols; 771b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 772b9ed4604SStefano Zampini PetscBool isseqdense,isseqsbaij,isseqaij,isseqbaij; 773b9ed4604SStefano Zampini #if defined (PETSC_USE_DEBUG) 774b9ed4604SStefano Zampini PetscBool lb[4],bb[4]; 775b9ed4604SStefano Zampini #endif 7767c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 777b7ce53b6SStefano Zampini /* values insertion */ 778b7ce53b6SStefano Zampini PetscScalar *array; 779b7ce53b6SStefano Zampini /* work */ 780b7ce53b6SStefano Zampini PetscErrorCode ierr; 781b7ce53b6SStefano Zampini 782b7ce53b6SStefano Zampini PetscFunctionBegin; 783b7ce53b6SStefano Zampini /* get info from mat */ 7847c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 7857c03b4e8SStefano Zampini if (nsubdomains == 1) { 7867c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 7877c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 7887c03b4e8SStefano Zampini } else { 7897c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 7907c03b4e8SStefano Zampini } 7917c03b4e8SStefano Zampini PetscFunctionReturn(0); 7927c03b4e8SStefano Zampini } 793b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 794b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 7953cfa4ea4SStefano Zampini ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 796b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 797b9ed4604SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isseqdense);CHKERRQ(ierr); 798686e3a49SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 799b9ed4604SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQBAIJ,&isseqbaij);CHKERRQ(ierr); 800b9ed4604SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr); 801b9ed4604SStefano Zampini if (!isseqdense && !isseqaij && !isseqbaij && !isseqsbaij) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)(matis->A))->type_name); 802b9ed4604SStefano Zampini #if defined (PETSC_USE_DEBUG) 803b9ed4604SStefano Zampini lb[0] = isseqdense; 804b9ed4604SStefano Zampini lb[1] = isseqaij; 805b9ed4604SStefano Zampini lb[2] = isseqbaij; 806b9ed4604SStefano Zampini lb[3] = isseqsbaij; 807b9ed4604SStefano Zampini ierr = MPIU_Allreduce(lb,bb,4,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 808b9ed4604SStefano Zampini if (!bb[0] && !bb[1] && !bb[2] && !bb[3]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local matrices must have the same type"); 809b9ed4604SStefano Zampini #endif 810b7ce53b6SStefano Zampini 811b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 8123927de2eSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 8133cfa4ea4SStefano Zampini ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 8143927de2eSStefano Zampini ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 815b9ed4604SStefano Zampini if (!isseqsbaij) { 816b9ed4604SStefano Zampini ierr = MatSetType(*M,MATAIJ);CHKERRQ(ierr); 817b9ed4604SStefano Zampini } else { 818b9ed4604SStefano Zampini ierr = MatSetType(*M,MATSBAIJ);CHKERRQ(ierr); 819b9ed4604SStefano Zampini } 8203927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 821b7ce53b6SStefano Zampini } else { 8223cfa4ea4SStefano Zampini PetscInt mbs,mrows,mcols,mlrows,mlcols; 823b7ce53b6SStefano Zampini /* some checks */ 824b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 825b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 8263cfa4ea4SStefano Zampini ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 8276c4ed002SBarry Smith if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 8286c4ed002SBarry Smith if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 8296c4ed002SBarry Smith if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 8306c4ed002SBarry Smith if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 8316c4ed002SBarry Smith if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 832b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 833b7ce53b6SStefano Zampini } 834d9a9e74cSStefano Zampini 835b9ed4604SStefano Zampini if (isseqsbaij) { 836d9a9e74cSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 837d9a9e74cSStefano Zampini } else { 838d9a9e74cSStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 839d9a9e74cSStefano Zampini local_mat = matis->A; 840d9a9e74cSStefano Zampini } 841686e3a49SStefano Zampini 842b7ce53b6SStefano Zampini /* Set values */ 843ecf5a873SStefano Zampini ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 844b9ed4604SStefano Zampini if (isseqdense) { /* special case for dense local matrices */ 845ecf5a873SStefano Zampini PetscInt i,*dummy_rows,*dummy_cols; 846ecf5a873SStefano Zampini 847ecf5a873SStefano Zampini if (local_rows != local_cols) { 848ecf5a873SStefano Zampini ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 849ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 850ecf5a873SStefano Zampini for (i=0;i<local_cols;i++) dummy_cols[i] = i; 851ecf5a873SStefano Zampini } else { 852ecf5a873SStefano Zampini ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 853ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 854ecf5a873SStefano Zampini dummy_cols = dummy_rows; 855ecf5a873SStefano Zampini } 856b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 857d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 858ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 859d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 860ecf5a873SStefano Zampini if (dummy_rows != dummy_cols) { 861ecf5a873SStefano Zampini ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 862ecf5a873SStefano Zampini } else { 863ecf5a873SStefano Zampini ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 864ecf5a873SStefano Zampini } 865686e3a49SStefano Zampini } else if (isseqaij) { 866ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 867686e3a49SStefano Zampini PetscBool done; 868686e3a49SStefano Zampini 869d9a9e74cSStefano Zampini ierr = MatGetRowIJ(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 MatGetRowIJ called from %s\n",__FUNCT__); 871d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 872686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 873ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 874686e3a49SStefano Zampini } 875d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 876bb1015c3SStefano Zampini if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__); 877d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 878686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 879ecf5a873SStefano Zampini PetscInt i; 880c0962df8SStefano Zampini 881686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 882686e3a49SStefano Zampini PetscInt j; 883ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 884686e3a49SStefano Zampini 885ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 886ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 887ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 888686e3a49SStefano Zampini } 889b7ce53b6SStefano Zampini } 890b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 891d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 892b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 893b9ed4604SStefano Zampini if (isseqdense) { 894b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 895b7ce53b6SStefano Zampini } 896b7ce53b6SStefano Zampini PetscFunctionReturn(0); 897b7ce53b6SStefano Zampini } 898b7ce53b6SStefano Zampini 899b7ce53b6SStefano Zampini #undef __FUNCT__ 900b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 901b7ce53b6SStefano Zampini /*@ 902b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 903b7ce53b6SStefano Zampini 904b7ce53b6SStefano Zampini Input Parameter: 905b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 906b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 907b7ce53b6SStefano Zampini 908b7ce53b6SStefano Zampini Output Parameter: 909b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 910b7ce53b6SStefano Zampini 911b7ce53b6SStefano Zampini Level: developer 912b7ce53b6SStefano Zampini 913eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 914b7ce53b6SStefano Zampini 915b7ce53b6SStefano Zampini .seealso: MATIS 916b7ce53b6SStefano Zampini @*/ 917b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 918b7ce53b6SStefano Zampini { 919b7ce53b6SStefano Zampini PetscErrorCode ierr; 920b7ce53b6SStefano Zampini 921b7ce53b6SStefano Zampini PetscFunctionBegin; 922b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 923b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 924b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 925b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 926b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 927b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 9286c4ed002SBarry Smith if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 929b7ce53b6SStefano Zampini } 930b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 931b7ce53b6SStefano Zampini PetscFunctionReturn(0); 932b7ce53b6SStefano Zampini } 933b7ce53b6SStefano Zampini 934b7ce53b6SStefano Zampini #undef __FUNCT__ 935ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 936ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 937ad6194a2SStefano Zampini { 938ad6194a2SStefano Zampini PetscErrorCode ierr; 939ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 940ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 941ad6194a2SStefano Zampini Mat B,localmat; 942ad6194a2SStefano Zampini 943ad6194a2SStefano Zampini PetscFunctionBegin; 944ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 945ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 946ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 947e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 948ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 949ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 950b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 951ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 952ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 953ad6194a2SStefano Zampini *newmat = B; 954ad6194a2SStefano Zampini PetscFunctionReturn(0); 955ad6194a2SStefano Zampini } 956ad6194a2SStefano Zampini 957ad6194a2SStefano Zampini #undef __FUNCT__ 95869796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 959a8116848SStefano Zampini static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 96069796d55SStefano Zampini { 96169796d55SStefano Zampini PetscErrorCode ierr; 96269796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 96369796d55SStefano Zampini PetscBool local_sym; 96469796d55SStefano Zampini 96569796d55SStefano Zampini PetscFunctionBegin; 96669796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 967b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 96869796d55SStefano Zampini PetscFunctionReturn(0); 96969796d55SStefano Zampini } 97069796d55SStefano Zampini 97169796d55SStefano Zampini #undef __FUNCT__ 97269796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 973a8116848SStefano Zampini static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 97469796d55SStefano Zampini { 97569796d55SStefano Zampini PetscErrorCode ierr; 97669796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 97769796d55SStefano Zampini PetscBool local_sym; 97869796d55SStefano Zampini 97969796d55SStefano Zampini PetscFunctionBegin; 98069796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 981b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 98269796d55SStefano Zampini PetscFunctionReturn(0); 98369796d55SStefano Zampini } 98469796d55SStefano Zampini 98569796d55SStefano Zampini #undef __FUNCT__ 98645471136SStefano Zampini #define __FUNCT__ "MatIsStructurallySymmetric_IS" 98745471136SStefano Zampini static PetscErrorCode MatIsStructurallySymmetric_IS(Mat A,PetscBool *flg) 98845471136SStefano Zampini { 98945471136SStefano Zampini PetscErrorCode ierr; 99045471136SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 99145471136SStefano Zampini PetscBool local_sym; 99245471136SStefano Zampini 99345471136SStefano Zampini PetscFunctionBegin; 99445471136SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 99545471136SStefano Zampini *flg = PETSC_FALSE; 99645471136SStefano Zampini PetscFunctionReturn(0); 99745471136SStefano Zampini } 99845471136SStefano Zampini ierr = MatIsStructurallySymmetric(matis->A,&local_sym);CHKERRQ(ierr); 99945471136SStefano Zampini ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 100045471136SStefano Zampini PetscFunctionReturn(0); 100145471136SStefano Zampini } 100245471136SStefano Zampini 100345471136SStefano Zampini #undef __FUNCT__ 1004b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 1005a8116848SStefano Zampini static PetscErrorCode MatDestroy_IS(Mat A) 1006b4319ba4SBarry Smith { 1007dfbe8321SBarry Smith PetscErrorCode ierr; 1008b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 1009b4319ba4SBarry Smith 1010b4319ba4SBarry Smith PetscFunctionBegin; 10116bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 1012e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 1013e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 10146bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 10156bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 10163fd1c9e7SStefano Zampini ierr = VecDestroy(&b->counter);CHKERRQ(ierr); 1017a8116848SStefano Zampini ierr = ISDestroy(&b->getsub_ris);CHKERRQ(ierr); 1018a8116848SStefano Zampini ierr = ISDestroy(&b->getsub_cis);CHKERRQ(ierr); 1019a8116848SStefano Zampini if (b->sf != b->csf) { 1020a8116848SStefano Zampini ierr = PetscSFDestroy(&b->csf);CHKERRQ(ierr); 1021a8116848SStefano Zampini ierr = PetscFree2(b->csf_rootdata,b->csf_leafdata);CHKERRQ(ierr); 1022a8116848SStefano Zampini } 102328f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 102428f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 1025bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 1026dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 1027bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 1028b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 1029b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 10302e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 1031cf0a3239SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",NULL);CHKERRQ(ierr); 1032b4319ba4SBarry Smith PetscFunctionReturn(0); 1033b4319ba4SBarry Smith } 1034b4319ba4SBarry Smith 1035b4319ba4SBarry Smith #undef __FUNCT__ 1036b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 1037a8116848SStefano Zampini static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 1038b4319ba4SBarry Smith { 1039dfbe8321SBarry Smith PetscErrorCode ierr; 1040b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1041b4319ba4SBarry Smith PetscScalar zero = 0.0; 1042b4319ba4SBarry Smith 1043b4319ba4SBarry Smith PetscFunctionBegin; 1044b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 1045e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1046e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1047b4319ba4SBarry Smith 1048b4319ba4SBarry Smith /* multiply the local matrix */ 1049b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 1050b4319ba4SBarry Smith 1051b4319ba4SBarry Smith /* scatter product back into global memory */ 10522dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 1053e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1054e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1055b4319ba4SBarry Smith PetscFunctionReturn(0); 1056b4319ba4SBarry Smith } 1057b4319ba4SBarry Smith 1058b4319ba4SBarry Smith #undef __FUNCT__ 10592e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 1060a8116848SStefano Zampini static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 10612e74eeadSLisandro Dalcin { 1062650997f4SStefano Zampini Vec temp_vec; 10632e74eeadSLisandro Dalcin PetscErrorCode ierr; 10642e74eeadSLisandro Dalcin 10652e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 1066650997f4SStefano Zampini if (v3 != v2) { 1067650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 1068650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 1069650997f4SStefano Zampini } else { 1070650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 1071650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 1072650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 1073650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 1074650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 1075650997f4SStefano Zampini } 10762e74eeadSLisandro Dalcin PetscFunctionReturn(0); 10772e74eeadSLisandro Dalcin } 10782e74eeadSLisandro Dalcin 10792e74eeadSLisandro Dalcin #undef __FUNCT__ 10802e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 1081a8116848SStefano Zampini static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 10822e74eeadSLisandro Dalcin { 10832e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 10842e74eeadSLisandro Dalcin PetscErrorCode ierr; 10852e74eeadSLisandro Dalcin 1086e176bc59SStefano Zampini PetscFunctionBegin; 10872e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 1088e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1089e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 10902e74eeadSLisandro Dalcin 10912e74eeadSLisandro Dalcin /* multiply the local matrix */ 1092e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 10932e74eeadSLisandro Dalcin 10942e74eeadSLisandro Dalcin /* scatter product back into global vector */ 1095e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 1096e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1097e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 10982e74eeadSLisandro Dalcin PetscFunctionReturn(0); 10992e74eeadSLisandro Dalcin } 11002e74eeadSLisandro Dalcin 11012e74eeadSLisandro Dalcin #undef __FUNCT__ 11022e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 1103a8116848SStefano Zampini static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 11042e74eeadSLisandro Dalcin { 1105650997f4SStefano Zampini Vec temp_vec; 11062e74eeadSLisandro Dalcin PetscErrorCode ierr; 11072e74eeadSLisandro Dalcin 11082e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 1109650997f4SStefano Zampini if (v3 != v2) { 1110650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 1111650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 1112650997f4SStefano Zampini } else { 1113650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 1114650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 1115650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 1116650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 1117650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 1118650997f4SStefano Zampini } 11192e74eeadSLisandro Dalcin PetscFunctionReturn(0); 11202e74eeadSLisandro Dalcin } 11212e74eeadSLisandro Dalcin 11222e74eeadSLisandro Dalcin #undef __FUNCT__ 1123b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 1124a8116848SStefano Zampini static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 1125b4319ba4SBarry Smith { 1126b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 1127dfbe8321SBarry Smith PetscErrorCode ierr; 1128b4319ba4SBarry Smith PetscViewer sviewer; 1129ee2491ecSStefano Zampini PetscBool isascii,view = PETSC_TRUE; 1130b4319ba4SBarry Smith 1131b4319ba4SBarry Smith PetscFunctionBegin; 1132ee2491ecSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 1133ee2491ecSStefano Zampini if (isascii) { 1134ee2491ecSStefano Zampini PetscViewerFormat format; 1135ee2491ecSStefano Zampini 1136ee2491ecSStefano Zampini ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1137ee2491ecSStefano Zampini if (format == PETSC_VIEWER_ASCII_INFO) view = PETSC_FALSE; 1138ee2491ecSStefano Zampini } 1139ee2491ecSStefano Zampini if (!view) PetscFunctionReturn(0); 11403f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 1141b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 11423f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 11436e520ac8SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1144b4319ba4SBarry Smith PetscFunctionReturn(0); 1145b4319ba4SBarry Smith } 1146b4319ba4SBarry Smith 1147b4319ba4SBarry Smith #undef __FUNCT__ 1148b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 1149a8116848SStefano Zampini static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 1150b4319ba4SBarry Smith { 1151dfbe8321SBarry Smith PetscErrorCode ierr; 1152e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 1153b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1154b4319ba4SBarry Smith IS from,to; 1155e176bc59SStefano Zampini Vec cglobal,rglobal; 1156b4319ba4SBarry Smith 1157b4319ba4SBarry Smith PetscFunctionBegin; 1158784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 1159e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 11603bbff08aSStefano Zampini /* Destroy any previous data */ 116170cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 116270cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 11633fd1c9e7SStefano Zampini ierr = VecDestroy(&is->counter);CHKERRQ(ierr); 1164e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 1165e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 11661c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 116728f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 116828f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 11693bbff08aSStefano Zampini 11703bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 1171fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1172fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1173fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 1174fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 1175b4319ba4SBarry Smith 1176b4319ba4SBarry Smith /* Create the local matrix A */ 1177e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 1178e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 1179e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 1180e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 1181f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 1182e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 1183e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 1184e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 1185ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 1186ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 1187b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 1188b4319ba4SBarry Smith 1189f26d0771SStefano Zampini if (!is->islocalref) { /* setup scatters and local vectors for MatMult */ 1190b4319ba4SBarry Smith /* Create the local work vectors */ 11913bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 1192b4319ba4SBarry Smith 1193e176bc59SStefano Zampini /* setup the global to local scatters */ 1194e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 1195e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 1196784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 1197e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 1198e176bc59SStefano Zampini if (rmapping != cmapping) { 1199e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 1200e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 1201e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 1202e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 1203e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 1204e176bc59SStefano Zampini } else { 1205e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 1206e176bc59SStefano Zampini is->cctx = is->rctx; 1207e176bc59SStefano Zampini } 12083fd1c9e7SStefano Zampini 12093fd1c9e7SStefano Zampini /* interface counter vector (local) */ 12103fd1c9e7SStefano Zampini ierr = VecDuplicate(is->y,&is->counter);CHKERRQ(ierr); 12113fd1c9e7SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 12123fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12133fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 12143fd1c9e7SStefano Zampini ierr = VecScatterBegin(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12153fd1c9e7SStefano Zampini ierr = VecScatterEnd(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 12163fd1c9e7SStefano Zampini 12173fd1c9e7SStefano Zampini /* free workspace */ 1218e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 1219e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 12206bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 12216bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1222f26d0771SStefano Zampini } 12239c0b00dbSStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 1224b4319ba4SBarry Smith PetscFunctionReturn(0); 1225b4319ba4SBarry Smith } 1226b4319ba4SBarry Smith 12272e74eeadSLisandro Dalcin #undef __FUNCT__ 12282e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 1229a8116848SStefano Zampini static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 12302e74eeadSLisandro Dalcin { 12312e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 12322e74eeadSLisandro Dalcin PetscErrorCode ierr; 123397563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 123497563a80SStefano Zampini PetscInt i,zm,zn; 123597563a80SStefano Zampini #endif 1236f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 12372e74eeadSLisandro Dalcin 12382e74eeadSLisandro Dalcin PetscFunctionBegin; 12392e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 1240f26d0771SStefano 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); 124197563a80SStefano Zampini /* count negative indices */ 124297563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 124397563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 12442e74eeadSLisandro Dalcin #endif 124597563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 124697563a80SStefano Zampini ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 124797563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 124897563a80SStefano Zampini /* count negative indices (should be the same as before) */ 124997563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 125097563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 125197563a80SStefano 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"); 125297563a80SStefano 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"); 125397563a80SStefano Zampini #endif 12542e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 12552e74eeadSLisandro Dalcin PetscFunctionReturn(0); 12562e74eeadSLisandro Dalcin } 12572e74eeadSLisandro Dalcin 1258b4319ba4SBarry Smith #undef __FUNCT__ 125997563a80SStefano Zampini #define __FUNCT__ "MatSetValuesBlocked_IS" 1260a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 126197563a80SStefano Zampini { 126297563a80SStefano Zampini Mat_IS *is = (Mat_IS*)mat->data; 126397563a80SStefano Zampini PetscErrorCode ierr; 126497563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 126597563a80SStefano Zampini PetscInt i,zm,zn; 126697563a80SStefano Zampini #endif 1267f26d0771SStefano Zampini PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 126897563a80SStefano Zampini 126997563a80SStefano Zampini PetscFunctionBegin; 127097563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 1271f26d0771SStefano 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); 127297563a80SStefano Zampini /* count negative indices */ 127397563a80SStefano Zampini for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 127497563a80SStefano Zampini for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 127597563a80SStefano Zampini #endif 127697563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 127797563a80SStefano Zampini ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 127897563a80SStefano Zampini #if defined(PETSC_USE_DEBUG) 127997563a80SStefano Zampini /* count negative indices (should be the same as before) */ 128097563a80SStefano Zampini for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 128197563a80SStefano Zampini for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 128297563a80SStefano 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"); 128397563a80SStefano 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"); 128497563a80SStefano Zampini #endif 128597563a80SStefano Zampini ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 128697563a80SStefano Zampini PetscFunctionReturn(0); 128797563a80SStefano Zampini } 128897563a80SStefano Zampini 128997563a80SStefano Zampini #undef __FUNCT__ 1290b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 1291a8116848SStefano Zampini static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1292b4319ba4SBarry Smith { 1293dfbe8321SBarry Smith PetscErrorCode ierr; 1294b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1295b4319ba4SBarry Smith 1296b4319ba4SBarry Smith PetscFunctionBegin; 1297b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1298b4319ba4SBarry Smith PetscFunctionReturn(0); 1299b4319ba4SBarry Smith } 1300b4319ba4SBarry Smith 1301b4319ba4SBarry Smith #undef __FUNCT__ 1302f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 1303a8116848SStefano Zampini static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1304f0006bf2SLisandro Dalcin { 1305f0006bf2SLisandro Dalcin PetscErrorCode ierr; 1306f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 1307f0006bf2SLisandro Dalcin 1308f0006bf2SLisandro Dalcin PetscFunctionBegin; 1309f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1310f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 1311f0006bf2SLisandro Dalcin } 1312f0006bf2SLisandro Dalcin 1313f0006bf2SLisandro Dalcin #undef __FUNCT__ 1314f0ae7da4SStefano Zampini #define __FUNCT__ "MatISZeroRowsColumnsLocal_Private" 1315f0ae7da4SStefano Zampini static PetscErrorCode MatISZeroRowsColumnsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,PetscBool columns) 13162e74eeadSLisandro Dalcin { 1317f0ae7da4SStefano Zampini Mat_IS *is = (Mat_IS*)A->data; 13182e74eeadSLisandro Dalcin PetscErrorCode ierr; 13192e74eeadSLisandro Dalcin 13202e74eeadSLisandro Dalcin PetscFunctionBegin; 1321f0ae7da4SStefano Zampini if (!n) { 1322f0ae7da4SStefano Zampini is->pure_neumann = PETSC_TRUE; 1323f0ae7da4SStefano Zampini } else { 1324f0ae7da4SStefano Zampini PetscInt i; 1325f0ae7da4SStefano Zampini is->pure_neumann = PETSC_FALSE; 1326f0ae7da4SStefano Zampini 1327f0ae7da4SStefano Zampini if (columns) { 1328f0ae7da4SStefano Zampini ierr = MatZeroRowsColumns(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 1329f0ae7da4SStefano Zampini } else { 1330f0ae7da4SStefano Zampini ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 13312e74eeadSLisandro Dalcin } 1332f0ae7da4SStefano Zampini if (diag != 0.) { 1333f0ae7da4SStefano Zampini const PetscScalar *array; 1334f0ae7da4SStefano Zampini ierr = VecGetArrayRead(is->counter,&array);CHKERRQ(ierr); 1335f0ae7da4SStefano Zampini for (i=0; i<n; i++) { 1336f0ae7da4SStefano Zampini ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 1337f0ae7da4SStefano Zampini } 1338f0ae7da4SStefano Zampini ierr = VecRestoreArrayRead(is->counter,&array);CHKERRQ(ierr); 1339f0ae7da4SStefano Zampini } 1340f0ae7da4SStefano Zampini ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1341f0ae7da4SStefano Zampini ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1342f0ae7da4SStefano Zampini } 13432e74eeadSLisandro Dalcin PetscFunctionReturn(0); 13442e74eeadSLisandro Dalcin } 13452e74eeadSLisandro Dalcin 13462e74eeadSLisandro Dalcin #undef __FUNCT__ 1347f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_Private_IS" 1348f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_Private_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b,PetscBool columns) 1349b4319ba4SBarry Smith { 13506e520ac8SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 13516e520ac8SStefano Zampini PetscInt nr,nl,len,i; 13526e520ac8SStefano Zampini PetscInt *lrows; 1353dfbe8321SBarry Smith PetscErrorCode ierr; 1354b4319ba4SBarry Smith 1355b4319ba4SBarry Smith PetscFunctionBegin; 1356f0ae7da4SStefano Zampini #if defined(PETSC_USE_DEBUG) 1357f0ae7da4SStefano Zampini if (columns || diag != 0. || (x && b)) { 1358f0ae7da4SStefano Zampini PetscBool cong; 1359f0ae7da4SStefano Zampini ierr = PetscLayoutCompare(A->rmap,A->cmap,&cong);CHKERRQ(ierr); 1360f0ae7da4SStefano 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"); 1361f0ae7da4SStefano 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"); 1362f0ae7da4SStefano Zampini if (!cong && x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"A->rmap and A->cmap need to be congruent"); 1363b4319ba4SBarry Smith } 1364f0ae7da4SStefano Zampini #endif 13656e520ac8SStefano Zampini /* get locally owned rows */ 1366f0ae7da4SStefano Zampini ierr = PetscLayoutMapLocal_Private(A->rmap,n,rows,&len,&lrows,NULL);CHKERRQ(ierr); 13676e520ac8SStefano Zampini /* fix right hand side if needed */ 13686e520ac8SStefano Zampini if (x && b) { 13696e520ac8SStefano Zampini const PetscScalar *xx; 13706e520ac8SStefano Zampini PetscScalar *bb; 13712205254eSKarl Rupp 13726e520ac8SStefano Zampini ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr); 13736e520ac8SStefano Zampini ierr = VecGetArray(b, &bb);CHKERRQ(ierr); 13746e520ac8SStefano Zampini for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]]; 13756e520ac8SStefano Zampini ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr); 13766e520ac8SStefano Zampini ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr); 1377b4319ba4SBarry Smith } 13786e520ac8SStefano Zampini /* get rows associated to the local matrices */ 13793d996552SStefano Zampini ierr = MatISSetUpSF(A);CHKERRQ(ierr); 13806e520ac8SStefano Zampini ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr); 13816e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr); 13826e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 13836e520ac8SStefano Zampini for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1; 13846e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 13856e520ac8SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 13866e520ac8SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 13876e520ac8SStefano Zampini ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr); 13886e520ac8SStefano Zampini for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i; 1389f0ae7da4SStefano Zampini ierr = MatISZeroRowsColumnsLocal_Private(A,nr,lrows,diag,columns);CHKERRQ(ierr); 13906e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 1391b4319ba4SBarry Smith PetscFunctionReturn(0); 1392b4319ba4SBarry Smith } 1393b4319ba4SBarry Smith 1394b4319ba4SBarry Smith #undef __FUNCT__ 1395f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRows_IS" 1396f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 1397b4319ba4SBarry Smith { 1398b4319ba4SBarry Smith PetscErrorCode ierr; 1399b4319ba4SBarry Smith 1400b4319ba4SBarry Smith PetscFunctionBegin; 1401f0ae7da4SStefano Zampini ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_FALSE);CHKERRQ(ierr); 1402f0ae7da4SStefano Zampini PetscFunctionReturn(0); 1403f0ae7da4SStefano Zampini } 1404b4319ba4SBarry Smith 1405f0ae7da4SStefano Zampini #undef __FUNCT__ 1406f0ae7da4SStefano Zampini #define __FUNCT__ "MatZeroRowsColumns_IS" 1407f0ae7da4SStefano Zampini static PetscErrorCode MatZeroRowsColumns_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 1408f0ae7da4SStefano Zampini { 1409f0ae7da4SStefano Zampini PetscErrorCode ierr; 1410f0ae7da4SStefano Zampini 1411f0ae7da4SStefano Zampini PetscFunctionBegin; 1412f0ae7da4SStefano Zampini ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_TRUE);CHKERRQ(ierr); 1413b4319ba4SBarry Smith PetscFunctionReturn(0); 1414b4319ba4SBarry Smith } 1415b4319ba4SBarry Smith 1416b4319ba4SBarry Smith #undef __FUNCT__ 1417b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 1418a8116848SStefano Zampini static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 1419b4319ba4SBarry Smith { 1420b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1421dfbe8321SBarry Smith PetscErrorCode ierr; 1422b4319ba4SBarry Smith 1423b4319ba4SBarry Smith PetscFunctionBegin; 1424b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 1425b4319ba4SBarry Smith PetscFunctionReturn(0); 1426b4319ba4SBarry Smith } 1427b4319ba4SBarry Smith 1428b4319ba4SBarry Smith #undef __FUNCT__ 1429b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 1430a8116848SStefano Zampini static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 1431b4319ba4SBarry Smith { 1432b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 1433dfbe8321SBarry Smith PetscErrorCode ierr; 1434b4319ba4SBarry Smith 1435b4319ba4SBarry Smith PetscFunctionBegin; 1436b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 1437b4319ba4SBarry Smith PetscFunctionReturn(0); 1438b4319ba4SBarry Smith } 1439b4319ba4SBarry Smith 1440b4319ba4SBarry Smith #undef __FUNCT__ 1441b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 1442a8116848SStefano Zampini static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 1443b4319ba4SBarry Smith { 1444b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 1445b4319ba4SBarry Smith 1446b4319ba4SBarry Smith PetscFunctionBegin; 1447b4319ba4SBarry Smith *local = is->A; 1448b4319ba4SBarry Smith PetscFunctionReturn(0); 1449b4319ba4SBarry Smith } 1450b4319ba4SBarry Smith 1451b4319ba4SBarry Smith #undef __FUNCT__ 1452b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 1453b4319ba4SBarry Smith /*@ 1454b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 1455b4319ba4SBarry Smith 1456b4319ba4SBarry Smith Input Parameter: 1457b4319ba4SBarry Smith . mat - the matrix 1458b4319ba4SBarry Smith 1459b4319ba4SBarry Smith Output Parameter: 1460eb82efa4SStefano Zampini . local - the local matrix 1461b4319ba4SBarry Smith 1462b4319ba4SBarry Smith Level: advanced 1463b4319ba4SBarry Smith 1464b4319ba4SBarry Smith Notes: 1465b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 1466b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 1467b4319ba4SBarry Smith of the MatSetValues() operation. 1468b4319ba4SBarry Smith 1469b4319ba4SBarry Smith .seealso: MATIS 1470b4319ba4SBarry Smith @*/ 14717087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 1472b4319ba4SBarry Smith { 14734ac538c5SBarry Smith PetscErrorCode ierr; 1474b4319ba4SBarry Smith 1475b4319ba4SBarry Smith PetscFunctionBegin; 14760700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1477b4319ba4SBarry Smith PetscValidPointer(local,2); 14784ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 1479b4319ba4SBarry Smith PetscFunctionReturn(0); 1480b4319ba4SBarry Smith } 1481b4319ba4SBarry Smith 14823b03a366Sstefano_zampini #undef __FUNCT__ 14833b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 1484a8116848SStefano Zampini static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 14853b03a366Sstefano_zampini { 14863b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 14873b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 14883b03a366Sstefano_zampini PetscErrorCode ierr; 14893b03a366Sstefano_zampini 14903b03a366Sstefano_zampini PetscFunctionBegin; 14914e4c7dbeSStefano Zampini if (is->A) { 14923b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 14933b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 1494f0ae7da4SStefano 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); 14954e4c7dbeSStefano Zampini } 14963b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 14973b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 14983b03a366Sstefano_zampini is->A = local; 14993b03a366Sstefano_zampini PetscFunctionReturn(0); 15003b03a366Sstefano_zampini } 15013b03a366Sstefano_zampini 15023b03a366Sstefano_zampini #undef __FUNCT__ 15033b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 15043b03a366Sstefano_zampini /*@ 1505eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 15063b03a366Sstefano_zampini 15073b03a366Sstefano_zampini Input Parameter: 15083b03a366Sstefano_zampini . mat - the matrix 1509eb82efa4SStefano Zampini . local - the local matrix 15103b03a366Sstefano_zampini 15113b03a366Sstefano_zampini Output Parameter: 15123b03a366Sstefano_zampini 15133b03a366Sstefano_zampini Level: advanced 15143b03a366Sstefano_zampini 15153b03a366Sstefano_zampini Notes: 15163b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 15173b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 15183b03a366Sstefano_zampini 15193b03a366Sstefano_zampini .seealso: MATIS 15203b03a366Sstefano_zampini @*/ 15213b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 15223b03a366Sstefano_zampini { 15233b03a366Sstefano_zampini PetscErrorCode ierr; 15243b03a366Sstefano_zampini 15253b03a366Sstefano_zampini PetscFunctionBegin; 15263b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1527b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 15283b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 15293b03a366Sstefano_zampini PetscFunctionReturn(0); 15303b03a366Sstefano_zampini } 15313b03a366Sstefano_zampini 15326726f965SBarry Smith #undef __FUNCT__ 15336726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 1534a8116848SStefano Zampini static PetscErrorCode MatZeroEntries_IS(Mat A) 15356726f965SBarry Smith { 15366726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 15376726f965SBarry Smith PetscErrorCode ierr; 15386726f965SBarry Smith 15396726f965SBarry Smith PetscFunctionBegin; 15406726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 15416726f965SBarry Smith PetscFunctionReturn(0); 15426726f965SBarry Smith } 15436726f965SBarry Smith 15446726f965SBarry Smith #undef __FUNCT__ 15452e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 1546a8116848SStefano Zampini static PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 15472e74eeadSLisandro Dalcin { 15482e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 15492e74eeadSLisandro Dalcin PetscErrorCode ierr; 15502e74eeadSLisandro Dalcin 15512e74eeadSLisandro Dalcin PetscFunctionBegin; 15522e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 15532e74eeadSLisandro Dalcin PetscFunctionReturn(0); 15542e74eeadSLisandro Dalcin } 15552e74eeadSLisandro Dalcin 15562e74eeadSLisandro Dalcin #undef __FUNCT__ 15572e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 1558a8116848SStefano Zampini static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 15592e74eeadSLisandro Dalcin { 15602e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 15612e74eeadSLisandro Dalcin PetscErrorCode ierr; 15622e74eeadSLisandro Dalcin 15632e74eeadSLisandro Dalcin PetscFunctionBegin; 15642e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 1565e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 15662e74eeadSLisandro Dalcin 15672e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 15682e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 1569e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1570e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 15712e74eeadSLisandro Dalcin PetscFunctionReturn(0); 15722e74eeadSLisandro Dalcin } 15732e74eeadSLisandro Dalcin 15742e74eeadSLisandro Dalcin #undef __FUNCT__ 15756726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 1576a8116848SStefano Zampini static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 15776726f965SBarry Smith { 15786726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 15796726f965SBarry Smith PetscErrorCode ierr; 15806726f965SBarry Smith 15816726f965SBarry Smith PetscFunctionBegin; 15824e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 15836726f965SBarry Smith PetscFunctionReturn(0); 15846726f965SBarry Smith } 15856726f965SBarry Smith 1586284134d9SBarry Smith #undef __FUNCT__ 1587f26d0771SStefano Zampini #define __FUNCT__ "MatAXPY_IS" 1588f26d0771SStefano Zampini static PetscErrorCode MatAXPY_IS(Mat Y,PetscScalar a,Mat X,MatStructure str) 1589f26d0771SStefano Zampini { 1590f26d0771SStefano Zampini Mat_IS *y = (Mat_IS*)Y->data; 1591f26d0771SStefano Zampini Mat_IS *x; 1592f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1593f26d0771SStefano Zampini PetscBool ismatis; 1594f26d0771SStefano Zampini #endif 1595f26d0771SStefano Zampini PetscErrorCode ierr; 1596f26d0771SStefano Zampini 1597f26d0771SStefano Zampini PetscFunctionBegin; 1598f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1599f26d0771SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)X,MATIS,&ismatis);CHKERRQ(ierr); 1600f26d0771SStefano Zampini if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS"); 1601f26d0771SStefano Zampini #endif 1602f26d0771SStefano Zampini x = (Mat_IS*)X->data; 1603f26d0771SStefano Zampini ierr = MatAXPY(y->A,a,x->A,str);CHKERRQ(ierr); 1604f26d0771SStefano Zampini PetscFunctionReturn(0); 1605f26d0771SStefano Zampini } 1606f26d0771SStefano Zampini 1607f26d0771SStefano Zampini #undef __FUNCT__ 1608f26d0771SStefano Zampini #define __FUNCT__ "MatGetLocalSubMatrix_IS" 1609f26d0771SStefano Zampini static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat) 1610f26d0771SStefano Zampini { 1611f26d0771SStefano Zampini Mat lA; 1612f26d0771SStefano Zampini Mat_IS *matis; 1613f26d0771SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 1614f26d0771SStefano Zampini IS is; 1615f26d0771SStefano Zampini const PetscInt *rg,*rl; 1616f26d0771SStefano Zampini PetscInt nrg; 1617f26d0771SStefano Zampini PetscInt N,M,nrl,i,*idxs; 1618f26d0771SStefano Zampini PetscErrorCode ierr; 1619f26d0771SStefano Zampini 1620f26d0771SStefano Zampini PetscFunctionBegin; 1621f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 1622f26d0771SStefano Zampini ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr); 1623f26d0771SStefano Zampini ierr = ISGetIndices(row,&rl);CHKERRQ(ierr); 1624f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr); 1625f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1626f0ae7da4SStefano 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); 1627f26d0771SStefano Zampini #endif 1628f26d0771SStefano Zampini ierr = PetscMalloc1(nrg,&idxs);CHKERRQ(ierr); 1629f26d0771SStefano Zampini /* map from [0,nrl) to row */ 1630f26d0771SStefano Zampini for (i=0;i<nrl;i++) idxs[i] = rl[i]; 1631f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1632f26d0771SStefano Zampini for (i=nrl;i<nrg;i++) idxs[i] = nrg; 1633f26d0771SStefano Zampini #else 1634f26d0771SStefano Zampini for (i=nrl;i<nrg;i++) idxs[i] = -1; 1635f26d0771SStefano Zampini #endif 1636f26d0771SStefano Zampini ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr); 1637f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 1638f26d0771SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 1639f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 1640f26d0771SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 1641f26d0771SStefano Zampini /* compute new l2g map for columns */ 1642f26d0771SStefano Zampini if (col != row || A->rmap->mapping != A->cmap->mapping) { 1643f26d0771SStefano Zampini const PetscInt *cg,*cl; 1644f26d0771SStefano Zampini PetscInt ncg; 1645f26d0771SStefano Zampini PetscInt ncl; 1646f26d0771SStefano Zampini 1647f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 1648f26d0771SStefano Zampini ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr); 1649f26d0771SStefano Zampini ierr = ISGetIndices(col,&cl);CHKERRQ(ierr); 1650f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr); 1651f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1652f0ae7da4SStefano 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); 1653f26d0771SStefano Zampini #endif 1654f26d0771SStefano Zampini ierr = PetscMalloc1(ncg,&idxs);CHKERRQ(ierr); 1655f26d0771SStefano Zampini /* map from [0,ncl) to col */ 1656f26d0771SStefano Zampini for (i=0;i<ncl;i++) idxs[i] = cl[i]; 1657f26d0771SStefano Zampini #if defined(PETSC_USE_DEBUG) 1658f26d0771SStefano Zampini for (i=ncl;i<ncg;i++) idxs[i] = ncg; 1659f26d0771SStefano Zampini #else 1660f26d0771SStefano Zampini for (i=ncl;i<ncg;i++) idxs[i] = -1; 1661f26d0771SStefano Zampini #endif 1662f26d0771SStefano Zampini ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr); 1663f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 1664f26d0771SStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 1665f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 1666f26d0771SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 1667f26d0771SStefano Zampini } else { 1668f26d0771SStefano Zampini ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr); 1669f26d0771SStefano Zampini cl2g = rl2g; 1670f26d0771SStefano Zampini } 1671f26d0771SStefano Zampini /* create the MATIS submatrix */ 1672f26d0771SStefano Zampini ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 1673f26d0771SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr); 1674f26d0771SStefano Zampini ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1675f26d0771SStefano Zampini ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr); 1676b0aa3428SStefano Zampini matis = (Mat_IS*)((*submat)->data); 1677f26d0771SStefano Zampini matis->islocalref = PETSC_TRUE; 1678f26d0771SStefano Zampini ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr); 1679f26d0771SStefano Zampini ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr); 1680f26d0771SStefano Zampini ierr = MatISSetLocalMat(*submat,lA);CHKERRQ(ierr); 1681f26d0771SStefano Zampini ierr = MatSetUp(*submat);CHKERRQ(ierr); 1682f26d0771SStefano Zampini ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1683f26d0771SStefano Zampini ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1684f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 1685f26d0771SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 1686f26d0771SStefano Zampini /* remove unsupported ops */ 1687f26d0771SStefano Zampini ierr = PetscMemzero((*submat)->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1688f26d0771SStefano Zampini (*submat)->ops->destroy = MatDestroy_IS; 1689f26d0771SStefano Zampini (*submat)->ops->setvalueslocal = MatSetValuesLocal_SubMat_IS; 1690f26d0771SStefano Zampini (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS; 1691f26d0771SStefano Zampini (*submat)->ops->assemblybegin = MatAssemblyBegin_IS; 1692f26d0771SStefano Zampini (*submat)->ops->assemblyend = MatAssemblyEnd_IS; 1693f26d0771SStefano Zampini PetscFunctionReturn(0); 1694f26d0771SStefano Zampini } 1695f26d0771SStefano Zampini 1696f26d0771SStefano Zampini #undef __FUNCT__ 1697284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 1698284134d9SBarry Smith /*@ 16993c212e90SHong Zhang MatCreateIS - Creates a "process" unassembled matrix, assembled on each 1700284134d9SBarry Smith process but not across processes. 1701284134d9SBarry Smith 1702284134d9SBarry Smith Input Parameters: 1703284134d9SBarry Smith + comm - MPI communicator that will share the matrix 1704e176bc59SStefano Zampini . bs - block size of the matrix 1705df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 1706e176bc59SStefano Zampini . rmap - local to global map for rows 1707e176bc59SStefano Zampini - cmap - local to global map for cols 1708284134d9SBarry Smith 1709284134d9SBarry Smith Output Parameter: 1710284134d9SBarry Smith . A - the resulting matrix 1711284134d9SBarry Smith 17128e6c10adSSatish Balay Level: advanced 17138e6c10adSSatish Balay 17143c212e90SHong Zhang Notes: See MATIS for more details. 17153c212e90SHong Zhang m and n are NOT related to the size of the map; they are the size of the part of the vector owned 1716e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 17173c212e90SHong Zhang If either rmap or cmap are NULL, then the matrix is assumed to be square. 1718284134d9SBarry Smith 1719284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 1720284134d9SBarry Smith @*/ 1721e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 1722284134d9SBarry Smith { 1723284134d9SBarry Smith PetscErrorCode ierr; 1724284134d9SBarry Smith 1725284134d9SBarry Smith PetscFunctionBegin; 1726e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 1727284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1728d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 1729284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1730284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 1731e176bc59SStefano Zampini if (rmap && cmap) { 1732e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 1733e176bc59SStefano Zampini } else if (!rmap) { 1734e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 1735e176bc59SStefano Zampini } else { 1736e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 1737e176bc59SStefano Zampini } 1738284134d9SBarry Smith PetscFunctionReturn(0); 1739284134d9SBarry Smith } 1740284134d9SBarry Smith 1741b4319ba4SBarry Smith /*MC 1742f26d0771SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition methods (e.g. PCBDDC or KSPFETIDP). 1743b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 1744b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 1745b4319ba4SBarry Smith product is handled "implicitly". 1746b4319ba4SBarry Smith 1747b4319ba4SBarry Smith Operations Provided: 17486726f965SBarry Smith + MatMult() 17492e74eeadSLisandro Dalcin . MatMultAdd() 17502e74eeadSLisandro Dalcin . MatMultTranspose() 17512e74eeadSLisandro Dalcin . MatMultTransposeAdd() 17526726f965SBarry Smith . MatZeroEntries() 17536726f965SBarry Smith . MatSetOption() 17542e74eeadSLisandro Dalcin . MatZeroRows() 17552e74eeadSLisandro Dalcin . MatSetValues() 175697563a80SStefano Zampini . MatSetValuesBlocked() 17576726f965SBarry Smith . MatSetValuesLocal() 175897563a80SStefano Zampini . MatSetValuesBlockedLocal() 17592e74eeadSLisandro Dalcin . MatScale() 17602e74eeadSLisandro Dalcin . MatGetDiagonal() 17612b404112SStefano Zampini . MatMissingDiagonal() 17622b404112SStefano Zampini . MatDuplicate() 17632b404112SStefano Zampini . MatCopy() 1764f26d0771SStefano Zampini . MatAXPY() 1765f26d0771SStefano Zampini . MatGetSubMatrix() 1766f26d0771SStefano Zampini . MatGetLocalSubMatrix() 17676726f965SBarry Smith - MatSetLocalToGlobalMapping() 1768b4319ba4SBarry Smith 1769b4319ba4SBarry Smith Options Database Keys: 1770b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1771b4319ba4SBarry Smith 1772b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1773b4319ba4SBarry Smith 1774b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1775b4319ba4SBarry Smith 1776b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1777eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1778b4319ba4SBarry Smith 1779b4319ba4SBarry Smith Level: advanced 1780b4319ba4SBarry Smith 1781f26d0771SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), MatCreateIS(), PCBDDC, KSPFETIDP 1782b4319ba4SBarry Smith 1783b4319ba4SBarry Smith M*/ 1784b4319ba4SBarry Smith 1785b4319ba4SBarry Smith #undef __FUNCT__ 1786b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 17878cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1788b4319ba4SBarry Smith { 1789dfbe8321SBarry Smith PetscErrorCode ierr; 1790b4319ba4SBarry Smith Mat_IS *b; 1791b4319ba4SBarry Smith 1792b4319ba4SBarry Smith PetscFunctionBegin; 1793b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1794b4319ba4SBarry Smith A->data = (void*)b; 1795b4319ba4SBarry Smith 1796e176bc59SStefano Zampini /* matrix ops */ 1797e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1798b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 17992e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 18002e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 18012e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1802b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1803b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 18042e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 180598921651SStefano Zampini A->ops->setvaluesblocked = MatSetValuesBlocked_IS; 1806b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1807f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 18082e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1809f0ae7da4SStefano Zampini A->ops->zerorowscolumns = MatZeroRowsColumns_IS; 1810b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1811b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1812b4319ba4SBarry Smith A->ops->view = MatView_IS; 18136726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 18142e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 18152e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 18166726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 181769796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 181869796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 181945471136SStefano Zampini A->ops->isstructurallysymmetric = MatIsStructurallySymmetric_IS; 1820ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 18216bd84002SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_IS; 18222b404112SStefano Zampini A->ops->copy = MatCopy_IS; 1823659959c5SStefano Zampini A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_IS; 1824a8116848SStefano Zampini A->ops->getsubmatrix = MatGetSubMatrix_IS; 1825f26d0771SStefano Zampini A->ops->axpy = MatAXPY_IS; 18263fd1c9e7SStefano Zampini A->ops->diagonalset = MatDiagonalSet_IS; 18273fd1c9e7SStefano Zampini A->ops->shift = MatShift_IS; 1828b4319ba4SBarry Smith 1829b7ce53b6SStefano Zampini /* special MATIS functions */ 1830bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1831bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1832b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 18332e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 1834cf0a3239SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",MatISSetUpSF_IS);CHKERRQ(ierr); 183517667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1836b4319ba4SBarry Smith PetscFunctionReturn(0); 1837b4319ba4SBarry Smith } 1838