1 /* 2 Creates a matrix class for using the Neumann-Neumann type preconditioners. 3 This stores the matrices in globally unassembled form. Each processor 4 assembles only its local Neumann problem and the parallel matrix vector 5 product is handled "implicitly". 6 7 Currently this allows for only one subdomain per processor. 8 */ 9 10 #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/ 11 #include <petsc/private/sfimpl.h> 12 13 #define MATIS_MAX_ENTRIES_INSERTION 2048 14 15 #undef __FUNCT__ 16 #define __FUNCT__ "MatISSetUpSF" 17 /*@ 18 MatISSetUpSF - Setup star forest objects used by MatIS. 19 20 Collective on MPI_Comm 21 22 Input Parameters: 23 + A - the matrix 24 25 Level: advanced 26 27 Notes: This function does not need to be called by the user. 28 29 .keywords: matrix 30 31 .seealso: MatCreate(), MatCreateIS(), MatISSetPreallocation(), MatISGetLocalMat() 32 @*/ 33 PetscErrorCode MatISSetUpSF(Mat A) 34 { 35 PetscErrorCode ierr; 36 37 PetscFunctionBegin; 38 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 39 PetscValidType(A,1); 40 ierr = PetscTryMethod(A,"MatISSetUpSF_C",(Mat),(A));CHKERRQ(ierr); 41 PetscFunctionReturn(0); 42 } 43 44 #undef __FUNCT__ 45 #define __FUNCT__ "MatDiagonalSet_IS" 46 PetscErrorCode MatDiagonalSet_IS(Mat A,Vec D,InsertMode insmode) 47 { 48 Mat_IS *is = (Mat_IS*)A->data; 49 PetscErrorCode ierr; 50 51 PetscFunctionBegin; 52 if (!D) { /* this code branch is used by MatShift_IS */ 53 ierr = VecSet(is->y,1.);CHKERRQ(ierr); 54 } else { 55 ierr = VecScatterBegin(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 56 ierr = VecScatterEnd(is->rctx,D,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 57 } 58 ierr = VecPointwiseDivide(is->y,is->y,is->counter);CHKERRQ(ierr); 59 ierr = MatDiagonalSet(is->A,is->y,insmode);CHKERRQ(ierr); 60 PetscFunctionReturn(0); 61 } 62 63 #undef __FUNCT__ 64 #define __FUNCT__ "MatShift_IS" 65 PetscErrorCode MatShift_IS(Mat A,PetscScalar a) 66 { 67 PetscErrorCode ierr; 68 69 PetscFunctionBegin; 70 ierr = MatDiagonalSet_IS(A,NULL,ADD_VALUES);CHKERRQ(ierr); 71 PetscFunctionReturn(0); 72 } 73 74 #undef __FUNCT__ 75 #define __FUNCT__ "MatSetValuesLocal_SubMat_IS" 76 static PetscErrorCode MatSetValuesLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 77 { 78 PetscErrorCode ierr; 79 Mat_IS *is = (Mat_IS*)A->data; 80 PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 81 82 PetscFunctionBegin; 83 #if defined(PETSC_USE_DEBUG) 84 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); 85 #endif 86 ierr = ISLocalToGlobalMappingApply(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 87 ierr = ISLocalToGlobalMappingApply(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 88 ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 89 PetscFunctionReturn(0); 90 } 91 92 #undef __FUNCT__ 93 #define __FUNCT__ "MatSetValuesBlockedLocal_SubMat_IS" 94 static PetscErrorCode MatSetValuesBlockedLocal_SubMat_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 95 { 96 PetscErrorCode ierr; 97 Mat_IS *is = (Mat_IS*)A->data; 98 PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 99 100 PetscFunctionBegin; 101 #if defined(PETSC_USE_DEBUG) 102 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); 103 #endif 104 ierr = ISLocalToGlobalMappingApplyBlock(A->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 105 ierr = ISLocalToGlobalMappingApplyBlock(A->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 106 ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 107 PetscFunctionReturn(0); 108 } 109 110 #undef __FUNCT__ 111 #define __FUNCT__ "PetscLayoutMapLocal_Private" 112 static PetscErrorCode PetscLayoutMapLocal_Private(PetscLayout map,PetscInt N,const PetscInt idxs[], PetscInt *on,PetscInt **oidxs,PetscInt **ogidxs) 113 { 114 PetscInt *owners = map->range; 115 PetscInt n = map->n; 116 PetscSF sf; 117 PetscInt *lidxs,*work = NULL; 118 PetscSFNode *ridxs; 119 PetscMPIInt rank; 120 PetscInt r, p = 0, len = 0; 121 PetscErrorCode ierr; 122 123 PetscFunctionBegin; 124 /* Create SF where leaves are input idxs and roots are owned idxs (code adapted from MatZeroRowsMapLocal_Private) */ 125 ierr = MPI_Comm_rank(map->comm,&rank);CHKERRQ(ierr); 126 ierr = PetscMalloc1(n,&lidxs);CHKERRQ(ierr); 127 for (r = 0; r < n; ++r) lidxs[r] = -1; 128 ierr = PetscMalloc1(N,&ridxs);CHKERRQ(ierr); 129 for (r = 0; r < N; ++r) { 130 const PetscInt idx = idxs[r]; 131 if (idx < 0 || map->N <= idx) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Index %D out of range [0,%D)",idx,map->N); 132 if (idx < owners[p] || owners[p+1] <= idx) { /* short-circuit the search if the last p owns this idx too */ 133 ierr = PetscLayoutFindOwner(map,idx,&p);CHKERRQ(ierr); 134 } 135 ridxs[r].rank = p; 136 ridxs[r].index = idxs[r] - owners[p]; 137 } 138 ierr = PetscSFCreate(map->comm,&sf);CHKERRQ(ierr); 139 ierr = PetscSFSetGraph(sf,n,N,NULL,PETSC_OWN_POINTER,ridxs,PETSC_OWN_POINTER);CHKERRQ(ierr); 140 ierr = PetscSFReduceBegin(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 141 ierr = PetscSFReduceEnd(sf,MPIU_INT,(PetscInt*)idxs,lidxs,MPI_LOR);CHKERRQ(ierr); 142 if (ogidxs) { /* communicate global idxs */ 143 PetscInt cum = 0,start,*work2; 144 145 ierr = PetscMalloc1(n,&work);CHKERRQ(ierr); 146 ierr = PetscCalloc1(N,&work2);CHKERRQ(ierr); 147 for (r = 0; r < N; ++r) if (idxs[r] >=0) cum++; 148 ierr = MPI_Scan(&cum,&start,1,MPIU_INT,MPI_SUM,map->comm);CHKERRQ(ierr); 149 start -= cum; 150 cum = 0; 151 for (r = 0; r < N; ++r) if (idxs[r] >=0) work2[r] = start+cum++; 152 ierr = PetscSFReduceBegin(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 153 ierr = PetscSFReduceEnd(sf,MPIU_INT,work2,work,MPIU_REPLACE);CHKERRQ(ierr); 154 ierr = PetscFree(work2);CHKERRQ(ierr); 155 } 156 ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 157 /* Compress and put in indices */ 158 for (r = 0; r < n; ++r) 159 if (lidxs[r] >= 0) { 160 if (work) work[len] = work[r]; 161 lidxs[len++] = r; 162 } 163 if (on) *on = len; 164 if (oidxs) *oidxs = lidxs; 165 if (ogidxs) *ogidxs = work; 166 PetscFunctionReturn(0); 167 } 168 169 #undef __FUNCT__ 170 #define __FUNCT__ "MatGetSubMatrix_IS" 171 static PetscErrorCode MatGetSubMatrix_IS(Mat mat,IS irow,IS icol,MatReuse scall,Mat *newmat) 172 { 173 Mat locmat,newlocmat; 174 Mat_IS *newmatis; 175 #if defined(PETSC_USE_DEBUG) 176 Vec rtest,ltest; 177 const PetscScalar *array; 178 #endif 179 const PetscInt *idxs; 180 PetscInt i,m,n; 181 PetscErrorCode ierr; 182 183 PetscFunctionBegin; 184 if (scall == MAT_REUSE_MATRIX) { 185 PetscBool ismatis; 186 187 ierr = PetscObjectTypeCompare((PetscObject)*newmat,MATIS,&ismatis);CHKERRQ(ierr); 188 if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Not of MATIS type"); 189 newmatis = (Mat_IS*)(*newmat)->data; 190 if (!newmatis->getsub_ris) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local row IS"); 191 if (!newmatis->getsub_cis) SETERRQ(PetscObjectComm((PetscObject)*newmat),PETSC_ERR_ARG_WRONG,"Cannot reuse matrix! Misses local col IS"); 192 } 193 /* irow and icol may not have duplicate entries */ 194 #if defined(PETSC_USE_DEBUG) 195 ierr = MatCreateVecs(mat,<est,&rtest);CHKERRQ(ierr); 196 ierr = ISGetLocalSize(irow,&n);CHKERRQ(ierr); 197 ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 198 for (i=0;i<n;i++) { 199 ierr = VecSetValue(rtest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 200 } 201 ierr = VecAssemblyBegin(rtest);CHKERRQ(ierr); 202 ierr = VecAssemblyEnd(rtest);CHKERRQ(ierr); 203 ierr = VecGetLocalSize(rtest,&n);CHKERRQ(ierr); 204 ierr = VecGetOwnershipRange(rtest,&m,NULL);CHKERRQ(ierr); 205 ierr = VecGetArrayRead(rtest,&array);CHKERRQ(ierr); 206 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])); 207 ierr = VecRestoreArrayRead(rtest,&array);CHKERRQ(ierr); 208 ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 209 ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 210 ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 211 for (i=0;i<n;i++) { 212 ierr = VecSetValue(ltest,idxs[i],1.0,ADD_VALUES);CHKERRQ(ierr); 213 } 214 ierr = VecAssemblyBegin(ltest);CHKERRQ(ierr); 215 ierr = VecAssemblyEnd(ltest);CHKERRQ(ierr); 216 ierr = VecGetLocalSize(ltest,&n);CHKERRQ(ierr); 217 ierr = VecGetOwnershipRange(ltest,&m,NULL);CHKERRQ(ierr); 218 ierr = VecGetArrayRead(ltest,&array);CHKERRQ(ierr); 219 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])); 220 ierr = VecRestoreArrayRead(ltest,&array);CHKERRQ(ierr); 221 ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 222 ierr = VecDestroy(&rtest);CHKERRQ(ierr); 223 ierr = VecDestroy(<est);CHKERRQ(ierr); 224 #endif 225 if (scall == MAT_INITIAL_MATRIX) { 226 Mat_IS *matis = (Mat_IS*)mat->data; 227 ISLocalToGlobalMapping rl2g; 228 IS is; 229 PetscInt *lidxs,*lgidxs,*newgidxs; 230 PetscInt ll,newloc; 231 MPI_Comm comm; 232 233 ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr); 234 ierr = ISGetLocalSize(irow,&m);CHKERRQ(ierr); 235 ierr = ISGetLocalSize(icol,&n);CHKERRQ(ierr); 236 ierr = MatCreate(comm,newmat);CHKERRQ(ierr); 237 ierr = MatSetType(*newmat,MATIS);CHKERRQ(ierr); 238 ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 239 /* communicate irow to their owners in the layout */ 240 ierr = ISGetIndices(irow,&idxs);CHKERRQ(ierr); 241 ierr = PetscLayoutMapLocal_Private(mat->rmap,m,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 242 ierr = ISRestoreIndices(irow,&idxs);CHKERRQ(ierr); 243 ierr = MatISSetUpSF(mat);CHKERRQ(ierr); 244 ierr = PetscMemzero(matis->sf_rootdata,matis->sf->nroots*sizeof(PetscInt));CHKERRQ(ierr); 245 for (i=0;i<ll;i++) matis->sf_rootdata[lidxs[i]] = lgidxs[i]+1; 246 ierr = PetscFree(lidxs);CHKERRQ(ierr); 247 ierr = PetscFree(lgidxs);CHKERRQ(ierr); 248 ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 249 ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 250 for (i=0,newloc=0;i<matis->sf->nleaves;i++) if (matis->sf_leafdata[i]) newloc++; 251 ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 252 ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 253 for (i=0,newloc=0;i<matis->sf->nleaves;i++) 254 if (matis->sf_leafdata[i]) { 255 lidxs[newloc] = i; 256 newgidxs[newloc++] = matis->sf_leafdata[i]-1; 257 } 258 ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 259 ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 260 ierr = ISDestroy(&is);CHKERRQ(ierr); 261 /* local is to extract local submatrix */ 262 newmatis = (Mat_IS*)(*newmat)->data; 263 ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_ris);CHKERRQ(ierr); 264 if (mat->congruentlayouts == -1) { /* first time we compare rows and cols layouts */ 265 PetscBool cong; 266 ierr = PetscLayoutCompare(mat->rmap,mat->cmap,&cong);CHKERRQ(ierr); 267 if (cong) mat->congruentlayouts = 1; 268 else mat->congruentlayouts = 0; 269 } 270 if (mat->congruentlayouts && irow == icol) { 271 ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,rl2g);CHKERRQ(ierr); 272 ierr = PetscObjectReference((PetscObject)newmatis->getsub_ris);CHKERRQ(ierr); 273 newmatis->getsub_cis = newmatis->getsub_ris; 274 } else { 275 ISLocalToGlobalMapping cl2g; 276 277 /* communicate icol to their owners in the layout */ 278 ierr = ISGetIndices(icol,&idxs);CHKERRQ(ierr); 279 ierr = PetscLayoutMapLocal_Private(mat->cmap,n,idxs,&ll,&lidxs,&lgidxs);CHKERRQ(ierr); 280 ierr = ISRestoreIndices(icol,&idxs);CHKERRQ(ierr); 281 ierr = PetscMemzero(matis->csf_rootdata,matis->csf->nroots*sizeof(PetscInt));CHKERRQ(ierr); 282 for (i=0;i<ll;i++) matis->csf_rootdata[lidxs[i]] = lgidxs[i]+1; 283 ierr = PetscFree(lidxs);CHKERRQ(ierr); 284 ierr = PetscFree(lgidxs);CHKERRQ(ierr); 285 ierr = PetscSFBcastBegin(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 286 ierr = PetscSFBcastEnd(matis->csf,MPIU_INT,matis->csf_rootdata,matis->csf_leafdata);CHKERRQ(ierr); 287 for (i=0,newloc=0;i<matis->csf->nleaves;i++) if (matis->csf_leafdata[i]) newloc++; 288 ierr = PetscMalloc1(newloc,&newgidxs);CHKERRQ(ierr); 289 ierr = PetscMalloc1(newloc,&lidxs);CHKERRQ(ierr); 290 for (i=0,newloc=0;i<matis->csf->nleaves;i++) 291 if (matis->csf_leafdata[i]) { 292 lidxs[newloc] = i; 293 newgidxs[newloc++] = matis->csf_leafdata[i]-1; 294 } 295 ierr = ISCreateGeneral(comm,newloc,newgidxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 296 ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 297 ierr = ISDestroy(&is);CHKERRQ(ierr); 298 /* local is to extract local submatrix */ 299 ierr = ISCreateGeneral(comm,newloc,lidxs,PETSC_OWN_POINTER,&newmatis->getsub_cis);CHKERRQ(ierr); 300 ierr = MatSetLocalToGlobalMapping(*newmat,rl2g,cl2g);CHKERRQ(ierr); 301 ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 302 } 303 ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 304 } else { 305 ierr = MatISGetLocalMat(*newmat,&newlocmat);CHKERRQ(ierr); 306 } 307 ierr = MatISGetLocalMat(mat,&locmat);CHKERRQ(ierr); 308 newmatis = (Mat_IS*)(*newmat)->data; 309 ierr = MatGetSubMatrix(locmat,newmatis->getsub_ris,newmatis->getsub_cis,scall,&newlocmat);CHKERRQ(ierr); 310 if (scall == MAT_INITIAL_MATRIX) { 311 ierr = MatISSetLocalMat(*newmat,newlocmat);CHKERRQ(ierr); 312 ierr = MatDestroy(&newlocmat);CHKERRQ(ierr); 313 } 314 ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 315 ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 316 PetscFunctionReturn(0); 317 } 318 319 #undef __FUNCT__ 320 #define __FUNCT__ "MatCopy_IS" 321 static PetscErrorCode MatCopy_IS(Mat A,Mat B,MatStructure str) 322 { 323 Mat_IS *a = (Mat_IS*)A->data,*b; 324 PetscBool ismatis; 325 PetscErrorCode ierr; 326 327 PetscFunctionBegin; 328 ierr = PetscObjectTypeCompare((PetscObject)B,MATIS,&ismatis);CHKERRQ(ierr); 329 if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Need to be implemented"); 330 b = (Mat_IS*)B->data; 331 ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr); 332 PetscFunctionReturn(0); 333 } 334 335 #undef __FUNCT__ 336 #define __FUNCT__ "MatMissingDiagonal_IS" 337 static PetscErrorCode MatMissingDiagonal_IS(Mat A,PetscBool *missing,PetscInt *d) 338 { 339 Vec v; 340 const PetscScalar *array; 341 PetscInt i,n; 342 PetscErrorCode ierr; 343 344 PetscFunctionBegin; 345 *missing = PETSC_FALSE; 346 ierr = MatCreateVecs(A,NULL,&v);CHKERRQ(ierr); 347 ierr = MatGetDiagonal(A,v);CHKERRQ(ierr); 348 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 349 ierr = VecGetArrayRead(v,&array);CHKERRQ(ierr); 350 for (i=0;i<n;i++) if (array[i] == 0.) break; 351 ierr = VecRestoreArrayRead(v,&array);CHKERRQ(ierr); 352 ierr = VecDestroy(&v);CHKERRQ(ierr); 353 if (i != n) *missing = PETSC_TRUE; 354 if (d) { 355 *d = -1; 356 if (*missing) { 357 PetscInt rstart; 358 ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 359 *d = i+rstart; 360 } 361 } 362 PetscFunctionReturn(0); 363 } 364 365 #undef __FUNCT__ 366 #define __FUNCT__ "MatISSetUpSF_IS" 367 static PetscErrorCode MatISSetUpSF_IS(Mat B) 368 { 369 Mat_IS *matis = (Mat_IS*)(B->data); 370 const PetscInt *gidxs; 371 PetscInt nleaves; 372 PetscErrorCode ierr; 373 374 PetscFunctionBegin; 375 if (matis->sf) PetscFunctionReturn(0); 376 ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 377 ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 378 ierr = ISLocalToGlobalMappingGetSize(B->rmap->mapping,&nleaves);CHKERRQ(ierr); 379 ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 380 ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 381 ierr = PetscMalloc2(matis->sf->nroots,&matis->sf_rootdata,matis->sf->nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 382 if (B->rmap->mapping != B->cmap->mapping) { /* setup SF for columns */ 383 ierr = ISLocalToGlobalMappingGetSize(B->cmap->mapping,&nleaves);CHKERRQ(ierr); 384 ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->csf);CHKERRQ(ierr); 385 ierr = ISLocalToGlobalMappingGetIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 386 ierr = PetscSFSetGraphLayout(matis->csf,B->cmap,nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 387 ierr = ISLocalToGlobalMappingRestoreIndices(B->cmap->mapping,&gidxs);CHKERRQ(ierr); 388 ierr = PetscMalloc2(matis->csf->nroots,&matis->csf_rootdata,matis->csf->nleaves,&matis->csf_leafdata);CHKERRQ(ierr); 389 } else { 390 matis->csf = matis->sf; 391 matis->csf_leafdata = matis->sf_leafdata; 392 matis->csf_rootdata = matis->sf_rootdata; 393 } 394 PetscFunctionReturn(0); 395 } 396 397 #undef __FUNCT__ 398 #define __FUNCT__ "MatISSetPreallocation" 399 /*@ 400 MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 401 402 Collective on MPI_Comm 403 404 Input Parameters: 405 + B - the matrix 406 . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 407 (same value is used for all local rows) 408 . d_nnz - array containing the number of nonzeros in the various rows of the 409 DIAGONAL portion of the local submatrix (possibly different for each row) 410 or NULL, if d_nz is used to specify the nonzero structure. 411 The size of this array is equal to the number of local rows, i.e 'm'. 412 For matrices that will be factored, you must leave room for (and set) 413 the diagonal entry even if it is zero. 414 . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 415 submatrix (same value is used for all local rows). 416 - o_nnz - array containing the number of nonzeros in the various rows of the 417 OFF-DIAGONAL portion of the local submatrix (possibly different for 418 each row) or NULL, if o_nz is used to specify the nonzero 419 structure. The size of this array is equal to the number 420 of local rows, i.e 'm'. 421 422 If the *_nnz parameter is given then the *_nz parameter is ignored 423 424 Level: intermediate 425 426 Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 427 from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 428 matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 429 430 .keywords: matrix 431 432 .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS 433 @*/ 434 PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 435 { 436 PetscErrorCode ierr; 437 438 PetscFunctionBegin; 439 PetscValidHeaderSpecific(B,MAT_CLASSID,1); 440 PetscValidType(B,1); 441 ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 442 PetscFunctionReturn(0); 443 } 444 445 #undef __FUNCT__ 446 #define __FUNCT__ "MatISSetPreallocation_IS" 447 static PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 448 { 449 Mat_IS *matis = (Mat_IS*)(B->data); 450 PetscInt bs,i,nlocalcols; 451 PetscErrorCode ierr; 452 453 PetscFunctionBegin; 454 if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 455 ierr = MatISSetUpSF(B);CHKERRQ(ierr); 456 457 if (!d_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nz; 458 else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 459 460 if (!o_nnz) for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nz; 461 else for (i=0;i<matis->sf->nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 462 463 ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 464 ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 465 ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 466 ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 467 468 for (i=0;i<matis->sf->nleaves;i++) matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 469 ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 470 471 for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 472 ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 473 474 for (i=0;i<matis->sf->nleaves/bs;i++) matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 475 ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 476 PetscFunctionReturn(0); 477 } 478 479 #undef __FUNCT__ 480 #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private" 481 PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 482 { 483 Mat_IS *matis = (Mat_IS*)(A->data); 484 PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 485 const PetscInt *global_indices_r,*global_indices_c; 486 PetscInt i,j,bs,rows,cols; 487 PetscInt lrows,lcols; 488 PetscInt local_rows,local_cols; 489 PetscMPIInt nsubdomains; 490 PetscBool isdense,issbaij; 491 PetscErrorCode ierr; 492 493 PetscFunctionBegin; 494 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr); 495 ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr); 496 ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 497 ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 498 ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 499 ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 500 ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 501 if (A->rmap->mapping != A->cmap->mapping) { 502 ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 503 } else { 504 global_indices_c = global_indices_r; 505 } 506 507 if (issbaij) { 508 ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 509 } 510 /* 511 An SF reduce is needed to sum up properly on shared rows. 512 Note that generally preallocation is not exact, since it overestimates nonzeros 513 */ 514 ierr = MatISSetUpSF(A);CHKERRQ(ierr); 515 ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr); 516 ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr); 517 /* All processes need to compute entire row ownership */ 518 ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 519 ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 520 for (i=0;i<nsubdomains;i++) { 521 for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 522 row_ownership[j] = i; 523 } 524 } 525 ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 526 527 /* 528 my_dnz and my_onz contains exact contribution to preallocation from each local mat 529 then, they will be summed up properly. This way, preallocation is always sufficient 530 */ 531 ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 532 /* preallocation as a MATAIJ */ 533 if (isdense) { /* special case for dense local matrices */ 534 for (i=0;i<local_rows;i++) { 535 PetscInt index_row = global_indices_r[i]; 536 for (j=i;j<local_rows;j++) { 537 PetscInt owner = row_ownership[index_row]; 538 PetscInt index_col = global_indices_c[j]; 539 if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 540 my_dnz[i] += 1; 541 } else { /* offdiag block */ 542 my_onz[i] += 1; 543 } 544 /* same as before, interchanging rows and cols */ 545 if (i != j) { 546 owner = row_ownership[index_col]; 547 if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 548 my_dnz[j] += 1; 549 } else { 550 my_onz[j] += 1; 551 } 552 } 553 } 554 } 555 } else if (matis->A->ops->getrowij) { 556 const PetscInt *ii,*jj,*jptr; 557 PetscBool done; 558 ierr = MatGetRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr); 559 if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__); 560 jptr = jj; 561 for (i=0;i<local_rows;i++) { 562 PetscInt index_row = global_indices_r[i]; 563 for (j=0;j<ii[i+1]-ii[i];j++,jptr++) { 564 PetscInt owner = row_ownership[index_row]; 565 PetscInt index_col = global_indices_c[*jptr]; 566 if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 567 my_dnz[i] += 1; 568 } else { /* offdiag block */ 569 my_onz[i] += 1; 570 } 571 /* same as before, interchanging rows and cols */ 572 if (issbaij && index_col != index_row) { 573 owner = row_ownership[index_col]; 574 if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 575 my_dnz[*jptr] += 1; 576 } else { 577 my_onz[*jptr] += 1; 578 } 579 } 580 } 581 } 582 ierr = MatRestoreRowIJ(matis->A,0,PETSC_FALSE,PETSC_FALSE,&local_rows,&ii,&jj,&done);CHKERRQ(ierr); 583 if (!done) SETERRQ1(PetscObjectComm((PetscObject)(matis->A)),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__); 584 } else { /* loop over rows and use MatGetRow */ 585 for (i=0;i<local_rows;i++) { 586 const PetscInt *cols; 587 PetscInt ncols,index_row = global_indices_r[i]; 588 ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 589 for (j=0;j<ncols;j++) { 590 PetscInt owner = row_ownership[index_row]; 591 PetscInt index_col = global_indices_c[cols[j]]; 592 if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 593 my_dnz[i] += 1; 594 } else { /* offdiag block */ 595 my_onz[i] += 1; 596 } 597 /* same as before, interchanging rows and cols */ 598 if (issbaij && index_col != index_row) { 599 owner = row_ownership[index_col]; 600 if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 601 my_dnz[cols[j]] += 1; 602 } else { 603 my_onz[cols[j]] += 1; 604 } 605 } 606 } 607 ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 608 } 609 } 610 ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 611 if (global_indices_c != global_indices_r) { 612 ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 613 } 614 ierr = PetscFree(row_ownership);CHKERRQ(ierr); 615 616 /* Reduce my_dnz and my_onz */ 617 if (maxreduce) { 618 ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 619 ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 620 ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 621 ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 622 } else { 623 ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 624 ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 625 ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 626 ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 627 } 628 ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 629 630 /* Resize preallocation if overestimated */ 631 for (i=0;i<lrows;i++) { 632 dnz[i] = PetscMin(dnz[i],lcols); 633 onz[i] = PetscMin(onz[i],cols-lcols); 634 } 635 /* set preallocation */ 636 ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 637 for (i=0;i<lrows/bs;i++) { 638 dnz[i] = dnz[i*bs]/bs; 639 onz[i] = onz[i*bs]/bs; 640 } 641 ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 642 ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 643 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 644 if (issbaij) { 645 ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 646 } 647 PetscFunctionReturn(0); 648 } 649 650 #undef __FUNCT__ 651 #define __FUNCT__ "MatISGetMPIXAIJ_IS" 652 static PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 653 { 654 Mat_IS *matis = (Mat_IS*)(mat->data); 655 Mat local_mat; 656 /* info on mat */ 657 PetscInt bs,rows,cols,lrows,lcols; 658 PetscInt local_rows,local_cols; 659 PetscBool isseqdense,isseqsbaij,isseqaij,isseqbaij; 660 #if defined (PETSC_USE_DEBUG) 661 PetscBool lb[4],bb[4]; 662 #endif 663 PetscMPIInt nsubdomains; 664 /* values insertion */ 665 PetscScalar *array; 666 /* work */ 667 PetscErrorCode ierr; 668 669 PetscFunctionBegin; 670 /* get info from mat */ 671 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 672 if (nsubdomains == 1) { 673 if (reuse == MAT_INITIAL_MATRIX) { 674 ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 675 } else { 676 ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 677 } 678 PetscFunctionReturn(0); 679 } 680 ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 681 ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 682 ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 683 ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 684 ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isseqdense);CHKERRQ(ierr); 685 ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 686 ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQBAIJ,&isseqbaij);CHKERRQ(ierr); 687 ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr); 688 if (!isseqdense && !isseqaij && !isseqbaij && !isseqsbaij) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not for matrix type %s",((PetscObject)(matis->A))->type_name); 689 #if defined (PETSC_USE_DEBUG) 690 lb[0] = isseqdense; 691 lb[1] = isseqaij; 692 lb[2] = isseqbaij; 693 lb[3] = isseqsbaij; 694 ierr = MPIU_Allreduce(lb,bb,4,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 695 if (!bb[0] && !bb[1] && !bb[2] && !bb[3]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Local matrices must have the same type"); 696 #endif 697 698 if (reuse == MAT_INITIAL_MATRIX) { 699 ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 700 ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 701 ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 702 if (!isseqsbaij) { 703 ierr = MatSetType(*M,MATAIJ);CHKERRQ(ierr); 704 } else { 705 ierr = MatSetType(*M,MATSBAIJ);CHKERRQ(ierr); 706 } 707 ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 708 } else { 709 PetscInt mbs,mrows,mcols,mlrows,mlcols; 710 /* some checks */ 711 ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 712 ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 713 ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 714 if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 715 if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 716 if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 717 if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 718 if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 719 ierr = MatZeroEntries(*M);CHKERRQ(ierr); 720 } 721 722 if (isseqsbaij) { 723 ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 724 } else { 725 ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 726 local_mat = matis->A; 727 } 728 729 /* Set values */ 730 ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 731 if (isseqdense) { /* special case for dense local matrices */ 732 PetscInt i,*dummy_rows,*dummy_cols; 733 734 if (local_rows != local_cols) { 735 ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 736 for (i=0;i<local_rows;i++) dummy_rows[i] = i; 737 for (i=0;i<local_cols;i++) dummy_cols[i] = i; 738 } else { 739 ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 740 for (i=0;i<local_rows;i++) dummy_rows[i] = i; 741 dummy_cols = dummy_rows; 742 } 743 ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 744 ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 745 ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 746 ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 747 if (dummy_rows != dummy_cols) { 748 ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 749 } else { 750 ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 751 } 752 } else if (isseqaij) { 753 PetscInt i,nvtxs,*xadj,*adjncy; 754 PetscBool done; 755 756 ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 757 if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatGetRowIJ called from %s\n",__FUNCT__); 758 ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 759 for (i=0;i<nvtxs;i++) { 760 ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 761 } 762 ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 763 if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called from %s\n",__FUNCT__); 764 ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 765 } else { /* very basic values insertion for all other matrix types */ 766 PetscInt i; 767 768 for (i=0;i<local_rows;i++) { 769 PetscInt j; 770 const PetscInt *local_indices_cols; 771 772 ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 773 ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 774 ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 775 } 776 } 777 ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 778 ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 779 ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 780 if (isseqdense) { 781 ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 782 } 783 PetscFunctionReturn(0); 784 } 785 786 #undef __FUNCT__ 787 #define __FUNCT__ "MatISGetMPIXAIJ" 788 /*@ 789 MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 790 791 Input Parameter: 792 . mat - the matrix (should be of type MATIS) 793 . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 794 795 Output Parameter: 796 . newmat - the matrix in AIJ format 797 798 Level: developer 799 800 Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 801 802 .seealso: MATIS 803 @*/ 804 PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 805 { 806 PetscErrorCode ierr; 807 808 PetscFunctionBegin; 809 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 810 PetscValidLogicalCollectiveEnum(mat,reuse,2); 811 PetscValidPointer(newmat,3); 812 if (reuse != MAT_INITIAL_MATRIX) { 813 PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 814 PetscCheckSameComm(mat,1,*newmat,3); 815 if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 816 } 817 ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 818 PetscFunctionReturn(0); 819 } 820 821 #undef __FUNCT__ 822 #define __FUNCT__ "MatDuplicate_IS" 823 PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 824 { 825 PetscErrorCode ierr; 826 Mat_IS *matis = (Mat_IS*)(mat->data); 827 PetscInt bs,m,n,M,N; 828 Mat B,localmat; 829 830 PetscFunctionBegin; 831 ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 832 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 833 ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 834 ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 835 ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 836 ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 837 ierr = MatDestroy(&localmat);CHKERRQ(ierr); 838 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 839 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 840 *newmat = B; 841 PetscFunctionReturn(0); 842 } 843 844 #undef __FUNCT__ 845 #define __FUNCT__ "MatIsHermitian_IS" 846 static PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 847 { 848 PetscErrorCode ierr; 849 Mat_IS *matis = (Mat_IS*)A->data; 850 PetscBool local_sym; 851 852 PetscFunctionBegin; 853 ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 854 ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 855 PetscFunctionReturn(0); 856 } 857 858 #undef __FUNCT__ 859 #define __FUNCT__ "MatIsSymmetric_IS" 860 static PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 861 { 862 PetscErrorCode ierr; 863 Mat_IS *matis = (Mat_IS*)A->data; 864 PetscBool local_sym; 865 866 PetscFunctionBegin; 867 ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 868 ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 869 PetscFunctionReturn(0); 870 } 871 872 #undef __FUNCT__ 873 #define __FUNCT__ "MatIsStructurallySymmetric_IS" 874 static PetscErrorCode MatIsStructurallySymmetric_IS(Mat A,PetscBool *flg) 875 { 876 PetscErrorCode ierr; 877 Mat_IS *matis = (Mat_IS*)A->data; 878 PetscBool local_sym; 879 880 PetscFunctionBegin; 881 if (A->rmap->mapping != A->cmap->mapping) { 882 *flg = PETSC_FALSE; 883 PetscFunctionReturn(0); 884 } 885 ierr = MatIsStructurallySymmetric(matis->A,&local_sym);CHKERRQ(ierr); 886 ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 887 PetscFunctionReturn(0); 888 } 889 890 #undef __FUNCT__ 891 #define __FUNCT__ "MatDestroy_IS" 892 static PetscErrorCode MatDestroy_IS(Mat A) 893 { 894 PetscErrorCode ierr; 895 Mat_IS *b = (Mat_IS*)A->data; 896 897 PetscFunctionBegin; 898 ierr = MatDestroy(&b->A);CHKERRQ(ierr); 899 ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 900 ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 901 ierr = VecDestroy(&b->x);CHKERRQ(ierr); 902 ierr = VecDestroy(&b->y);CHKERRQ(ierr); 903 ierr = VecDestroy(&b->counter);CHKERRQ(ierr); 904 ierr = ISDestroy(&b->getsub_ris);CHKERRQ(ierr); 905 ierr = ISDestroy(&b->getsub_cis);CHKERRQ(ierr); 906 if (b->sf != b->csf) { 907 ierr = PetscSFDestroy(&b->csf);CHKERRQ(ierr); 908 ierr = PetscFree2(b->csf_rootdata,b->csf_leafdata);CHKERRQ(ierr); 909 } 910 ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 911 ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 912 ierr = PetscFree(A->data);CHKERRQ(ierr); 913 ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 914 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 915 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 916 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 917 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 918 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",NULL);CHKERRQ(ierr); 919 PetscFunctionReturn(0); 920 } 921 922 #undef __FUNCT__ 923 #define __FUNCT__ "MatMult_IS" 924 static PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 925 { 926 PetscErrorCode ierr; 927 Mat_IS *is = (Mat_IS*)A->data; 928 PetscScalar zero = 0.0; 929 930 PetscFunctionBegin; 931 /* scatter the global vector x into the local work vector */ 932 ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 933 ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 934 935 /* multiply the local matrix */ 936 ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 937 938 /* scatter product back into global memory */ 939 ierr = VecSet(y,zero);CHKERRQ(ierr); 940 ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 941 ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 942 PetscFunctionReturn(0); 943 } 944 945 #undef __FUNCT__ 946 #define __FUNCT__ "MatMultAdd_IS" 947 static PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 948 { 949 Vec temp_vec; 950 PetscErrorCode ierr; 951 952 PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 953 if (v3 != v2) { 954 ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 955 ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 956 } else { 957 ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 958 ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 959 ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 960 ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 961 ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 962 } 963 PetscFunctionReturn(0); 964 } 965 966 #undef __FUNCT__ 967 #define __FUNCT__ "MatMultTranspose_IS" 968 static PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 969 { 970 Mat_IS *is = (Mat_IS*)A->data; 971 PetscErrorCode ierr; 972 973 PetscFunctionBegin; 974 /* scatter the global vector x into the local work vector */ 975 ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 976 ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 977 978 /* multiply the local matrix */ 979 ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 980 981 /* scatter product back into global vector */ 982 ierr = VecSet(x,0);CHKERRQ(ierr); 983 ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 984 ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 985 PetscFunctionReturn(0); 986 } 987 988 #undef __FUNCT__ 989 #define __FUNCT__ "MatMultTransposeAdd_IS" 990 static PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 991 { 992 Vec temp_vec; 993 PetscErrorCode ierr; 994 995 PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 996 if (v3 != v2) { 997 ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 998 ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 999 } else { 1000 ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 1001 ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 1002 ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 1003 ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 1004 ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 1005 } 1006 PetscFunctionReturn(0); 1007 } 1008 1009 #undef __FUNCT__ 1010 #define __FUNCT__ "MatView_IS" 1011 static PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 1012 { 1013 Mat_IS *a = (Mat_IS*)A->data; 1014 PetscErrorCode ierr; 1015 PetscViewer sviewer; 1016 PetscBool isascii,view = PETSC_TRUE; 1017 1018 PetscFunctionBegin; 1019 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 1020 if (isascii) { 1021 PetscViewerFormat format; 1022 1023 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1024 if (format == PETSC_VIEWER_ASCII_INFO) view = PETSC_FALSE; 1025 } 1026 if (!view) PetscFunctionReturn(0); 1027 ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 1028 ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 1029 ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 1030 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1031 PetscFunctionReturn(0); 1032 } 1033 1034 #undef __FUNCT__ 1035 #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 1036 static PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 1037 { 1038 PetscErrorCode ierr; 1039 PetscInt nr,rbs,nc,cbs; 1040 Mat_IS *is = (Mat_IS*)A->data; 1041 IS from,to; 1042 Vec cglobal,rglobal; 1043 1044 PetscFunctionBegin; 1045 PetscCheckSameComm(A,1,rmapping,2); 1046 PetscCheckSameComm(A,1,cmapping,3); 1047 /* Destroy any previous data */ 1048 ierr = VecDestroy(&is->x);CHKERRQ(ierr); 1049 ierr = VecDestroy(&is->y);CHKERRQ(ierr); 1050 ierr = VecDestroy(&is->counter);CHKERRQ(ierr); 1051 ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 1052 ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 1053 ierr = MatDestroy(&is->A);CHKERRQ(ierr); 1054 ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 1055 ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 1056 1057 /* Setup Layout and set local to global maps */ 1058 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1059 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1060 ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 1061 ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 1062 1063 /* Create the local matrix A */ 1064 ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 1065 ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 1066 ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 1067 ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 1068 ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 1069 ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 1070 ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 1071 ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 1072 ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 1073 ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 1074 ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 1075 1076 if (!is->islocalref) { /* setup scatters and local vectors for MatMult */ 1077 /* Create the local work vectors */ 1078 ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 1079 1080 /* setup the global to local scatters */ 1081 ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 1082 ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 1083 ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 1084 ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 1085 if (rmapping != cmapping) { 1086 ierr = ISDestroy(&to);CHKERRQ(ierr); 1087 ierr = ISDestroy(&from);CHKERRQ(ierr); 1088 ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 1089 ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 1090 ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 1091 } else { 1092 ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 1093 is->cctx = is->rctx; 1094 } 1095 1096 /* interface counter vector (local) */ 1097 ierr = VecDuplicate(is->y,&is->counter);CHKERRQ(ierr); 1098 ierr = VecSet(is->y,1.);CHKERRQ(ierr); 1099 ierr = VecScatterBegin(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1100 ierr = VecScatterEnd(is->rctx,is->y,rglobal,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1101 ierr = VecScatterBegin(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1102 ierr = VecScatterEnd(is->rctx,rglobal,is->counter,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1103 1104 /* free workspace */ 1105 ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 1106 ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 1107 ierr = ISDestroy(&to);CHKERRQ(ierr); 1108 ierr = ISDestroy(&from);CHKERRQ(ierr); 1109 } 1110 ierr = MatSetUp(A);CHKERRQ(ierr); 1111 PetscFunctionReturn(0); 1112 } 1113 1114 #undef __FUNCT__ 1115 #define __FUNCT__ "MatSetValues_IS" 1116 static PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 1117 { 1118 Mat_IS *is = (Mat_IS*)mat->data; 1119 PetscErrorCode ierr; 1120 #if defined(PETSC_USE_DEBUG) 1121 PetscInt i,zm,zn; 1122 #endif 1123 PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 1124 1125 PetscFunctionBegin; 1126 #if defined(PETSC_USE_DEBUG) 1127 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); 1128 /* count negative indices */ 1129 for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 1130 for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 1131 #endif 1132 ierr = ISGlobalToLocalMappingApply(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 1133 ierr = ISGlobalToLocalMappingApply(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 1134 #if defined(PETSC_USE_DEBUG) 1135 /* count negative indices (should be the same as before) */ 1136 for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 1137 for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 1138 if (zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the row indices can not be mapped! Maybe you should not use MATIS"); 1139 if (zn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the column indices can not be mapped! Maybe you should not use MATIS"); 1140 #endif 1141 ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 1142 PetscFunctionReturn(0); 1143 } 1144 1145 #undef __FUNCT__ 1146 #define __FUNCT__ "MatSetValuesBlocked_IS" 1147 static PetscErrorCode MatSetValuesBlocked_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 1148 { 1149 Mat_IS *is = (Mat_IS*)mat->data; 1150 PetscErrorCode ierr; 1151 #if defined(PETSC_USE_DEBUG) 1152 PetscInt i,zm,zn; 1153 #endif 1154 PetscInt rows_l[MATIS_MAX_ENTRIES_INSERTION],cols_l[MATIS_MAX_ENTRIES_INSERTION]; 1155 1156 PetscFunctionBegin; 1157 #if defined(PETSC_USE_DEBUG) 1158 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); 1159 /* count negative indices */ 1160 for (i=0,zm=0;i<m;i++) if (rows[i] < 0) zm++; 1161 for (i=0,zn=0;i<n;i++) if (cols[i] < 0) zn++; 1162 #endif 1163 ierr = ISGlobalToLocalMappingApplyBlock(mat->rmap->mapping,IS_GTOLM_MASK,m,rows,&m,rows_l);CHKERRQ(ierr); 1164 ierr = ISGlobalToLocalMappingApplyBlock(mat->cmap->mapping,IS_GTOLM_MASK,n,cols,&n,cols_l);CHKERRQ(ierr); 1165 #if defined(PETSC_USE_DEBUG) 1166 /* count negative indices (should be the same as before) */ 1167 for (i=0;i<m;i++) if (rows_l[i] < 0) zm--; 1168 for (i=0;i<n;i++) if (cols_l[i] < 0) zn--; 1169 if (zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the row indices can not be mapped! Maybe you should not use MATIS"); 1170 if (zn) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Some of the column indices can not be mapped! Maybe you should not use MATIS"); 1171 #endif 1172 ierr = MatSetValuesBlocked(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 1173 PetscFunctionReturn(0); 1174 } 1175 1176 #undef __FUNCT__ 1177 #define __FUNCT__ "MatSetValuesLocal_IS" 1178 static PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1179 { 1180 PetscErrorCode ierr; 1181 Mat_IS *is = (Mat_IS*)A->data; 1182 1183 PetscFunctionBegin; 1184 ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1185 PetscFunctionReturn(0); 1186 } 1187 1188 #undef __FUNCT__ 1189 #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 1190 static PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 1191 { 1192 PetscErrorCode ierr; 1193 Mat_IS *is = (Mat_IS*)A->data; 1194 1195 PetscFunctionBegin; 1196 ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 1197 PetscFunctionReturn(0); 1198 } 1199 1200 #undef __FUNCT__ 1201 #define __FUNCT__ "MatISZeroRowsColumnsLocal_Private" 1202 static PetscErrorCode MatISZeroRowsColumnsLocal_Private(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,PetscBool columns) 1203 { 1204 Mat_IS *is = (Mat_IS*)A->data; 1205 PetscErrorCode ierr; 1206 1207 PetscFunctionBegin; 1208 if (!n) { 1209 is->pure_neumann = PETSC_TRUE; 1210 } else { 1211 PetscInt i; 1212 is->pure_neumann = PETSC_FALSE; 1213 1214 if (columns) { 1215 ierr = MatZeroRowsColumns(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 1216 } else { 1217 ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 1218 } 1219 if (diag != 0.) { 1220 const PetscScalar *array; 1221 ierr = VecGetArrayRead(is->counter,&array);CHKERRQ(ierr); 1222 for (i=0; i<n; i++) { 1223 ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 1224 } 1225 ierr = VecRestoreArrayRead(is->counter,&array);CHKERRQ(ierr); 1226 } 1227 ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1228 ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1229 } 1230 PetscFunctionReturn(0); 1231 } 1232 1233 #undef __FUNCT__ 1234 #define __FUNCT__ "MatZeroRowsColumns_Private_IS" 1235 static PetscErrorCode MatZeroRowsColumns_Private_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b,PetscBool columns) 1236 { 1237 Mat_IS *matis = (Mat_IS*)A->data; 1238 PetscInt nr,nl,len,i; 1239 PetscInt *lrows; 1240 PetscErrorCode ierr; 1241 1242 PetscFunctionBegin; 1243 #if defined(PETSC_USE_DEBUG) 1244 if (columns || diag != 0. || (x && b)) { 1245 PetscBool cong; 1246 ierr = PetscLayoutCompare(A->rmap,A->cmap,&cong);CHKERRQ(ierr); 1247 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"); 1248 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"); 1249 if (!cong && x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"A->rmap and A->cmap need to be congruent"); 1250 } 1251 #endif 1252 /* get locally owned rows */ 1253 ierr = PetscLayoutMapLocal_Private(A->rmap,n,rows,&len,&lrows,NULL);CHKERRQ(ierr); 1254 /* fix right hand side if needed */ 1255 if (x && b) { 1256 const PetscScalar *xx; 1257 PetscScalar *bb; 1258 1259 ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr); 1260 ierr = VecGetArray(b, &bb);CHKERRQ(ierr); 1261 for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]]; 1262 ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr); 1263 ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr); 1264 } 1265 /* get rows associated to the local matrices */ 1266 ierr = MatISSetUpSF(A);CHKERRQ(ierr); 1267 ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr); 1268 ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr); 1269 ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 1270 for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1; 1271 ierr = PetscFree(lrows);CHKERRQ(ierr); 1272 ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 1273 ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 1274 ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr); 1275 for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i; 1276 ierr = MatISZeroRowsColumnsLocal_Private(A,nr,lrows,diag,columns);CHKERRQ(ierr); 1277 ierr = PetscFree(lrows);CHKERRQ(ierr); 1278 PetscFunctionReturn(0); 1279 } 1280 1281 #undef __FUNCT__ 1282 #define __FUNCT__ "MatZeroRows_IS" 1283 static PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 1284 { 1285 PetscErrorCode ierr; 1286 1287 PetscFunctionBegin; 1288 ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_FALSE);CHKERRQ(ierr); 1289 PetscFunctionReturn(0); 1290 } 1291 1292 #undef __FUNCT__ 1293 #define __FUNCT__ "MatZeroRowsColumns_IS" 1294 static PetscErrorCode MatZeroRowsColumns_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 1295 { 1296 PetscErrorCode ierr; 1297 1298 PetscFunctionBegin; 1299 ierr = MatZeroRowsColumns_Private_IS(A,n,rows,diag,x,b,PETSC_TRUE);CHKERRQ(ierr); 1300 PetscFunctionReturn(0); 1301 } 1302 1303 #undef __FUNCT__ 1304 #define __FUNCT__ "MatAssemblyBegin_IS" 1305 static PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 1306 { 1307 Mat_IS *is = (Mat_IS*)A->data; 1308 PetscErrorCode ierr; 1309 1310 PetscFunctionBegin; 1311 ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 1312 PetscFunctionReturn(0); 1313 } 1314 1315 #undef __FUNCT__ 1316 #define __FUNCT__ "MatAssemblyEnd_IS" 1317 static PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 1318 { 1319 Mat_IS *is = (Mat_IS*)A->data; 1320 PetscErrorCode ierr; 1321 1322 PetscFunctionBegin; 1323 ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 1324 PetscFunctionReturn(0); 1325 } 1326 1327 #undef __FUNCT__ 1328 #define __FUNCT__ "MatISGetLocalMat_IS" 1329 static PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 1330 { 1331 Mat_IS *is = (Mat_IS*)mat->data; 1332 1333 PetscFunctionBegin; 1334 *local = is->A; 1335 PetscFunctionReturn(0); 1336 } 1337 1338 #undef __FUNCT__ 1339 #define __FUNCT__ "MatISGetLocalMat" 1340 /*@ 1341 MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 1342 1343 Input Parameter: 1344 . mat - the matrix 1345 1346 Output Parameter: 1347 . local - the local matrix 1348 1349 Level: advanced 1350 1351 Notes: 1352 This can be called if you have precomputed the nonzero structure of the 1353 matrix and want to provide it to the inner matrix object to improve the performance 1354 of the MatSetValues() operation. 1355 1356 .seealso: MATIS 1357 @*/ 1358 PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 1359 { 1360 PetscErrorCode ierr; 1361 1362 PetscFunctionBegin; 1363 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1364 PetscValidPointer(local,2); 1365 ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 1366 PetscFunctionReturn(0); 1367 } 1368 1369 #undef __FUNCT__ 1370 #define __FUNCT__ "MatISSetLocalMat_IS" 1371 static PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 1372 { 1373 Mat_IS *is = (Mat_IS*)mat->data; 1374 PetscInt nrows,ncols,orows,ocols; 1375 PetscErrorCode ierr; 1376 1377 PetscFunctionBegin; 1378 if (is->A) { 1379 ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 1380 ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 1381 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); 1382 } 1383 ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 1384 ierr = MatDestroy(&is->A);CHKERRQ(ierr); 1385 is->A = local; 1386 PetscFunctionReturn(0); 1387 } 1388 1389 #undef __FUNCT__ 1390 #define __FUNCT__ "MatISSetLocalMat" 1391 /*@ 1392 MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 1393 1394 Input Parameter: 1395 . mat - the matrix 1396 . local - the local matrix 1397 1398 Output Parameter: 1399 1400 Level: advanced 1401 1402 Notes: 1403 This can be called if you have precomputed the local matrix and 1404 want to provide it to the matrix object MATIS. 1405 1406 .seealso: MATIS 1407 @*/ 1408 PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 1409 { 1410 PetscErrorCode ierr; 1411 1412 PetscFunctionBegin; 1413 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 1414 PetscValidHeaderSpecific(local,MAT_CLASSID,2); 1415 ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 1416 PetscFunctionReturn(0); 1417 } 1418 1419 #undef __FUNCT__ 1420 #define __FUNCT__ "MatZeroEntries_IS" 1421 static PetscErrorCode MatZeroEntries_IS(Mat A) 1422 { 1423 Mat_IS *a = (Mat_IS*)A->data; 1424 PetscErrorCode ierr; 1425 1426 PetscFunctionBegin; 1427 ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 1428 PetscFunctionReturn(0); 1429 } 1430 1431 #undef __FUNCT__ 1432 #define __FUNCT__ "MatScale_IS" 1433 static PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 1434 { 1435 Mat_IS *is = (Mat_IS*)A->data; 1436 PetscErrorCode ierr; 1437 1438 PetscFunctionBegin; 1439 ierr = MatScale(is->A,a);CHKERRQ(ierr); 1440 PetscFunctionReturn(0); 1441 } 1442 1443 #undef __FUNCT__ 1444 #define __FUNCT__ "MatGetDiagonal_IS" 1445 static PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 1446 { 1447 Mat_IS *is = (Mat_IS*)A->data; 1448 PetscErrorCode ierr; 1449 1450 PetscFunctionBegin; 1451 /* get diagonal of the local matrix */ 1452 ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 1453 1454 /* scatter diagonal back into global vector */ 1455 ierr = VecSet(v,0);CHKERRQ(ierr); 1456 ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1457 ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1458 PetscFunctionReturn(0); 1459 } 1460 1461 #undef __FUNCT__ 1462 #define __FUNCT__ "MatSetOption_IS" 1463 static PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 1464 { 1465 Mat_IS *a = (Mat_IS*)A->data; 1466 PetscErrorCode ierr; 1467 1468 PetscFunctionBegin; 1469 ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 1470 PetscFunctionReturn(0); 1471 } 1472 1473 #undef __FUNCT__ 1474 #define __FUNCT__ "MatAXPY_IS" 1475 static PetscErrorCode MatAXPY_IS(Mat Y,PetscScalar a,Mat X,MatStructure str) 1476 { 1477 Mat_IS *y = (Mat_IS*)Y->data; 1478 Mat_IS *x; 1479 #if defined(PETSC_USE_DEBUG) 1480 PetscBool ismatis; 1481 #endif 1482 PetscErrorCode ierr; 1483 1484 PetscFunctionBegin; 1485 #if defined(PETSC_USE_DEBUG) 1486 ierr = PetscObjectTypeCompare((PetscObject)X,MATIS,&ismatis);CHKERRQ(ierr); 1487 if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_SUP,"Cannot call MatAXPY(Y,a,X,str) with X not of type MATIS"); 1488 #endif 1489 x = (Mat_IS*)X->data; 1490 ierr = MatAXPY(y->A,a,x->A,str);CHKERRQ(ierr); 1491 PetscFunctionReturn(0); 1492 } 1493 1494 #undef __FUNCT__ 1495 #define __FUNCT__ "MatGetLocalSubMatrix_IS" 1496 static PetscErrorCode MatGetLocalSubMatrix_IS(Mat A,IS row,IS col,Mat *submat) 1497 { 1498 Mat lA; 1499 Mat_IS *matis; 1500 ISLocalToGlobalMapping rl2g,cl2g; 1501 IS is; 1502 const PetscInt *rg,*rl; 1503 PetscInt nrg; 1504 PetscInt N,M,nrl,i,*idxs; 1505 PetscErrorCode ierr; 1506 1507 PetscFunctionBegin; 1508 ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 1509 ierr = ISGetLocalSize(row,&nrl);CHKERRQ(ierr); 1510 ierr = ISGetIndices(row,&rl);CHKERRQ(ierr); 1511 ierr = ISLocalToGlobalMappingGetSize(A->rmap->mapping,&nrg);CHKERRQ(ierr); 1512 #if defined(PETSC_USE_DEBUG) 1513 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); 1514 #endif 1515 ierr = PetscMalloc1(nrg,&idxs);CHKERRQ(ierr); 1516 /* map from [0,nrl) to row */ 1517 for (i=0;i<nrl;i++) idxs[i] = rl[i]; 1518 #if defined(PETSC_USE_DEBUG) 1519 for (i=nrl;i<nrg;i++) idxs[i] = nrg; 1520 #else 1521 for (i=nrl;i<nrg;i++) idxs[i] = -1; 1522 #endif 1523 ierr = ISRestoreIndices(row,&rl);CHKERRQ(ierr); 1524 ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&rg);CHKERRQ(ierr); 1525 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),nrg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 1526 ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 1527 ierr = ISDestroy(&is);CHKERRQ(ierr); 1528 /* compute new l2g map for columns */ 1529 if (col != row || A->rmap->mapping != A->cmap->mapping) { 1530 const PetscInt *cg,*cl; 1531 PetscInt ncg; 1532 PetscInt ncl; 1533 1534 ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 1535 ierr = ISGetLocalSize(col,&ncl);CHKERRQ(ierr); 1536 ierr = ISGetIndices(col,&cl);CHKERRQ(ierr); 1537 ierr = ISLocalToGlobalMappingGetSize(A->cmap->mapping,&ncg);CHKERRQ(ierr); 1538 #if defined(PETSC_USE_DEBUG) 1539 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); 1540 #endif 1541 ierr = PetscMalloc1(ncg,&idxs);CHKERRQ(ierr); 1542 /* map from [0,ncl) to col */ 1543 for (i=0;i<ncl;i++) idxs[i] = cl[i]; 1544 #if defined(PETSC_USE_DEBUG) 1545 for (i=ncl;i<ncg;i++) idxs[i] = ncg; 1546 #else 1547 for (i=ncl;i<ncg;i++) idxs[i] = -1; 1548 #endif 1549 ierr = ISRestoreIndices(col,&cl);CHKERRQ(ierr); 1550 ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&cg);CHKERRQ(ierr); 1551 ierr = ISCreateGeneral(PetscObjectComm((PetscObject)A),ncg,idxs,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 1552 ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 1553 ierr = ISDestroy(&is);CHKERRQ(ierr); 1554 } else { 1555 ierr = PetscObjectReference((PetscObject)rl2g);CHKERRQ(ierr); 1556 cl2g = rl2g; 1557 } 1558 /* create the MATIS submatrix */ 1559 ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 1560 ierr = MatCreate(PetscObjectComm((PetscObject)A),submat);CHKERRQ(ierr); 1561 ierr = MatSetSizes(*submat,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 1562 ierr = MatSetType(*submat,MATIS);CHKERRQ(ierr); 1563 matis = (Mat_IS*)((*submat)->data); 1564 matis->islocalref = PETSC_TRUE; 1565 ierr = MatSetLocalToGlobalMapping(*submat,rl2g,cl2g);CHKERRQ(ierr); 1566 ierr = MatISGetLocalMat(A,&lA);CHKERRQ(ierr); 1567 ierr = MatISSetLocalMat(*submat,lA);CHKERRQ(ierr); 1568 ierr = MatSetUp(*submat);CHKERRQ(ierr); 1569 ierr = MatAssemblyBegin(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1570 ierr = MatAssemblyEnd(*submat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1571 ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 1572 ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 1573 /* remove unsupported ops */ 1574 ierr = PetscMemzero((*submat)->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1575 (*submat)->ops->destroy = MatDestroy_IS; 1576 (*submat)->ops->setvalueslocal = MatSetValuesLocal_SubMat_IS; 1577 (*submat)->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_SubMat_IS; 1578 (*submat)->ops->assemblybegin = MatAssemblyBegin_IS; 1579 (*submat)->ops->assemblyend = MatAssemblyEnd_IS; 1580 PetscFunctionReturn(0); 1581 } 1582 1583 #undef __FUNCT__ 1584 #define __FUNCT__ "MatCreateIS" 1585 /*@ 1586 MatCreateIS - Creates a "process" unassembled matrix, assembled on each 1587 process but not across processes. 1588 1589 Input Parameters: 1590 + comm - MPI communicator that will share the matrix 1591 . bs - block size of the matrix 1592 . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 1593 . rmap - local to global map for rows 1594 - cmap - local to global map for cols 1595 1596 Output Parameter: 1597 . A - the resulting matrix 1598 1599 Level: advanced 1600 1601 Notes: See MATIS for more details. 1602 m and n are NOT related to the size of the map; they are the size of the part of the vector owned 1603 by that process. The sizes of rmap and cmap define the size of the local matrices. 1604 If either rmap or cmap are NULL, then the matrix is assumed to be square. 1605 1606 .seealso: MATIS, MatSetLocalToGlobalMapping() 1607 @*/ 1608 PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 1609 { 1610 PetscErrorCode ierr; 1611 1612 PetscFunctionBegin; 1613 if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 1614 ierr = MatCreate(comm,A);CHKERRQ(ierr); 1615 ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 1616 ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1617 ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 1618 if (rmap && cmap) { 1619 ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 1620 } else if (!rmap) { 1621 ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 1622 } else { 1623 ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 1624 } 1625 PetscFunctionReturn(0); 1626 } 1627 1628 /*MC 1629 MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition methods (e.g. PCBDDC or KSPFETIDP). 1630 This stores the matrices in globally unassembled form. Each processor 1631 assembles only its local Neumann problem and the parallel matrix vector 1632 product is handled "implicitly". 1633 1634 Operations Provided: 1635 + MatMult() 1636 . MatMultAdd() 1637 . MatMultTranspose() 1638 . MatMultTransposeAdd() 1639 . MatZeroEntries() 1640 . MatSetOption() 1641 . MatZeroRows() 1642 . MatSetValues() 1643 . MatSetValuesBlocked() 1644 . MatSetValuesLocal() 1645 . MatSetValuesBlockedLocal() 1646 . MatScale() 1647 . MatGetDiagonal() 1648 . MatMissingDiagonal() 1649 . MatDuplicate() 1650 . MatCopy() 1651 . MatAXPY() 1652 . MatGetSubMatrix() 1653 . MatGetLocalSubMatrix() 1654 - MatSetLocalToGlobalMapping() 1655 1656 Options Database Keys: 1657 . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1658 1659 Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1660 1661 You must call MatSetLocalToGlobalMapping() before using this matrix type. 1662 1663 You can do matrix preallocation on the local matrix after you obtain it with 1664 MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1665 1666 Level: advanced 1667 1668 .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), MatCreateIS(), PCBDDC, KSPFETIDP 1669 1670 M*/ 1671 1672 #undef __FUNCT__ 1673 #define __FUNCT__ "MatCreate_IS" 1674 PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1675 { 1676 PetscErrorCode ierr; 1677 Mat_IS *b; 1678 1679 PetscFunctionBegin; 1680 ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1681 A->data = (void*)b; 1682 1683 /* matrix ops */ 1684 ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1685 A->ops->mult = MatMult_IS; 1686 A->ops->multadd = MatMultAdd_IS; 1687 A->ops->multtranspose = MatMultTranspose_IS; 1688 A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1689 A->ops->destroy = MatDestroy_IS; 1690 A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 1691 A->ops->setvalues = MatSetValues_IS; 1692 A->ops->setvaluesblocked = MatSetValuesBlocked_IS; 1693 A->ops->setvalueslocal = MatSetValuesLocal_IS; 1694 A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 1695 A->ops->zerorows = MatZeroRows_IS; 1696 A->ops->zerorowscolumns = MatZeroRowsColumns_IS; 1697 A->ops->assemblybegin = MatAssemblyBegin_IS; 1698 A->ops->assemblyend = MatAssemblyEnd_IS; 1699 A->ops->view = MatView_IS; 1700 A->ops->zeroentries = MatZeroEntries_IS; 1701 A->ops->scale = MatScale_IS; 1702 A->ops->getdiagonal = MatGetDiagonal_IS; 1703 A->ops->setoption = MatSetOption_IS; 1704 A->ops->ishermitian = MatIsHermitian_IS; 1705 A->ops->issymmetric = MatIsSymmetric_IS; 1706 A->ops->isstructurallysymmetric = MatIsStructurallySymmetric_IS; 1707 A->ops->duplicate = MatDuplicate_IS; 1708 A->ops->missingdiagonal = MatMissingDiagonal_IS; 1709 A->ops->copy = MatCopy_IS; 1710 A->ops->getlocalsubmatrix = MatGetLocalSubMatrix_IS; 1711 A->ops->getsubmatrix = MatGetSubMatrix_IS; 1712 A->ops->axpy = MatAXPY_IS; 1713 A->ops->diagonalset = MatDiagonalSet_IS; 1714 A->ops->shift = MatShift_IS; 1715 1716 /* special MATIS functions */ 1717 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1718 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1719 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 1720 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 1721 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetUpSF_C",MatISSetUpSF_IS);CHKERRQ(ierr); 1722 ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1723 PetscFunctionReturn(0); 1724 } 1725